xref: /petsc/src/dm/impls/da/dagetelem.c (revision 517f4e3460cd8c7e73c19b8bf533f1efe47f30a7)
1 #include <petsc/private/dmdaimpl.h> /*I  "petscdmda.h"   I*/
2 
3 static PetscErrorCode DMDAGetElements_1D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[])
4 {
5   DM_DA   *da = (DM_DA *)dm->data;
6   PetscInt i, xs, xe, Xs, Xe;
7   PetscInt cnt = 0;
8 
9   PetscFunctionBegin;
10   if (!da->e) {
11     PetscInt corners[2];
12 
13     PetscCheck(da->s, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot get elements for DMDA with zero stencil width");
14     PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xe, NULL, NULL));
15     PetscCall(DMDAGetGhostCorners(dm, &Xs, NULL, NULL, &Xe, NULL, NULL));
16     xe += xs;
17     Xe += Xs;
18     if (xs != Xs) xs -= 1;
19     da->ne = 1 * (xe - xs - 1);
20     PetscCall(PetscMalloc1(1 + 2 * da->ne, &da->e));
21     for (i = xs; i < xe - 1; i++) {
22       da->e[cnt++] = (i - Xs);
23       da->e[cnt++] = (i - Xs + 1);
24     }
25     da->nen = 2;
26 
27     corners[0] = (xs - Xs);
28     corners[1] = (xe - 1 - Xs);
29     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2, corners, PETSC_COPY_VALUES, &da->ecorners));
30   }
31   *nel = da->ne;
32   *nen = da->nen;
33   *e   = da->e;
34   PetscFunctionReturn(PETSC_SUCCESS);
35 }
36 
37 static PetscErrorCode DMDAGetElements_2D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[])
38 {
39   DM_DA   *da = (DM_DA *)dm->data;
40   PetscInt i, xs, xe, Xs, Xe;
41   PetscInt j, ys, ye, Ys, Ye;
42   PetscInt cnt = 0, cell[4], ns = 2;
43   PetscInt c, split[] = {0, 1, 3, 2, 3, 1};
44 
45   PetscFunctionBegin;
46   if (!da->e) {
47     PetscInt corners[4], nn = 0;
48 
49     PetscCheck(da->s, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot get elements for DMDA with zero stencil width");
50 
51     switch (da->elementtype) {
52     case DMDA_ELEMENT_Q1:
53       da->nen = 4;
54       break;
55     case DMDA_ELEMENT_P1:
56       da->nen = 3;
57       break;
58     default:
59       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown element type %d", da->elementtype);
60     }
61     nn = da->nen;
62 
63     if (da->elementtype == DMDA_ELEMENT_P1) ns = 2;
64     if (da->elementtype == DMDA_ELEMENT_Q1) ns = 1;
65     PetscCall(DMDAGetCorners(dm, &xs, &ys, NULL, &xe, &ye, NULL));
66     PetscCall(DMDAGetGhostCorners(dm, &Xs, &Ys, NULL, &Xe, &Ye, NULL));
67     xe += xs;
68     Xe += Xs;
69     if (xs != Xs) xs -= 1;
70     ye += ys;
71     Ye += Ys;
72     if (ys != Ys) ys -= 1;
73     da->ne = ns * (xe - xs - 1) * (ye - ys - 1);
74     PetscCall(PetscMalloc1(1 + nn * da->ne, &da->e));
75     for (j = ys; j < ye - 1; j++) {
76       for (i = xs; i < xe - 1; i++) {
77         cell[0] = (i - Xs) + (j - Ys) * (Xe - Xs);
78         cell[1] = (i - Xs + 1) + (j - Ys) * (Xe - Xs);
79         cell[2] = (i - Xs + 1) + (j - Ys + 1) * (Xe - Xs);
80         cell[3] = (i - Xs) + (j - Ys + 1) * (Xe - Xs);
81         if (da->elementtype == DMDA_ELEMENT_P1) {
82           for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[split[c]];
83         }
84         if (da->elementtype == DMDA_ELEMENT_Q1) {
85           for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[c];
86         }
87       }
88     }
89 
90     corners[0] = (xs - Xs) + (ys - Ys) * (Xe - Xs);
91     corners[1] = (xe - 1 - Xs) + (ys - Ys) * (Xe - Xs);
92     corners[2] = (xs - Xs) + (ye - 1 - Ys) * (Xe - Xs);
93     corners[3] = (xe - 1 - Xs) + (ye - 1 - Ys) * (Xe - Xs);
94     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 4, corners, PETSC_COPY_VALUES, &da->ecorners));
95   }
96   *nel = da->ne;
97   *nen = da->nen;
98   *e   = da->e;
99   PetscFunctionReturn(PETSC_SUCCESS);
100 }
101 
102 static PetscErrorCode DMDAGetElements_3D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[])
103 {
104   DM_DA   *da = (DM_DA *)dm->data;
105   PetscInt i, xs, xe, Xs, Xe;
106   PetscInt j, ys, ye, Ys, Ye;
107   PetscInt k, zs, ze, Zs, Ze;
108   PetscInt cnt = 0, cell[8], ns = 6;
109   PetscInt c, split[] = {0, 1, 3, 7, 0, 1, 7, 4, 1, 2, 3, 7, 1, 2, 7, 6, 1, 4, 5, 7, 1, 5, 6, 7};
110 
111   PetscFunctionBegin;
112   if (!da->e) {
113     PetscInt corners[8], nn = 0;
114 
115     PetscCheck(da->s, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot get elements for DMDA with zero stencil width");
116 
117     switch (da->elementtype) {
118     case DMDA_ELEMENT_Q1:
119       da->nen = 8;
120       break;
121     case DMDA_ELEMENT_P1:
122       da->nen = 4;
123       break;
124     default:
125       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown element type %d", da->elementtype);
126     }
127     nn = da->nen;
128 
129     if (da->elementtype == DMDA_ELEMENT_P1) ns = 6;
130     if (da->elementtype == DMDA_ELEMENT_Q1) ns = 1;
131     PetscCall(DMDAGetCorners(dm, &xs, &ys, &zs, &xe, &ye, &ze));
132     PetscCall(DMDAGetGhostCorners(dm, &Xs, &Ys, &Zs, &Xe, &Ye, &Ze));
133     xe += xs;
134     Xe += Xs;
135     if (xs != Xs) xs -= 1;
136     ye += ys;
137     Ye += Ys;
138     if (ys != Ys) ys -= 1;
139     ze += zs;
140     Ze += Zs;
141     if (zs != Zs) zs -= 1;
142     da->ne = ns * (xe - xs - 1) * (ye - ys - 1) * (ze - zs - 1);
143     PetscCall(PetscMalloc1(1 + nn * da->ne, &da->e));
144     for (k = zs; k < ze - 1; k++) {
145       for (j = ys; j < ye - 1; j++) {
146         for (i = xs; i < xe - 1; i++) {
147           cell[0] = (i - Xs) + (j - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys);
148           cell[1] = (i + 1 - Xs) + (j - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys);
149           cell[2] = (i + 1 - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys);
150           cell[3] = (i - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys);
151           cell[4] = (i - Xs) + (j - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys);
152           cell[5] = (i + 1 - Xs) + (j - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys);
153           cell[6] = (i + 1 - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys);
154           cell[7] = (i - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys);
155           if (da->elementtype == DMDA_ELEMENT_P1) {
156             for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[split[c]];
157           }
158           if (da->elementtype == DMDA_ELEMENT_Q1) {
159             for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[c];
160           }
161         }
162       }
163     }
164 
165     corners[0] = (xs - Xs) + (ys - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys);
166     corners[1] = (xe - 1 - Xs) + (ys - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys);
167     corners[2] = (xs - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys);
168     corners[3] = (xe - 1 - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys);
169     corners[4] = (xs - Xs) + (ys - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys);
170     corners[5] = (xe - 1 - Xs) + (ys - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys);
171     corners[6] = (xs - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys);
172     corners[7] = (xe - 1 - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys);
173     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 8, corners, PETSC_COPY_VALUES, &da->ecorners));
174   }
175   *nel = da->ne;
176   *nen = da->nen;
177   *e   = da->e;
178   PetscFunctionReturn(PETSC_SUCCESS);
179 }
180 
181 /*@
182   DMDAGetElementsCorners - Returns the global (x,y,z) indices of the lower left
183   corner of the non-overlapping decomposition identified by `DMDAGetElements()`
184 
185   Not Collective
186 
187   Input Parameter:
188 . da - the `DMDA` object
189 
190   Output Parameters:
191 + gx - the x index
192 . gy - the y index
193 - gz - the z index
194 
195   Level: intermediate
196 
197 .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`
198 @*/
199 PetscErrorCode DMDAGetElementsCorners(DM da, PetscInt *gx, PetscInt *gy, PetscInt *gz)
200 {
201   PetscInt  xs, Xs;
202   PetscInt  ys, Ys;
203   PetscInt  zs, Zs;
204   PetscBool isda;
205 
206   PetscFunctionBegin;
207   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
208   if (gx) PetscAssertPointer(gx, 2);
209   if (gy) PetscAssertPointer(gy, 3);
210   if (gz) PetscAssertPointer(gz, 4);
211   PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda));
212   PetscCheck(isda, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)da)->type_name);
213   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, NULL, NULL, NULL));
214   PetscCall(DMDAGetGhostCorners(da, &Xs, &Ys, &Zs, NULL, NULL, NULL));
215   if (xs != Xs) xs -= 1;
216   if (ys != Ys) ys -= 1;
217   if (zs != Zs) zs -= 1;
218   if (gx) *gx = xs;
219   if (gy) *gy = ys;
220   if (gz) *gz = zs;
221   PetscFunctionReturn(PETSC_SUCCESS);
222 }
223 
224 /*@
225   DMDAGetElementsSizes - Gets the local number of elements per direction for the non-overlapping decomposition identified by `DMDAGetElements()`
226 
227   Not Collective
228 
229   Input Parameter:
230 . da - the `DMDA` object
231 
232   Output Parameters:
233 + mx - number of local elements in x-direction
234 . my - number of local elements in y-direction
235 - mz - number of local elements in z-direction
236 
237   Level: intermediate
238 
239   Note:
240   It returns the same number of elements, irrespective of the `DMDAElementType`
241 
242 .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements`
243 @*/
244 PetscErrorCode DMDAGetElementsSizes(DM da, PetscInt *mx, PetscInt *my, PetscInt *mz)
245 {
246   PetscInt  xs, xe, Xs;
247   PetscInt  ys, ye, Ys;
248   PetscInt  zs, ze, Zs;
249   PetscInt  dim;
250   PetscBool isda;
251 
252   PetscFunctionBegin;
253   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
254   if (mx) PetscAssertPointer(mx, 2);
255   if (my) PetscAssertPointer(my, 3);
256   if (mz) PetscAssertPointer(mz, 4);
257   PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda));
258   PetscCheck(isda, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)da)->type_name);
259   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xe, &ye, &ze));
260   PetscCall(DMDAGetGhostCorners(da, &Xs, &Ys, &Zs, NULL, NULL, NULL));
261   xe += xs;
262   if (xs != Xs) xs -= 1;
263   ye += ys;
264   if (ys != Ys) ys -= 1;
265   ze += zs;
266   if (zs != Zs) zs -= 1;
267   if (mx) *mx = 0;
268   if (my) *my = 0;
269   if (mz) *mz = 0;
270   PetscCall(DMGetDimension(da, &dim));
271   switch (dim) {
272   case 3:
273     if (mz) *mz = ze - zs - 1; /* fall through */
274   case 2:
275     if (my) *my = ye - ys - 1; /* fall through */
276   case 1:
277     if (mx) *mx = xe - xs - 1;
278     break;
279   }
280   PetscFunctionReturn(PETSC_SUCCESS);
281 }
282 
283 /*@
284   DMDASetElementType - Sets the element type to be returned by `DMDAGetElements()`
285 
286   Not Collective
287 
288   Input Parameter:
289 . da - the `DMDA` object
290 
291   Output Parameter:
292 . etype - the element type, currently either `DMDA_ELEMENT_P1` or `DMDA_ELEMENT_Q1`
293 
294   Level: intermediate
295 
296 .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDAGetElementType()`, `DMDAGetElements()`, `DMDARestoreElements()`
297 @*/
298 PetscErrorCode DMDASetElementType(DM da, DMDAElementType etype)
299 {
300   DM_DA    *dd = (DM_DA *)da->data;
301   PetscBool isda;
302 
303   PetscFunctionBegin;
304   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
305   PetscValidLogicalCollectiveEnum(da, etype, 2);
306   PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda));
307   if (!isda) PetscFunctionReturn(PETSC_SUCCESS);
308   if (dd->elementtype != etype) {
309     PetscCall(PetscFree(dd->e));
310     PetscCall(ISDestroy(&dd->ecorners));
311 
312     dd->elementtype = etype;
313     dd->ne          = 0;
314     dd->nen         = 0;
315     dd->e           = NULL;
316   }
317   PetscFunctionReturn(PETSC_SUCCESS);
318 }
319 
320 /*@
321   DMDAGetElementType - Gets the element type to be returned by `DMDAGetElements()`
322 
323   Not Collective
324 
325   Input Parameter:
326 . da - the `DMDA` object
327 
328   Output Parameter:
329 . etype - the element type, currently either `DMDA_ELEMENT_P1` or `DMDA_ELEMENT_Q1`
330 
331   Level: intermediate
332 
333 .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`, `DMDARestoreElements()`
334 @*/
335 PetscErrorCode DMDAGetElementType(DM da, DMDAElementType *etype)
336 {
337   DM_DA    *dd = (DM_DA *)da->data;
338   PetscBool isda;
339 
340   PetscFunctionBegin;
341   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
342   PetscAssertPointer(etype, 2);
343   PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda));
344   PetscCheck(isda, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)da)->type_name);
345   *etype = dd->elementtype;
346   PetscFunctionReturn(PETSC_SUCCESS);
347 }
348 
349 /*@C
350   DMDAGetElements - Gets an array containing the indices (in local coordinates)
351   of all the local elements
352 
353   Not Collective
354 
355   Input Parameter:
356 . dm - the `DMDA` object
357 
358   Output Parameters:
359 + nel - number of local elements
360 . nen - number of nodes in each element (for example in one dimension it is 2, in two dimensions it is 3 or 4)
361 - e   - the local indices of the elements' vertices
362 
363   Level: intermediate
364 
365   Notes:
366   Call `DMDARestoreElements()` once you have finished accessing the elements.
367 
368   Each process uniquely owns a subset of the elements. That is no element is owned by two or more processes.
369 
370   If on each process you integrate over its owned elements and use `ADD_VALUES` in `Vec`/`MatSetValuesLocal()` then you'll obtain the correct result.
371 
372   Fortran Notes:
373   Use
374 .vb
375    PetscScalar, pointer :: e(:)
376 .ve
377   to declare the element array
378 
379 .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `VecSetValuesLocal()`, `MatSetValuesLocal()`, `DMGlobalToLocalBegin()`, `DMLocalToGlobalBegin()`, `DMDARestoreElements()`
380 @*/
381 PetscErrorCode DMDAGetElements(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[])
382 {
383   PetscInt  dim;
384   DM_DA    *dd = (DM_DA *)dm->data;
385   PetscBool isda;
386 
387   PetscFunctionBegin;
388   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
389   PetscAssertPointer(nel, 2);
390   PetscAssertPointer(nen, 3);
391   PetscAssertPointer(e, 4);
392   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda));
393   PetscCheck(isda, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)dm)->type_name);
394   PetscCheck(dd->stencil_type != DMDA_STENCIL_STAR, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMDAGetElements() requires you use a stencil type of DMDA_STENCIL_BOX");
395   PetscCall(DMGetDimension(dm, &dim));
396   if (dd->e) {
397     *nel = dd->ne;
398     *nen = dd->nen;
399     *e   = dd->e;
400     PetscFunctionReturn(PETSC_SUCCESS);
401   }
402   if (dim == -1) {
403     *nel = 0;
404     *nen = 0;
405     *e   = NULL;
406   } else if (dim == 1) {
407     PetscCall(DMDAGetElements_1D(dm, nel, nen, e));
408   } else if (dim == 2) {
409     PetscCall(DMDAGetElements_2D(dm, nel, nen, e));
410   } else if (dim == 3) {
411     PetscCall(DMDAGetElements_3D(dm, nel, nen, e));
412   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
413   PetscFunctionReturn(PETSC_SUCCESS);
414 }
415 
416 /*@
417   DMDAGetSubdomainCornersIS - Gets an index set containing the corner indices (in local coordinates)
418   of the non-overlapping decomposition identified by `DMDAGetElements()`
419 
420   Not Collective
421 
422   Input Parameter:
423 . dm - the `DMDA` object
424 
425   Output Parameter:
426 . is - the index set
427 
428   Level: intermediate
429 
430   Note:
431   Call `DMDARestoreSubdomainCornersIS()` once you have finished accessing the index set.
432 
433 .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`, `DMDARestoreElementsCornersIS()`
434 @*/
435 PetscErrorCode DMDAGetSubdomainCornersIS(DM dm, IS *is)
436 {
437   DM_DA    *dd = (DM_DA *)dm->data;
438   PetscBool isda;
439 
440   PetscFunctionBegin;
441   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
442   PetscAssertPointer(is, 2);
443   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda));
444   PetscCheck(isda, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)dm)->type_name);
445   PetscCheck(dd->stencil_type != DMDA_STENCIL_STAR, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMDAGetElement() requires you use a stencil type of DMDA_STENCIL_BOX");
446   if (!dd->ecorners) { /* compute elements if not yet done */
447     const PetscInt *e;
448     PetscInt        nel, nen;
449 
450     PetscCall(DMDAGetElements(dm, &nel, &nen, &e));
451     PetscCall(DMDARestoreElements(dm, &nel, &nen, &e));
452   }
453   *is = dd->ecorners;
454   PetscFunctionReturn(PETSC_SUCCESS);
455 }
456 
457 /*@C
458   DMDARestoreElements - Restores the array obtained with `DMDAGetElements()`
459 
460   Not Collective
461 
462   Input Parameters:
463 + dm  - the `DM` object
464 . nel - number of local elements
465 . nen - number of nodes in each element
466 - e   - the local indices of the elements' vertices
467 
468   Level: intermediate
469 
470   Note:
471   This restore signals the `DMDA` object that you no longer need access to the array information.
472 
473 .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`
474 @*/
475 PetscErrorCode DMDARestoreElements(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[])
476 {
477   PetscFunctionBegin;
478   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
479   PetscAssertPointer(nel, 2);
480   PetscAssertPointer(nen, 3);
481   PetscAssertPointer(e, 4);
482   if (nel) *nel = 0;
483   if (nen) *nen = -1;
484   if (e) *e = NULL;
485   PetscFunctionReturn(PETSC_SUCCESS);
486 }
487 
488 /*@
489   DMDARestoreSubdomainCornersIS - Restores the `IS` obtained with `DMDAGetSubdomainCornersIS()`
490 
491   Not Collective
492 
493   Input Parameters:
494 + dm - the `DM` object
495 - is - the index set
496 
497   Level: intermediate
498 
499 .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetSubdomainCornersIS()`
500 @*/
501 PetscErrorCode DMDARestoreSubdomainCornersIS(DM dm, IS *is)
502 {
503   PetscFunctionBegin;
504   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
505   PetscValidHeaderSpecific(*is, IS_CLASSID, 2);
506   *is = NULL;
507   PetscFunctionReturn(PETSC_SUCCESS);
508 }
509