xref: /petsc/src/dm/impls/da/daghost.c (revision 63f3c55c12ae2f190c391f6df6c540efe018a44a)
1 
2 /*
3   Code for manipulating distributed regular arrays in parallel.
4 */
5 
6 #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
7 
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.
11 
12    Not Collective
13 
14    Input Parameter:
15 .  da - the distributed array
16 
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)
22 
23    Level: beginner
24 
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.
31 
32 .seealso: DMDAGetCorners(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDAGetOwnershipRanges(), DMStagGetGhostCorners()
33 
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;
39 
40   PetscFunctionBegin;
41   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
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(0);
52 }
53 
54