xref: /petsc/src/dm/impls/da/daghost.c (revision 09573ac72a50d3e7ecd55a2b7f0ef28450cd0a8b)
1 #define PETSCDM_DLL
2 
3 /*
4   Code for manipulating distributed regular arrays in parallel.
5 */
6 
7 #include "private/daimpl.h"    /*I   "petscdm.h"   I*/
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "DMDAGetGhostCorners"
11 /*@
12    DMDAGetGhostCorners - Returns the global (x,y,z) indices of the lower left
13    corner of the local region, including ghost points.
14 
15    Not Collective
16 
17    Input Parameter:
18 .  da - the distributed array
19 
20    Output Parameters:
21 +  x,y,z - the corner indices (where y and z are optional; these are used
22            for 2D and 3D problems)
23 -  m,n,p - widths in the corresponding directions (where n and p are optional;
24            these are used for 2D and 3D problems)
25 
26    Level: beginner
27 
28    Note:
29    The corner information is independent of the number of degrees of
30    freedom per node set with the DMDACreateXX() routine. Thus the x, y, z, and
31    m, n, p can be thought of as coordinates on a logical grid, where each
32    grid point has (potentially) several degrees of freedom.
33    Any of y, z, n, and p can be passed in as PETSC_NULL if not needed.
34 
35 .keywords: distributed array, get, ghost, corners, nodes, local indices
36 
37 .seealso: DMDAGetCorners(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDAGetOwnershipRanges()
38 
39 @*/
40 PetscErrorCode PETSCDM_DLLEXPORT DMDAGetGhostCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p)
41 {
42   PetscInt w;
43   DM_DA    *dd = (DM_DA*)da->data;
44 
45   PetscFunctionBegin;
46   PetscValidHeaderSpecific(da,DM_CLASSID,1);
47   /* since the xs, xe ... have all been multiplied by the number of degrees
48      of freedom per cell, w = dd->w, we divide that out before returning.*/
49   w = dd->w;
50   if (x) *x = dd->Xs/w; if (m) *m = (dd->Xe - dd->Xs)/w;
51   if (y) *y = dd->Ys;   if (n) *n = (dd->Ye - dd->Ys);
52   if (z) *z = dd->Zs;   if (p) *p = (dd->Ze - dd->Zs);
53   PetscFunctionReturn(0);
54 }
55 
56