petsc-3.9.4 2018-09-11
Report Typos and Errors

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
csize - The size of the input values array, or NULL
values - An array to use for the values, or NULL to have it allocated automatically

Output Parameters

csize - The number of values in the closure
values - The array of values. If the user provided NULL, it is a borrowed array and should not be freed

Note that 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;
 ierr = DMPlexVecGetClosure(dm, NULL, v, p, &clSize, &values);CHKERRQ(ierr);
 for (cl = 0; cl < clSize; ++cl) {
   <Compute on closure>
 }
 ierr = DMPlexVecRestoreClosure(dm, NULL, v, p, &clSize, &values);CHKERRQ(ierr);

or

 PetscMalloc1(clMaxSize, &values);
 for (p = pStart; p < pEnd; ++p) {
   clSize = clMaxSize;
   ierr = DMPlexVecGetClosure(dm, NULL, v, p, &clSize, &values);CHKERRQ(ierr);
   for (cl = 0; cl < clSize; ++cl) {
     <Compute on closure>
   }
 }
 PetscFree(values);

Fortran Notes

Since it returns an array, this routine is only available in Fortran 90, and you must include petsc.h90 in your code.

The csize argument is not present in the Fortran 90 binding since it is internal to the array.

See Also

DMPlexVecRestoreClosure(), DMPlexVecSetClosure(), DMPlexMatSetClosure()

Level

intermediate

Location

src/dm/impls/plex/plex.c

Examples

src/snes/examples/tutorials/ex56.c.html
src/snes/examples/tutorials/ex62.c.html
src/snes/examples/tutorials/ex77.c.html

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