xref: /petsc/src/dm/impls/da/daghost.c (revision 6d8694c4fbab79f9439f1ad13c0386ba7ee1ca4b)
147c6ae99SBarry Smith /*
247c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
347c6ae99SBarry Smith */
447c6ae99SBarry Smith 
5af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
647c6ae99SBarry Smith 
75d83a8b1SBarry Smith /*@
872fd0fbdSBarry Smith   DMDAGetGhostCorners - Returns the global (`i`,`j`,`k`) indices of the lower left
959f3ab6dSMatthew G. Knepley   corner and size of the local region, including ghost points.
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith   Not Collective
1247c6ae99SBarry Smith 
1347c6ae99SBarry Smith   Input Parameter:
1472fd0fbdSBarry Smith . da - the `DMDA`
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith   Output Parameters:
176b867d5aSJose E. Roman + x - the corner index for the first dimension
186b867d5aSJose E. Roman . y - the corner index for the second dimension (only used in 2D and 3D problems)
196b867d5aSJose E. Roman . z - the corner index for the third dimension (only used in 3D problems)
206b867d5aSJose E. Roman . m - the width in the first dimension
216b867d5aSJose E. Roman . n - the width in the second dimension (only used in 2D and 3D problems)
226b867d5aSJose E. Roman - p - the width in the third dimension (only used in 3D problems)
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   Level: beginner
2547c6ae99SBarry Smith 
26f13dfd9eSBarry Smith   Notes:
27f13dfd9eSBarry Smith   Any of `y`, `z`, `n`, and `p` can be passed in as `NULL` if not needed.
28f13dfd9eSBarry Smith 
2947c6ae99SBarry Smith   The corner information is independent of the number of degrees of
3072fd0fbdSBarry Smith   freedom per node set with the `DMDACreateXX()` routine. Thus the `x`, `y`, and `z`
3172fd0fbdSBarry Smith   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
3272fd0fbdSBarry Smith   extent of the patch. Where
3347c6ae99SBarry Smith   grid point has (potentially) several degrees of freedom.
3447c6ae99SBarry Smith 
3572fd0fbdSBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetCorners()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDAGetOwnershipRanges()`, `DMStagGetGhostCorners()`, `DMSTAG`
3647c6ae99SBarry Smith @*/
DMDAGetGhostCorners(DM da,PeOp PetscInt * x,PeOp PetscInt * y,PeOp PetscInt * z,PeOp PetscInt * m,PeOp PetscInt * n,PeOp PetscInt * p)37*ce78bad3SBarry Smith PetscErrorCode DMDAGetGhostCorners(DM da, PeOp PetscInt *x, PeOp PetscInt *y, PeOp PetscInt *z, PeOp PetscInt *m, PeOp PetscInt *n, PeOp PetscInt *p)
38d71ae5a4SJacob Faibussowitsch {
3947c6ae99SBarry Smith   PetscInt w;
4047c6ae99SBarry Smith   DM_DA   *dd = (DM_DA *)da->data;
4147c6ae99SBarry Smith 
4247c6ae99SBarry Smith   PetscFunctionBegin;
43a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4447c6ae99SBarry Smith   /* since the xs, xe ... have all been multiplied by the number of degrees
4547c6ae99SBarry Smith      of freedom per cell, w = dd->w, we divide that out before returning.*/
4647c6ae99SBarry Smith   w = dd->w;
4759bc5b24SSatish Balay   if (x) *x = dd->Xs / w + dd->xo;
4859bc5b24SSatish Balay   if (y) *y = dd->Ys + dd->yo;
4959bc5b24SSatish Balay   if (z) *z = dd->Zs + dd->zo;
5059bc5b24SSatish Balay   if (m) *m = (dd->Xe - dd->Xs) / w;
5159bc5b24SSatish Balay   if (n) *n = (dd->Ye - dd->Ys);
5259bc5b24SSatish Balay   if (p) *p = (dd->Ze - dd->Zs);
533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5447c6ae99SBarry Smith }
55