# DMPlexComputeCellGeometryFEM Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell ## Synopsis ``` #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 ## Input Parameters - ***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 ## Output Parameters - ***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 ## Fortran Notes Since it returns arrays, this routine is only available in Fortran 90, and you must include petsc.h90 in your code. ## See Also `DMGetCoordinateSection()`, `DMGetCoordinates()` ## Level advanced ## Location src/dm/impls/plex/plexgeometry.c ## Examples src/tao/tutorials/ex3.c.html
## Implementations DMPlexComputeCellGeometryFEM_Implicit in src/dm/impls/plex/plexgeometry.c
DMPlexComputeCellGeometryFEM_FE in src/dm/impls/plex/plexgeometry.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/impls/plex/plexgeometry.c) [Index of all DMPlex routines](index.md) [Table of Contents for all manual pages](/docs/manualpages/index.md) [Index of all manual pages](/docs/manualpages/singleindex.md)