DMPlexTransferVecTree#
transfer a vector between two meshes that differ from each other by refinement/coarsening that can be represented by a common reference tree used by both. This routine can be used for a combination of coarsening and refinement at the same time.
Synopsis#
Collective
Input Parameters#
dmIn - The
DMPLEX
mesh for the input vectordmOut - The second
DMPLEX
meshvecIn - The input vector
sfRefine - A star forest indicating points in the mesh
dmIn
(roots in the star forest) that are parents to points in the meshdmOut
(leaves in the star forest), i.e. wheredmOut
is more refined thandmIn
sfCoarsen - A star forest indicating points in the mesh
dmOut
(roots in the star forest) that are parents to points in the meshdmIn
(leaves in the star forest), i.e. wheredmOut
is more coarsened thandmIn
cidsRefine - The childIds of the points in
dmOut
. These childIds relate back to the reference tree: childid[j] = k implies that mesh point j ofdmOut
was refined from a point indmIn
just as the mesh point k in the reference tree was refined from its parent. childid[j] = -1 indicates that the point j indmOut
is exactly equivalent to its root indmIn
, so no interpolation is necessary. childid[j] = -2 indicates that this point j indmOut
is not a leaf ofsfRefine
.cidsCoarsen - The childIds of the points in
dmIn
. These childIds relate back to the reference tree: childid[j] = k implies that mesh point j of dmIn coarsens to a point indmOut
just as the mesh point k in the reference tree coarsens to its parent. childid[j] = -2 indicates that point j indmOut
is not a leaf insfCoarsen
.useBCs -
PETSC_TRUE
indicates that boundary values should be inserted intovecIn
before transfer.time - Used if boundary values are time dependent.
Output Parameter#
vecOut - Using interpolation and injection operators calculated on the reference tree, the transferred projection of
vecIn
fromdmIn
todmOut
. Note that any field discretized with aPetscFV
finite volume method that uses gradient reconstruction will use reconstructed gradients when interpolating from coarse points to fine points.
See Also#
DMPlex: Unstructured Grids, DM
, DMPLEX
, PetscSF
, Vec
, PetscFV
, DMPlexSetReferenceTree()
, DMPlexGetReferenceTree()
, PetscFVGetComputeGradients()
Level#
developer
Location#
Implementations#
DMPlexTransferVecTree_Interpolate() in src/dm/impls/plex/plextree.c
DMPlexTransferVecTree_Inject() in src/dm/impls/plex/plextree.c
Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages