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