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#

#include "petscdmplex.h"   
PetscErrorCode DMPlexTransferVecTree(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscSF sfRefine, PetscSF sfCoarsen, PetscInt *cidsRefine, PetscInt *cidsCoarsen, PetscBool useBCs, PetscReal time)

Collective

Input Parameters#

  • dmIn - The DMPLEX mesh for the input vector

  • dmOut - The second DMPLEX mesh

  • vecIn - 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 mesh dmOut (leaves in the star forest), i.e. where dmOut is more refined than dmIn

  • sfCoarsen - A star forest indicating points in the mesh dmOut (roots in the star forest) that are parents to points in the mesh dmIn (leaves in the star forest), i.e. where dmOut is more coarsened than dmIn

  • cidsRefine - The childIds of the points in dmOut. These childIds relate back to the reference tree: childid[j] = k implies that mesh point j of dmOut was refined from a point in dmIn just as the mesh point k in the reference tree was refined from its parent. childid[j] = -1 indicates that the point j in dmOut is exactly equivalent to its root in dmIn, so no interpolation is necessary. childid[j] = -2 indicates that this point j in dmOut is not a leaf of sfRefine.

  • 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 in dmOut just as the mesh point k in the reference tree coarsens to its parent. childid[j] = -2 indicates that point j in dmOut is not a leaf in sfCoarsen.

  • useBCs - PETSC_TRUE indicates that boundary values should be inserted into vecIn 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 from dmIn to dmOut. Note that any field discretized with a PetscFV 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#

src/dm/impls/plex/plextree.c

Implementations#

DMPlexTransferVecTree_Interpolate in src/dm/impls/plex/plextree.c
DMPlexTransferVecTree_Inject in src/dm/impls/plex/plextree.c


Edit on GitLab

Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages