petsc-3.14.6 2021-03-30
Report Typos and Errors

DMPlexGetConeRecursive

Expand each given point into its cone points and do that recursively until we end up just with vertices (DAG points of depth 0, i.e. without cones).

Synopsis

#include "petscdmplex.h"   
PetscErrorCode DMPlexGetConeRecursive(DM dm, IS points, PetscInt *depth, IS *expandedPoints[], PetscSection *sections[])
Not collective

Input Parameters

dm - The DMPlex
points - The IS of points, which must lie in the chart set with DMPlexSetChart()

Output Parameter

depth - (optional) Size of the output arrays, equal to DMPlex depth, returned by DMPlexGetDepth()
expandedPoints - (optional) An array of index sets with recursively expanded cones
sections - (optional) An array of sections which describe mappings from points to their cone points

Notes

Like DMPlexGetConeTuple() but recursive.

Array expandedPoints has size equal to depth. Each expandedPoints[d] contains DAG points with maximum depth d, recursively cone-wise expanded from the input points. For example, for d=0 it contains only vertices, for d=1 it can contain vertices and edges, etc.

Array section has size equal to depth. Each PetscSection sections[d] realizes mapping from expandedPoints[d+1] (section points) to expandedPoints[d] (section dofs) as follows

(1) DAG points in expandedPoints[d+1] with depth d+1 to their cone points in expandedPoints[d]; (2) DAG points in expandedPoints[d+1] with depth in [0,d] to the same points in expandedPoints[d].

See Also

DMPlexCreate(), DMPlexGetCone(), DMPlexGetConeTuple(), DMPlexRestoreConeRecursive(), DMPlexGetConeRecursiveVertices(), DMPlexGetDepth()

Level

advanced

Location

src/dm/impls/plex/plex.c
Index of all DMPLEX routines
Table of Contents for all manual pages
Index of all manual pages