DMDAGetCorners#

Returns the global (x,y,z) indices of the lower left corner and size of the local region, excluding ghost points.

Synopsis#

#include "petscdmda.h"   
PetscErrorCode DMDAGetCorners(DM da, PetscInt *x, PetscInt *y, PetscInt *z, PetscInt *m, PetscInt *n, PetscInt *p)

Not Collective

Input Parameter#

Output Parameters#

  • x - the corner index for the first dimension

  • y - the corner index for the second dimension (only used in 2D and 3D problems)

  • z - the corner index for the third dimension (only used in 3D problems)

  • m - the width in the first dimension

  • n - the width in the second dimension (only used in 2D and 3D problems)

  • p - the width in the third dimension (only used in 3D problems)

Notes#

Any of y, z, n, and p can be passed in as NULL if not needed.

The corner information is independent of the number of degrees of freedom per node set with the DMDACreateXX() routine. Thus the x, y, and z can be thought of as the lower left coordinates of the patch of values on process on a logical grid and m, n, and p as the extent of the patch, where each grid point has (potentially) several degrees of freedom.

See Also#

DMDA - Creating vectors for structured grids, DM, DMDA, DMDAGetGhostCorners(), DMDAGetOwnershipRanges(), DMStagGetCorners(), DMSTAG

Level#

beginner

Location#

src/dm/impls/da/dacorn.c

Examples#

src/tao/unconstrained/tutorials/eptorsion2f.F90
src/tao/unconstrained/tutorials/eptorsion2.c
src/tao/bound/tutorials/jbearing2.c
src/tao/unconstrained/tutorials/minsurf2.c
src/tao/bound/tutorials/plate2f.F90
src/tao/unconstrained/tutorials/spectraladjointassimilation.c
src/tao/unconstrained/tutorials/burgers_spectral.c
src/tao/bound/tutorials/plate2.c
src/tao/complementarity/tutorials/blackscholes.c
src/dm/tutorials/ex11f90.F90


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