petsc-3.12.5 2020-03-29
DMPlexComputeCellGeometryAffineFEM
Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
Synopsis
#include "petscdmplex.h"
#include "petscfe.h"
PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
Collective on dm
Input Arguments
| dm | - the DM
|
| cell | - the cell
|
Output Arguments
| v0 | - the translation part of this affine transform
|
| J | - the Jacobian of the transform from the reference element
|
| invJ | - the inverse of the Jacobian
|
| detJ | - the Jacobian determinant
|
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
DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinates()
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