[petsc-users] Visualizing higher order finite element output in ParaView

Matthew Knepley knepley at gmail.com
Thu Jan 30 09:19:13 CST 2025


On Thu, Jan 30, 2025 at 9:43 AM Anna Dalklint <anna.dalklint at solid.lth.se>
wrote:

> I looked deeper into the petsc codebase regarding HDF5. From what I
> understood (which of course can be wrong), the current version of petsc
> does not save edge degrees-of-freedom to HDF5? Is this something you plan
> to allow?
>

We write two different outputs (by default). One has all the data, and one
has only cell and vertex data because Paraview does not understand anything
else. This can be customized with options. What do you want to save?


> Otherwise I’m fine with using CGNS. But could you please explain how I
> could save timeseries that paraview recognizes using this format? Right now
> I’m saving files e.g. file0001.cgns, file0002.cgns, … where each .cgns file
> is written using VecView (i.e. it stores a discretized field). But paraview
> cannot load this as a timeseries.
>

Jed can explain how this works.


> Also, do you have any documentation regarding node (vertex, edge, face,
> cell) numbering? E.g. how would a 10 node tetrahedral be numbered? From the
> documentation on your webpage (https://urldefense.us/v3/__https://petsc.org/release/manual/dmplex/__;!!G_uCfscf7eWS!ejLt_rGcvE6uln2mYtMTHg6zBUMpobt0KviK7ZOeyB9vp_BRV_c7m_xMnDdSnqjEerY4ApgAANbxBmrzLVYe$ )
> it looks like cell dofs -> vertex dofs-> face dofs-> edge dofs. Is this
> correct?
>

When you call DMPlexVecGetClosure(), the closure follows the point
numbering, in that for each point, we lookup the dofs in the local Section,
and push them into the array in order. So then you need the point ordering.
For the closure, it goes by dimension, so cell dofs, face dofs, edge dofs,
vertex dofs. You can see the definition of faces (and edges) here:


https://urldefense.us/v3/__https://gitlab.com/petsc/petsc/-/blob/main/src/dm/impls/plex/plexinterpolate.c?ref_type=heads*L196__;Iw!!G_uCfscf7eWS!ejLt_rGcvE6uln2mYtMTHg6zBUMpobt0KviK7ZOeyB9vp_BRV_c7m_xMnDdSnqjEerY4ApgAANbxBonw9GhL$ 

and triangles are ordered here


https://urldefense.us/v3/__https://gitlab.com/petsc/petsc/-/blob/main/src/dm/impls/plex/plexinterpolate.c?ref_type=heads*L115__;Iw!!G_uCfscf7eWS!ejLt_rGcvE6uln2mYtMTHg6zBUMpobt0KviK7ZOeyB9vp_BRV_c7m_xMnDdSnqjEerY4ApgAANbxBmRfa-Ea$ 

The idea is that DMPlexVecGetClosure() delivers the dofs in a standard
order on the element, so that you can write
your residual function once. Also, for multiple fields, they are stacked
contiguously, so the numbering is [field, point, dof on point].

Let me know if that does not make sense.

  Thanks,

     Matt


> Thanks,
>
> Anna
>
>
>
> *From: *Matthew Knepley <knepley at gmail.com>
> *Date: *Thursday, 30 January 2025 at 00:39
> *To: *Jed Brown <jed at jedbrown.org>
> *Cc: *Anna Dalklint <anna.dalklint at solid.lth.se>, petsc-users at mcs.anl.gov
> <petsc-users at mcs.anl.gov>
> *Subject: *Re: [petsc-users] Visualizing higher order finite element
> output in ParaView
>
> That is all true. If you want lower level pieces to make it yourself, I
> have -dm_plex_high_order_view, which activates
>
> DMPlexCreateHighOrderSurrogate_Internal(). This is a simple function that
> refines the mesh lg(p) times to try and
>
> resolve the high order behavior.
>
>
>
>   Thanks,
>
>
>
>      Matt
>
>
>
> On Wed, Jan 29, 2025 at 4:55 PM Jed Brown <jed at jedbrown.org> wrote:
>
> I like the CGNS workflow for this, at least with quadratic and cubic
> elements. You can use options like -snes_view_solution cgns:solution.cgns
> (configure with --download-cgns). It can also monitor transient solves with
> flexible batch sizes (geometry and connectivity are stored only once within
> a batch of output frames).
>
> Anna Dalklint via petsc-users <petsc-users at mcs.anl.gov> writes:
>
> > Hello,
> >
> > We have created a finite element code in PETSc for unstructured meshes
> using DMPlex. The first order meshes are created in gmsh and loaded into
> PETSc. To introduce higher order elements, e.g. 10 node tetrahedral
> elements, we start from scratch using PetscSection and loop over the
> relevant points it the DM to introduce additional degrees-of-freedom
> (example; for 10 node tets we have 4 vertices “nodes” and 6 edge “nodes”).
> The coordinates of the new “nodes” are obtained by interpolation using the
> finite element basis functions.
> >
> > The simulations seem to run well, but we face issues when trying to
> visualize the results in ParaView. We have tried to use both CGNS and
> HDF5+XDMF file formats for e.g. VecView. CGNS works, but the edge
> degrees-of-freedom appear to not be interpolated correctly (we observe
> oscillations in the fields, don’t know if this is a PETSc och ParaView
> issue). Also, we would prefer to use another file format than CGNS since it
> does not appear to directly allow timeseries (at least ParaView doesn’t
> recognize it). We haven’t got the HDF5+XDMF file format to work at all when
> running on more than one core (the mesh is highly distorted when saving
> using VecView and DMView + running the “petsc_gen_xdmf.py” script on the
> .h5 output file).
> >
> > VTU format works but then only the vertices’ degrees-of-freedom are
> visualized. As far as we have understood it, this is because VTU/VTK only
> supports degrees-of-freedom on vertices/cell level.
> >
> > Does anyone have any idea of how to visualize fields generated from
> higher order elements in ParaView? Or understand what we might be doing
> wrong?
> >
> > Best regards,
> > Anna
>
>
>
>
> --
>
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
>
>
> https://urldefense.us/v3/__https://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!ejLt_rGcvE6uln2mYtMTHg6zBUMpobt0KviK7ZOeyB9vp_BRV_c7m_xMnDdSnqjEerY4ApgAANbxBqOUos6t$ 
> <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!ejLt_rGcvE6uln2mYtMTHg6zBUMpobt0KviK7ZOeyB9vp_BRV_c7m_xMnDdSnqjEerY4ApgAANbxBuwcaQjh$ >
>


-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener

https://urldefense.us/v3/__https://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!ejLt_rGcvE6uln2mYtMTHg6zBUMpobt0KviK7ZOeyB9vp_BRV_c7m_xMnDdSnqjEerY4ApgAANbxBqOUos6t$  <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!ejLt_rGcvE6uln2mYtMTHg6zBUMpobt0KviK7ZOeyB9vp_BRV_c7m_xMnDdSnqjEerY4ApgAANbxBuwcaQjh$ >
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20250130/def6b125/attachment.html>


More information about the petsc-users mailing list