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
NULL
to have it allocated automatically; if the user providedNULL
, it is a borrowed array and should not be freed
Notes#
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 not present in the Fortran binding since it is internal to the array.
See Also#
DMPlex: Unstructured Grids, DM
, DMPLEX
, DMPlexVecRestoreClosure()
, DMPlexVecSetClosure()
, DMPlexMatSetClosure()
Level#
intermediate
Location#
Examples#
src/snes/tutorials/ex77.c
src/dm/impls/plex/tutorials/ex11.c
src/dm/impls/plex/tutorials/ex6.c
src/snes/tutorials/ex56.c
Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages