Actual source code: daghost.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

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

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

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

 12:    Not Collective

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

 17:    Output Parameters:
 18: +  x,y,z - the corner indices (where y and z are optional; these are used
 19:            for 2D and 3D problems)
 20: -  m,n,p - widths in the corresponding directions (where n and p are optional;
 21:            these are used for 2D and 3D problems)

 23:    Level: beginner

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

 32: .keywords: distributed array, get, ghost, corners, nodes, local indices

 34: .seealso: DMDAGetCorners(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDAGetOwnershipRanges()

 36: @*/
 37: PetscErrorCode  DMDAGetGhostCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p)
 38: {
 39:   PetscInt w;
 40:   DM_DA    *dd = (DM_DA*)da->data;

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