DMPlexVecGetClosure#
Get an array of the values on the closure of ‘point’
Synopsis#
#include "petscdmplex.h"
PetscErrorCode DMPlexVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
Not collective
Input Parameters#
Input/Output Parameters#
csize - The size of the input values array, or
NULL; on output the number of values in the closurevalues - An array to use for the values, or *values =
NULLto have it allocated automatically; if the user providedNULL, it is a borrowed array and should not be freed, useDMPlexVecRestoreClosure()to return it
Notes#
This is used for getting the all values in a Vec in the closure of a mesh point.
To get only the values in the closure of a mesh point at a specific depth (for example, at mesh vertices), use DMPlexVecGetClosureAtDepth().
DMPlexVecGetClosure()/DMPlexVecRestoreClosure() only allocates the values array if it set to NULL in the
calling function. This is because DMPlexVecGetClosure() is typically called in the inner loop of a Vec or Mat
assembly function, and a user may already have allocated storage for this operation.
A typical use could be
values = NULL;
PetscCall(DMPlexVecGetClosure(dm, NULL, v, p, &clSize, &values));
for (cl = 0; cl < clSize; ++cl) {
<Compute on closure>
}
PetscCall(DMPlexVecRestoreClosure(dm, NULL, v, p, &clSize, &values));
or
PetscMalloc1(clMaxSize, &values);
for (p = pStart; p < pEnd; ++p) {
clSize = clMaxSize;
PetscCall(DMPlexVecGetClosure(dm, NULL, v, p, &clSize, &values));
for (cl = 0; cl < clSize; ++cl) {
<Compute on closure>
}
}
PetscFree(values);
Fortran Notes#
The csize argument is present in the Fortran binding. Since the Fortran values array contains its length information this argument may not be needed.
In that case one may pass PETSC_NULL_INTEGER for csize.
values must be declared with
PetscScalar,dimension(:),pointer :: values
and it will be allocated internally by PETSc to hold the values returned
See Also#
DMPlex: Unstructured Grids, DM, DMPLEX, DMPlexVecGetClosureAtDepth(), DMPlexVecRestoreClosure(), DMPlexVecSetClosure(), DMPlexMatSetClosure()
Level#
intermediate
Location#
Examples#
src/snes/tutorials/ex56.c
src/snes/tutorials/ex77.c
src/dm/impls/plex/tutorials/ex6.c
src/dm/impls/plex/tutorials/ex11.c
Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages