xref: /petsc/src/dm/impls/patch/patch.c (revision 6187e55b7be868de3dab288e5aecadaea4c1b258)
1 #include <petsc-private/dmpatchimpl.h>   /*I      "petscdmpatch.h"   I*/
2 #include <petscdmda.h>
3 #include <petscsf.h>
4 
5 /*
6 Solver loop to update \tau:
7 
8   DMZoom(dmc, &dmz)
9   DMRefine(dmz, &dmf),
10   Scatter Xcoarse -> Xzoom,
11   Interpolate Xzoom -> Xfine (note that this may be on subcomms),
12   Smooth Xfine using two-step smoother
13     normal smoother plus Kaczmarz---moves back and forth from dmzoom to dmfine
14   Compute residual Rfine
15   Restrict Rfine to Rzoom_restricted
16   Scatter Rzoom_restricted -> Rcoarse_restricted
17   Compute global residual Rcoarse
18   TauCoarse = Rcoarse - Rcoarse_restricted
19 */
20 
21 #undef __FUNCT__
22 #define __FUNCT__ "DMPatchZoom"
23 /*
24   DMPatchZoom - Create a version of the coarse patch (identified by rank) with halo on communicator commz
25 
26   Collective on DM
27 
28   Input Parameters:
29   + dm - the DM
30   . rank - the rank which holds the given patch
31   - commz - the new communicator for the patch
32 
33   Output Parameters:
34   + dmz  - the patch DM
35   . sfz  - the PetscSF mapping the patch+halo to the zoomed version
36   . sfzr - the PetscSF mapping the patch to the restricted zoomed version
37 
38   Level: intermediate
39 
40   Note: All processes in commz should have the same rank (could autosplit comm)
41 
42 .seealso: DMPatchSolve()
43 */
44 PetscErrorCode DMPatchZoom(DM dm, Vec X, MatStencil lower, MatStencil upper, MPI_Comm commz, DM *dmz, PetscSF *sfz, PetscSF *sfzr)
45 {
46   DMDAStencilType st;
47   MatStencil      blower, bupper, loclower, locupper;
48   IS              is;
49   const PetscInt  *ranges, *indices;
50   PetscInt        *localPoints  = NULL;
51   PetscSFNode     *remotePoints = NULL;
52   PetscInt        dim, dof;
53   PetscInt        M, N, P, rM, rN, rP, halo = 1, sxb, syb, szb, sxr, syr, szr, exr, eyr, ezr, mxb, myb, mzb, i, j, k, q;
54   PetscMPIInt     size;
55   PetscErrorCode  ierr;
56 
57   PetscFunctionBegin;
58   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);CHKERRQ(ierr);
59   /* Create patch DM */
60   ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0, &dof, 0,0,0,0, &st);CHKERRQ(ierr);
61 
62   /* Get piece for rank r, expanded by halo */
63   bupper.i = PetscMin(M, upper.i + halo); blower.i = PetscMax(lower.i - halo, 0);
64   bupper.j = PetscMin(N, upper.j + halo); blower.j = PetscMax(lower.j - halo, 0);
65   bupper.k = PetscMin(P, upper.k + halo); blower.k = PetscMax(lower.k - halo, 0);
66   rM       = bupper.i - blower.i;
67   rN       = bupper.j - blower.j;
68   rP       = bupper.k - blower.k;
69 
70   if (commz != MPI_COMM_NULL) {
71     ierr = DMDACreate(commz, dmz);CHKERRQ(ierr);
72     ierr = DMSetDimension(*dmz, dim);CHKERRQ(ierr);
73     ierr = DMDASetSizes(*dmz, rM, rN, rP);CHKERRQ(ierr);
74     ierr = DMDASetNumProcs(*dmz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE);CHKERRQ(ierr);
75     ierr = DMDASetBoundaryType(*dmz, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE);CHKERRQ(ierr);
76     ierr = DMDASetDof(*dmz, dof);CHKERRQ(ierr);
77     ierr = DMDASetStencilType(*dmz, st);CHKERRQ(ierr);
78     ierr = DMDASetStencilWidth(*dmz, 0);CHKERRQ(ierr);
79     ierr = DMDASetOwnershipRanges(*dmz, NULL, NULL, NULL);CHKERRQ(ierr);
80     ierr = DMSetFromOptions(*dmz);CHKERRQ(ierr);
81     ierr = DMSetUp(*dmz);CHKERRQ(ierr);
82     ierr = DMDAGetCorners(*dmz, &sxb, &syb, &szb, &mxb, &myb, &mzb);CHKERRQ(ierr);
83     sxr  = PetscMax(sxb,     lower.i - blower.i);
84     syr  = PetscMax(syb,     lower.j - blower.j);
85     szr  = PetscMax(szb,     lower.k - blower.k);
86     exr  = PetscMin(sxb+mxb, upper.i - blower.i);
87     eyr  = PetscMin(syb+myb, upper.j - blower.j);
88     ezr  = PetscMin(szb+mzb, upper.k - blower.k);
89     ierr = PetscMalloc2(rM*rN*rP,&localPoints,rM*rN*rP,&remotePoints);CHKERRQ(ierr);
90   } else {
91     sxr = syr = szr = exr = eyr = ezr = sxb = syb = szb = mxb = myb = mzb = 0;
92   }
93 
94   /* Create SF for restricted map */
95   ierr = VecGetOwnershipRanges(X,&ranges);CHKERRQ(ierr);
96 
97   loclower.i = blower.i + sxr; locupper.i = blower.i + exr;
98   loclower.j = blower.j + syr; locupper.j = blower.j + eyr;
99   loclower.k = blower.k + szr; locupper.k = blower.k + ezr;
100 
101   ierr = DMDACreatePatchIS(dm, &loclower, &locupper, &is);CHKERRQ(ierr);
102   ierr = ISGetIndices(is, &indices);CHKERRQ(ierr);
103 
104   q = 0;
105   for (k = szb; k < szb+mzb; ++k) {
106     if ((k < szr) || (k >= ezr)) continue;
107     for (j = syb; j < syb+myb; ++j) {
108       if ((j < syr) || (j >= eyr)) continue;
109       for (i = sxb; i < sxb+mxb; ++i) {
110         const PetscInt lp = ((k-szb)*rN + (j-syb))*rM + i-sxb;
111         PetscInt       r;
112 
113         if ((i < sxr) || (i >= exr)) continue;
114         localPoints[q]        = lp;
115         ierr = PetscFindInt(indices[q], size+1, ranges, &r);CHKERRQ(ierr);
116 
117         remotePoints[q].rank  = r < 0 ? -(r+1) - 1 : r;
118         remotePoints[q].index = indices[q] - ranges[remotePoints[q].rank];
119         ++q;
120       }
121     }
122   }
123   ierr = ISRestoreIndices(is, &indices);CHKERRQ(ierr);
124   ierr = ISDestroy(&is);CHKERRQ(ierr);
125   ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfzr);CHKERRQ(ierr);
126   ierr = PetscObjectSetName((PetscObject) *sfzr, "Restricted Map");CHKERRQ(ierr);
127   ierr = PetscSFSetGraph(*sfzr, M*N*P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES);CHKERRQ(ierr);
128 
129   /* Create SF for buffered map */
130   loclower.i = blower.i + sxb; locupper.i = blower.i + sxb+mxb;
131   loclower.j = blower.j + syb; locupper.j = blower.j + syb+myb;
132   loclower.k = blower.k + szb; locupper.k = blower.k + szb+mzb;
133 
134   ierr = DMDACreatePatchIS(dm, &loclower, &locupper, &is);CHKERRQ(ierr);
135   ierr = ISGetIndices(is, &indices);CHKERRQ(ierr);
136 
137   q = 0;
138   for (k = szb; k < szb+mzb; ++k) {
139     for (j = syb; j < syb+myb; ++j) {
140       for (i = sxb; i < sxb+mxb; ++i, ++q) {
141         PetscInt r;
142 
143         localPoints[q]        = q;
144         ierr = PetscFindInt(indices[q], size+1, ranges, &r);CHKERRQ(ierr);
145         remotePoints[q].rank  = r < 0 ? -(r+1) - 1 : r;
146         remotePoints[q].index = indices[q] - ranges[remotePoints[q].rank];
147       }
148     }
149   }
150   ierr = ISRestoreIndices(is, &indices);CHKERRQ(ierr);
151   ierr = ISDestroy(&is);CHKERRQ(ierr);
152   ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfz);CHKERRQ(ierr);
153   ierr = PetscObjectSetName((PetscObject) *sfz, "Buffered Map");CHKERRQ(ierr);
154   ierr = PetscSFSetGraph(*sfz, M*N*P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES);CHKERRQ(ierr);
155 
156   ierr = PetscFree2(localPoints, remotePoints);CHKERRQ(ierr);
157   PetscFunctionReturn(0);
158 }
159 
160 typedef enum {PATCH_COMM_TYPE_WORLD = 0, PATCH_COMM_TYPE_SELF = 1} PatchCommType;
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "DMPatchSolve"
164 PetscErrorCode DMPatchSolve(DM dm)
165 {
166   MPI_Comm       comm;
167   MPI_Comm       commz;
168   DM             dmc;
169   PetscSF        sfz, sfzr;
170   Vec            XC;
171   MatStencil     patchSize, commSize, gridRank, lower, upper;
172   PetscInt       M, N, P, i, j, k, l, m, n, p = 0;
173   PetscMPIInt    rank, size;
174   PetscInt       debug = 0;
175   PetscErrorCode ierr;
176 
177   PetscFunctionBegin;
178   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
179   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
180   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
181   ierr = DMPatchGetCoarse(dm, &dmc);CHKERRQ(ierr);
182   ierr = DMPatchGetPatchSize(dm, &patchSize);CHKERRQ(ierr);
183   ierr = DMPatchGetCommSize(dm, &commSize);CHKERRQ(ierr);
184   ierr = DMPatchGetCommSize(dm, &commSize);CHKERRQ(ierr);
185   ierr = DMGetGlobalVector(dmc, &XC);CHKERRQ(ierr);
186   ierr = DMDAGetInfo(dmc, 0, &M, &N, &P, &l, &m, &n, 0,0,0,0,0,0);CHKERRQ(ierr);
187   M    = PetscMax(M, 1); l = PetscMax(l, 1);
188   N    = PetscMax(N, 1); m = PetscMax(m, 1);
189   P    = PetscMax(P, 1); n = PetscMax(n, 1);
190 
191   gridRank.i = rank       % l;
192   gridRank.j = rank/l     % m;
193   gridRank.k = rank/(l*m) % n;
194 
195   if (commSize.i*commSize.j*commSize.k == size || commSize.i*commSize.j*commSize.k == 0) {
196     commSize.i = l; commSize.j = m; commSize.k = n;
197     commz      = comm;
198   } else if (commSize.i*commSize.j*commSize.k == 1) {
199     commz = PETSC_COMM_SELF;
200   } else {
201     const PetscMPIInt newComm = ((gridRank.k/commSize.k)*(m/commSize.j) + gridRank.j/commSize.j)*(l/commSize.i) + (gridRank.i/commSize.i);
202     const PetscMPIInt newRank = ((gridRank.k%commSize.k)*commSize.j     + gridRank.j%commSize.j)*commSize.i     + (gridRank.i%commSize.i);
203 
204     ierr = MPI_Comm_split(comm, newComm, newRank, &commz);CHKERRQ(ierr);
205     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "Rank %d color %d key %d commz %d\n", rank, newComm, newRank, *((PetscMPIInt*) &commz));CHKERRQ(ierr);}
206   }
207   /*
208    Assumptions:
209      - patchSize divides gridSize
210      - commSize divides gridSize
211      - commSize divides l,m,n
212    Ignore multiple patches per rank for now
213 
214    Multiple ranks per patch:
215      - l,m,n divides patchSize
216      - commSize divides patchSize
217    */
218   for (k = 0; k < P; k += PetscMax(patchSize.k, 1)) {
219     for (j = 0; j < N; j += PetscMax(patchSize.j, 1)) {
220       for (i = 0; i < M; i += PetscMax(patchSize.i, 1), ++p) {
221         MPI_Comm commp = MPI_COMM_NULL;
222         DM       dmz   = NULL;
223 #if 0
224         DM  dmf     = NULL;
225         Mat interpz = NULL;
226 #endif
227         Vec         XZ       = NULL;
228         PetscScalar *xcarray = NULL;
229         PetscScalar *xzarray = NULL;
230 
231         if ((gridRank.k/commSize.k == p/(l/commSize.i * m/commSize.j) % n/commSize.k) &&
232             (gridRank.j/commSize.j == p/(l/commSize.i)                % m/commSize.j) &&
233             (gridRank.i/commSize.i == p                               % l/commSize.i)) {
234           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "Rank %d is accepting Patch %d\n", rank, p);CHKERRQ(ierr);}
235           commp = commz;
236         }
237         /* Zoom to coarse patch */
238         lower.i = i; lower.j = j; lower.k = k;
239         upper.i = i + patchSize.i; upper.j = j + patchSize.j; upper.k = k + patchSize.k;
240         ierr    = DMPatchZoom(dmc, XC, lower, upper, commp, &dmz, &sfz, &sfzr);CHKERRQ(ierr);
241         lower.c = 0; /* initialize member, otherwise compiler issues warnings */
242         upper.c = 0; /* initialize member, otherwise compiler issues warnings */
243         /* Debug */
244         ierr = PetscPrintf(comm, "Patch %d: (%d, %d, %d)--(%d, %d, %d)\n", p, lower.i, lower.j, lower.k, upper.i, upper.j, upper.k);CHKERRQ(ierr);
245         if (dmz) {ierr = DMView(dmz, PETSC_VIEWER_STDOUT_(commz));CHKERRQ(ierr);}
246         ierr = PetscSFView(sfz,  PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr);
247         ierr = PetscSFView(sfzr, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr);
248         /* Scatter Xcoarse -> Xzoom */
249         if (dmz) {ierr = DMGetGlobalVector(dmz, &XZ);CHKERRQ(ierr);}
250         if (XZ)  {ierr = VecGetArray(XZ, &xzarray);CHKERRQ(ierr);}
251         ierr = VecGetArray(XC, &xcarray);CHKERRQ(ierr);
252         ierr = PetscSFBcastBegin(sfz, MPIU_SCALAR, xcarray, xzarray);CHKERRQ(ierr);
253         ierr = PetscSFBcastEnd(sfz, MPIU_SCALAR, xcarray, xzarray);CHKERRQ(ierr);
254         ierr = VecRestoreArray(XC, &xcarray);CHKERRQ(ierr);
255         if (XZ)  {ierr = VecRestoreArray(XZ, &xzarray);CHKERRQ(ierr);}
256 #if 0
257         /* Interpolate Xzoom -> Xfine, note that this may be on subcomms */
258         ierr = DMRefine(dmz, MPI_COMM_NULL, &dmf);CHKERRQ(ierr);
259         ierr = DMCreateInterpolation(dmz, dmf, &interpz, NULL);CHKERRQ(ierr);
260         ierr = DMInterpolate(dmz, interpz, dmf);CHKERRQ(ierr);
261         /* Smooth Xfine using two-step smoother, normal smoother plus Kaczmarz---moves back and forth from dmzoom to dmfine */
262         /* Compute residual Rfine */
263         /* Restrict Rfine to Rzoom_restricted */
264 #endif
265         /* Scatter Rzoom_restricted -> Rcoarse_restricted */
266         if (XZ)  {ierr = VecGetArray(XZ, &xzarray);CHKERRQ(ierr);}
267         ierr = VecGetArray(XC, &xcarray);CHKERRQ(ierr);
268         ierr = PetscSFReduceBegin(sfzr, MPIU_SCALAR, xzarray, xcarray, MPI_SUM);CHKERRQ(ierr);
269         ierr = PetscSFReduceEnd(sfzr, MPIU_SCALAR, xzarray, xcarray, MPI_SUM);CHKERRQ(ierr);
270         ierr = VecRestoreArray(XC, &xcarray);CHKERRQ(ierr);
271         if (XZ)  {ierr = VecRestoreArray(XZ, &xzarray);CHKERRQ(ierr);}
272         if (dmz) {ierr = DMRestoreGlobalVector(dmz, &XZ);CHKERRQ(ierr);}
273         /* Compute global residual Rcoarse */
274         /* TauCoarse = Rcoarse - Rcoarse_restricted */
275 
276         ierr = PetscSFDestroy(&sfz);CHKERRQ(ierr);
277         ierr = PetscSFDestroy(&sfzr);CHKERRQ(ierr);
278         ierr = DMDestroy(&dmz);CHKERRQ(ierr);
279       }
280     }
281   }
282   ierr = DMRestoreGlobalVector(dmc, &XC);CHKERRQ(ierr);
283   PetscFunctionReturn(0);
284 }
285 
286 #undef __FUNCT__
287 #define __FUNCT__ "DMPatchView_Ascii"
288 PetscErrorCode DMPatchView_Ascii(DM dm, PetscViewer viewer)
289 {
290   DM_Patch          *mesh = (DM_Patch*) dm->data;
291   PetscViewerFormat format;
292   const char        *name;
293   PetscErrorCode    ierr;
294 
295   PetscFunctionBegin;
296   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
297   /* if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) */
298   ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr);
299   ierr = PetscViewerASCIIPrintf(viewer, "Patch DM %s\n", name);CHKERRQ(ierr);
300   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
301   ierr = PetscViewerASCIIPrintf(viewer, "Coarse DM\n");CHKERRQ(ierr);
302   ierr = DMView(mesh->dmCoarse, viewer);CHKERRQ(ierr);
303   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
304   PetscFunctionReturn(0);
305 }
306 
307 #undef __FUNCT__
308 #define __FUNCT__ "DMView_Patch"
309 PetscErrorCode DMView_Patch(DM dm, PetscViewer viewer)
310 {
311   PetscBool      iascii, isbinary;
312   PetscErrorCode ierr;
313 
314   PetscFunctionBegin;
315   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
316   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
317   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
318   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);CHKERRQ(ierr);
319   if (iascii) {
320     ierr = DMPatchView_Ascii(dm, viewer);CHKERRQ(ierr);
321 #if 0
322   } else if (isbinary) {
323     ierr = DMPatchView_Binary(dm, viewer);CHKERRQ(ierr);
324 #endif
325   }
326   PetscFunctionReturn(0);
327 }
328 
329 #undef __FUNCT__
330 #define __FUNCT__ "DMDestroy_Patch"
331 PetscErrorCode DMDestroy_Patch(DM dm)
332 {
333   DM_Patch       *mesh = (DM_Patch*) dm->data;
334   PetscErrorCode ierr;
335 
336   PetscFunctionBegin;
337   if (--mesh->refct > 0) PetscFunctionReturn(0);
338   ierr = DMDestroy(&mesh->dmCoarse);CHKERRQ(ierr);
339   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
340   ierr = PetscFree(mesh);CHKERRQ(ierr);
341   PetscFunctionReturn(0);
342 }
343 
344 #undef __FUNCT__
345 #define __FUNCT__ "DMSetUp_Patch"
346 PetscErrorCode DMSetUp_Patch(DM dm)
347 {
348   DM_Patch       *mesh = (DM_Patch*) dm->data;
349   PetscErrorCode ierr;
350 
351   PetscFunctionBegin;
352   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
353   ierr = DMSetUp(mesh->dmCoarse);CHKERRQ(ierr);
354   PetscFunctionReturn(0);
355 }
356 
357 #undef __FUNCT__
358 #define __FUNCT__ "DMCreateGlobalVector_Patch"
359 PetscErrorCode DMCreateGlobalVector_Patch(DM dm, Vec *g)
360 {
361   DM_Patch       *mesh = (DM_Patch*) dm->data;
362   PetscErrorCode ierr;
363 
364   PetscFunctionBegin;
365   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
366   ierr = DMCreateGlobalVector(mesh->dmCoarse, g);CHKERRQ(ierr);
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "DMCreateLocalVector_Patch"
372 PetscErrorCode DMCreateLocalVector_Patch(DM dm, Vec *l)
373 {
374   DM_Patch       *mesh = (DM_Patch*) dm->data;
375   PetscErrorCode ierr;
376 
377   PetscFunctionBegin;
378   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
379   ierr = DMCreateLocalVector(mesh->dmCoarse, l);CHKERRQ(ierr);
380   PetscFunctionReturn(0);
381 }
382 
383 #undef __FUNCT__
384 #define __FUNCT__ "DMCreateSubDM_Patch"
385 PetscErrorCode DMCreateSubDM_Patch(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
386 {
387   PetscFunctionBegin;
388   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
389   SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tell me to code this");
390   PetscFunctionReturn(0);
391 }
392 
393 #undef __FUNCT__
394 #define __FUNCT__ "DMPatchGetCoarse"
395 PetscErrorCode DMPatchGetCoarse(DM dm, DM *dmCoarse)
396 {
397   DM_Patch *mesh = (DM_Patch*) dm->data;
398 
399   PetscFunctionBegin;
400   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
401   *dmCoarse = mesh->dmCoarse;
402   PetscFunctionReturn(0);
403 }
404 
405 #undef __FUNCT__
406 #define __FUNCT__ "DMPatchGetPatchSize"
407 PetscErrorCode DMPatchGetPatchSize(DM dm, MatStencil *patchSize)
408 {
409   DM_Patch *mesh = (DM_Patch*) dm->data;
410 
411   PetscFunctionBegin;
412   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
413   PetscValidPointer(patchSize, 2);
414   *patchSize = mesh->patchSize;
415   PetscFunctionReturn(0);
416 }
417 
418 #undef __FUNCT__
419 #define __FUNCT__ "DMPatchSetPatchSize"
420 PetscErrorCode DMPatchSetPatchSize(DM dm, MatStencil patchSize)
421 {
422   DM_Patch *mesh = (DM_Patch*) dm->data;
423 
424   PetscFunctionBegin;
425   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
426   mesh->patchSize = patchSize;
427   PetscFunctionReturn(0);
428 }
429 
430 #undef __FUNCT__
431 #define __FUNCT__ "DMPatchGetCommSize"
432 PetscErrorCode DMPatchGetCommSize(DM dm, MatStencil *commSize)
433 {
434   DM_Patch *mesh = (DM_Patch*) dm->data;
435 
436   PetscFunctionBegin;
437   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
438   PetscValidPointer(commSize, 2);
439   *commSize = mesh->commSize;
440   PetscFunctionReturn(0);
441 }
442 
443 #undef __FUNCT__
444 #define __FUNCT__ "DMPatchSetCommSize"
445 PetscErrorCode DMPatchSetCommSize(DM dm, MatStencil commSize)
446 {
447   DM_Patch *mesh = (DM_Patch*) dm->data;
448 
449   PetscFunctionBegin;
450   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
451   mesh->commSize = commSize;
452   PetscFunctionReturn(0);
453 }
454