petsc-3.13.6 2020-09-29
Report Typos and Errors

DMPlexGetClosureIndices

Get the global indices for all local points in the closure of the given point

Synopsis

#include "petscdmplex.h"   
PetscErrorCode DMPlexGetClosureIndices(DM dm, PetscSection section, PetscSection idxSection, PetscInt point, PetscInt *numIndices, PetscInt **indices, PetscInt *outOffsets)
Not collective

Input Parameters

dm - The DM
section - The section describing the points (a local section)
idxSection - The section on which to obtain indices (may be local or global)
point - The mesh point

Output parameters

numIndices - The number of indices
indices - The indices
outOffsets - Field offset if not NULL

Notes

Must call DMPlexRestoreClosureIndices() to free allocated memory

If idxSection is global, any constrained dofs (see DMAddBoundary(), for example) will get negative indices. The value of those indices is not significant. If idxSection is local, the constrained dofs will yield the involution -(idx+1) of their index in a local vector. A caller who does not wish to distinguish those points may recover the nonnegative indices via involution, -(-(idx+1)+1)==idx. Local indices are provided when idxSection == section, otherwise global indices (with the above semantics) are implied.

See Also

DMPlexRestoreClosureIndices(), DMPlexVecGetClosure(), DMPlexMatSetClosure(), DMGetLocalSection(), DMGetGlobalSection()

Level

advanced

Location

src/dm/impls/plex/plex.c

Examples

src/dm/impls/plex/tutorials/ex8.c.html

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