DMPlexVecSetClosure#
Set an array of the values on the closure of point
Synopsis#
#include "petscdmplex.h"
PetscErrorCode DMPlexVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt point, const PetscScalar values[], InsertMode mode)
Not collective
Input Parameters#
dm - The
DM
section - The section describing the layout in
v
, orNULL
to use the default sectionv - The local vector
point - The point in the
DM
values - The array of values
mode - The insert mode. One of
INSERT_ALL_VALUES
,ADD_ALL_VALUES
,INSERT_VALUES
,ADD_VALUES
,INSERT_BC_VALUES
, andADD_BC_VALUES
, whereINSERT_ALL_VALUES
andADD_ALL_VALUES
also overwrite boundary conditions.
Note#
Usually the input arrays were obtained with DMPlexVecGetClosure()
Fortran Note#
values
must be declared with
PetscScalar,dimension(:),pointer :: values
See Also#
DMPlex: Unstructured Grids, DM
, DMPLEX
, DMPlexVecGetClosure()
, DMPlexMatSetClosure()
Level#
intermediate
Location#
Examples#
src/dm/impls/plex/tutorials/ex6.c
Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages