2: /*
3: Code for manipulating distributed regular arrays in parallel.
4: */
6: #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/
10: /*@
11: DMDAGetGhostCorners - Returns the global (x,y,z) indices of the lower left
12: corner and size of the local region, including ghost points.
14: Not Collective
16: Input Parameter:
17: . da - the distributed array
19: Output Parameters:
20: + x,y,z - the corner indices (where y and z are optional; these are used
21: for 2D and 3D problems)
22: - m,n,p - widths in the corresponding directions (where n and p are optional;
23: these are used for 2D and 3D problems)
25: Level: beginner
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.
34: .keywords: distributed array, get, ghost, corners, nodes, local indices
36: .seealso: DMDAGetCorners(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDAGetOwnershipRanges()
38: @*/
39: PetscErrorCodeDMDAGetGhostCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p) 40: {
41: PetscInt w;
42: DM_DA *dd = (DM_DA*)da->data;
46: /* since the xs, xe ... have all been multiplied by the number of degrees
47: of freedom per cell, w = dd->w, we divide that out before returning.*/
48: w = dd->w;
49: if (x) *x = dd->Xs/w + dd->xo;
50: if (y) *y = dd->Ys + dd->yo;
51: if (z) *z = dd->Zs + dd->zo;
52: if (m) *m = (dd->Xe - dd->Xs)/w;
53: if (n) *n = (dd->Ye - dd->Ys);
54: if (p) *p = (dd->Ze - dd->Zs);
55: return(0);
56: }