DMPlexRefineRegularGetAffineTransforms#
Gets the affine map from the reference cell to each subcell
Synopsis#
#include "petscdmplextransform.h"
PetscErrorCode DMPlexRefineRegularGetAffineTransforms(DMPlexTransform tr, DMPolytopeType ct, PetscInt *Nc, PetscReal *v0[], PetscReal *J[], PetscReal *invJ[])
Input Parameters#
tr - The
DMPlexTransform
objectct - The cell type
Output Parameters#
Nc - The number of subcells produced from this cell type
v0 - The translation of the first vertex for each subcell, an array of length \(dim * Nc\). Pass
NULL
to ignore.J - The Jacobian for each subcell (map from reference cell to subcell), an array of length \(dim^2 * Nc\). Pass
NULL
to ignore.invJ - The inverse Jacobian for each subcell, an array of length \(dim^2 * Nc\). Pass
NULL
to ignore.
Note#
Do not free these output arrays
See Also#
DMPLEX
, DM
, DMPlexTransform
, DMPolytopeType
, DMPlexRefineRegularGetAffineFaceTransforms()
, DMPLEXREFINEREGULAR
Level#
developer
Location#
src/dm/impls/plex/transform/impls/refine/regular/plexrefregular.c
Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages