#include "petscdmplex.h" #include "petscfe.h" PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ)Collective on DM
| dm | - the DM | |
| cell | - the cell | |
| quad | - the quadrature containing the points in the reference element where the geometry will be evaluated. If quad == NULL, geometry will be evaluated at the first vertex of the reference element |
| v | - the image of the transformed quadrature points, otherwise the image of the first vertex in the closure of the reference element | |
| J | - the Jacobian of the transform from the reference element at each quadrature point | |
| invJ | - the inverse of the Jacobian at each quadrature point | |
| detJ | - the Jacobian determinant at each quadrature point |
Level:advanced
Location:src/dm/impls/plex/plexgeometry.c
Index of all DMPLEX routines
Table of Contents for all manual pages
Index of all manual pages