PetscFECompositeGetMapping#
Returns the mappings from the reference element to each subelement
Synopsis#
#include "petscfe.h"
#include "petscdt.h"
PetscErrorCode PetscFECompositeGetMapping(PetscFE fem, PetscInt *numSubelements, const PetscReal *v0[], const PetscReal *jac[], const PetscReal *invjac[])
Not Collective
Input Parameter#
fem - The
PetscFE
object
Output Parameters#
numSubelements - The number of sub elements
v0 - The affine transformation for each element, an array of length \(dim * Nc\). Pass
NULL
to ignore.jac - The Jacobian for each element, an array of length \(dim^2 * Nc\). Pass
NULL
to ignore.invjac - The inverse of the Jacobian, an array of length \(dim^2 * Nc\). Pass
NULL
to ignore.
Note#
Do not free the output arrays.
See Also#
Level#
intermediate
Location#
src/dm/dt/fe/impls/composite/fecomposite.c
Index of all FE routines
Table of Contents for all manual pages
Index of all manual pages