:orphan: # 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 - ***dm -*** The `DM` - ***section -*** The section describing the layout in `v`, or `NULL` to use the default section - ***v -*** The local vector - ***point -*** The point in the `DM` ## Input/Output Parameters - ***csize -*** The size of the input values array, or `NULL`; on output the number of values in the closure - ***values -*** An array to use for the values, or `NULL` to have it allocated automatically; if the user provided `NULL`, 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 ```none values = NULL; PetscCall(DMPlexVecGetClosure(dm, NULL, v, p, &clSize, &values)); for (cl = 0; cl < clSize; ++cl) { } PetscCall(DMPlexVecRestoreClosure(dm, NULL, v, p, &clSize, &values)); ``` or ```none 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) { } } PetscFree(values); ``` ## Fortran Note The `csize` argument is not present in the Fortran binding since it is internal to the array. ## See Also [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexVecRestoreClosure()`, `DMPlexVecSetClosure()`, `DMPlexMatSetClosure()` ## Level intermediate ## Location src/dm/impls/plex/plex.c ## Examples src/dm/impls/plex/tutorials/ex11.c
src/dm/impls/plex/tutorials/ex6.c
src/snes/tutorials/ex56.c
src/snes/tutorials/ex77.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/impls/plex/plex.c) [Index of all DMPlex routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)