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

DMPlexGetClosureIndices

Gets the global dof indices associated with the closure of the given point within the provided sections.

Synopsis

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

Input Parameters

dm - The DM
section - The PetscSection describing the points (a local section)
idxSection - The PetscSection from which to obtain indices (may be local or global)
point - The point defining the closure
useClPerm - Use the closure point permutation if available

Output Parameters

numIndices - The number of dof indices in the closure of point with the input sections
indices - The dof indices
outOffsets - Array to write the field offsets into, or NULL
values - The input values, which may be modified if sign flips are induced by the point symmetries, or 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