Actual source code: daghost.c

  1: /*
  2:   Code for manipulating distributed regular arrays in parallel.
  3: */

  5: #include <petsc/private/dmdaimpl.h>

  7: /*@C
  8:   DMDAGetGhostCorners - Returns the global (x,y,z) indices of the lower left
  9:   corner and size of the local region, including ghost points.

 11:   Not Collective

 13:   Input Parameter:
 14: . da - the distributed array

 16:   Output Parameters:
 17: + x - the corner index for the first dimension
 18: . y - the corner index for the second dimension (only used in 2D and 3D problems)
 19: . z - the corner index for the third dimension (only used in 3D problems)
 20: . m - the width in the first dimension
 21: . n - the width in the second dimension (only used in 2D and 3D problems)
 22: - p - the width in the third dimension (only used in 3D problems)

 24:   Level: beginner

 26:   Note:
 27:   The corner information is independent of the number of degrees of
 28:   freedom per node set with the `DMDACreateXX()` routine. Thus the `x`, `y`, `z`, and
 29:   `m`, `n`, `p` can be thought of as coordinates on a logical grid, where each
 30:   grid point has (potentially) several degrees of freedom.
 31:   Any of `y`, `z`, `n`, and `p` can be passed in as `NULL` if not needed.

 33: .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetCorners()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDAGetOwnershipRanges()`, `DMStagGetGhostCorners()`
 34: @*/
 35: PetscErrorCode DMDAGetGhostCorners(DM da, PetscInt *x, PetscInt *y, PetscInt *z, PetscInt *m, PetscInt *n, PetscInt *p)
 36: {
 37:   PetscInt w;
 38:   DM_DA   *dd = (DM_DA *)da->data;

 40:   PetscFunctionBegin;
 42:   /* since the xs, xe ... have all been multiplied by the number of degrees
 43:      of freedom per cell, w = dd->w, we divide that out before returning.*/
 44:   w = dd->w;
 45:   if (x) *x = dd->Xs / w + dd->xo;
 46:   if (y) *y = dd->Ys + dd->yo;
 47:   if (z) *z = dd->Zs + dd->zo;
 48:   if (m) *m = (dd->Xe - dd->Xs) / w;
 49:   if (n) *n = (dd->Ye - dd->Ys);
 50:   if (p) *p = (dd->Ze - dd->Zs);
 51:   PetscFunctionReturn(PETSC_SUCCESS);
 52: }