PetscFEGeomGetPoint#
Get the geometry for cell c
at point p
as a PetscFEGeom
Synopsis#
#include "petscfe.h"
PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *geom, PetscInt c, PetscInt p, const PetscReal pcoords[], PetscFEGeom *pgeom)
Input Parameters#
geom -
PetscFEGeom
objectc - The cell
p - The point
pcoords - The reference coordinates of point
p
, orNULL
Output Parameter#
pgeom - The geometry of cell
c
at pointp
Notes#
For affine geometries, this only copies to pgeom
at point 0. Since we copy pointers into pgeom
,
nothing needs to be done with it afterwards.
In the affine case, pgeom
must have storage for the integration point coordinates in pgeom->v if pcoords
is passed in.
See Also#
Level#
intermediate
Location#
src/dm/dt/fe/interface/fegeom.c
Index of all FE routines
Table of Contents for all manual pages
Index of all manual pages