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; No Fortran Support

Input Parameter#

  • dm - the DMDA 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.

See Also#

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

Level#

intermediate

Location#

src/dm/impls/da/dagetelem.c

Examples#

src/dm/tutorials/ex5.c
src/ksp/ksp/tutorials/ex70.c
src/ksp/ksp/tutorials/ex71.c

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


Edit on GitLab

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