petsc-3.11.4 2019-09-28
Report Typos and Errors

DMDAGetElements

Gets an array containing the indices (in local coordinates) of all the local elements

Synopsis

#include "petscdmda.h"   
PetscErrorCode  DMDAGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
Not Collective

Input Parameter

dm -the DM object

Output Parameters

nel - number of local elements
nen - number of element nodes
e - the local indices of the elements' vertices

Notes

Call DMDARestoreElements() once you have finished accessing the elements.

Each process uniquely owns a subset of the elements. That is no element is owned by two or more processes.

If on each process you integrate over its owned elements and use ADD_VALUES in Vec/MatSetValuesLocal() then you'll obtain the correct result.

Not supported in Fortran

See Also

DMDAElementType, DMDASetElementType(), VecSetValuesLocal(), MatSetValuesLocal(), DMGlobalToLocalBegin(), DMLocalToGlobalBegin()

Level

intermediate

Location

src/dm/impls/da/dagetelem.c

Examples

src/dm/examples/tutorials/ex5.c.html
src/ksp/ksp/examples/tutorials/ex70.c.html

Implementations

DMDAGetElements_1D in src/dm/impls/da/dagetelem.c
DMDAGetElements_2D in src/dm/impls/da/dagetelem.c
DMDAGetElements_3D in src/dm/impls/da/dagetelem.c

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