petsc-3.7.7 2017-09-25
Report Typos and Errors

DMPlexComputeCellGeometryFEM

Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell

Synopsis

#include "petscdmplex.h"   
PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscFE fe, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
Collective on DM

Input Arguments

dm - the DM
cell - the cell
fe - the finite element containing the quadrature

Output Arguments

v0 - the translation part of this transform
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(), DMGetCoordinateVec()

Level:advanced
Location:
src/dm/impls/plex/plexgeometry.c
Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages