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 - the corner index for the first dimension 19 . y - the corner index for the second dimension (only used in 2D and 3D problems) 20 . z - the corner index for the third dimension (only used in 3D problems) 21 . m - the width in the first dimension 22 . n - the width in the second dimension (only used in 2D and 3D problems) 23 - p - the width in the third dimension (only used in 3D problems) 24 25 Level: beginner 26 27 Note: 28 The corner information is independent of the number of degrees of 29 freedom per node set with the DMDACreateXX() routine. Thus the x, y, z, and 30 m, n, p can be thought of as coordinates on a logical grid, where each 31 grid point has (potentially) several degrees of freedom. 32 Any of y, z, n, and p can be passed in as NULL if not needed. 33 34 .seealso: `DMDAGetCorners()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDAGetOwnershipRanges()`, `DMStagGetGhostCorners()` 35 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; 41 42 PetscFunctionBegin; 43 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 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 PetscFunctionReturn(0); 54 } 55 56