DMPlexGetFaceGeometry#

Retrieve the geometric values for a chunk of faces

Synopsis#

#include "petscdmplex.h"   
PetscErrorCode DMPlexGetFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscInt *Nface, PetscFVFaceGeom **fgeom, PetscReal **vol)

Input Parameters#

  • dm - The DM

  • fStart - The first face to include

  • fEnd - The first face to exclude

  • faceGeometry - A local vector with face geometry

  • cellGeometry - A local vector with cell geometry

Output Parameters#

  • Nface - The number of faces with field values

  • fgeom - The extract the face centroid and normal

  • vol - The cell volume

See Also#

DMPlex: Unstructured Grids, DM, DMPLEX, DMPlexGetCellFields()

Level#

developer

Location#

src/dm/impls/plex/plexfem.c


Edit on GitLab

Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages