xref: /petsc/src/dm/impls/plex/plexgenerate.c (revision e6d0a238963c2a97dd04845ea512b529992c7cdb)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 
3 PetscErrorCode DMPlexInvertCell_Internal(PetscInt dim, PetscInt numCorners, PetscInt cone[])
4 {
5   int tmpc;
6 
7   PetscFunctionBegin;
8   if (dim != 3) PetscFunctionReturn(0);
9   switch (numCorners) {
10   case 4:
11     tmpc    = cone[0];
12     cone[0] = cone[1];
13     cone[1] = tmpc;
14     break;
15   case 8:
16     tmpc    = cone[1];
17     cone[1] = cone[3];
18     cone[3] = tmpc;
19     break;
20   default: break;
21   }
22   PetscFunctionReturn(0);
23 }
24 
25 /*@C
26   DMPlexInvertCell - This flips tetrahedron and hexahedron orientation since Plex stores them internally with outward normals. Other cells are left untouched.
27 
28   Input Parameters:
29 + numCorners - The number of vertices in a cell
30 - cone - The incoming cone
31 
32   Output Parameter:
33 . cone - The inverted cone (in-place)
34 
35   Level: developer
36 
37 .seealso: DMPlexGenerate()
38 @*/
39 PetscErrorCode DMPlexInvertCell(PetscInt dim, PetscInt numCorners, int cone[])
40 {
41   int tmpc;
42 
43   PetscFunctionBegin;
44   if (dim != 3) PetscFunctionReturn(0);
45   switch (numCorners) {
46   case 4:
47     tmpc    = cone[0];
48     cone[0] = cone[1];
49     cone[1] = tmpc;
50     break;
51   case 8:
52     tmpc    = cone[1];
53     cone[1] = cone[3];
54     cone[3] = tmpc;
55     break;
56   default: break;
57   }
58   PetscFunctionReturn(0);
59 }
60 
61 /* This is to fix the tetrahedron orientation from TetGen */
62 PETSC_UNUSED static PetscErrorCode DMPlexInvertCells_Internal(PetscInt dim, PetscInt numCells, PetscInt numCorners, int cells[])
63 {
64   PetscInt       bound = numCells*numCorners, coff;
65   PetscErrorCode ierr;
66 
67   PetscFunctionBegin;
68   for (coff = 0; coff < bound; coff += numCorners) {
69     ierr = DMPlexInvertCell(dim, numCorners, &cells[coff]);CHKERRQ(ierr);
70   }
71   PetscFunctionReturn(0);
72 }
73 
74 /*@C
75   DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator
76 
77   Not Collective
78 
79   Inputs Parameters:
80 + dm - The DMPlex object
81 - opts - The command line options
82 
83   Level: developer
84 
85 .keywords: mesh, points
86 .seealso: DMPlexTetgenSetOptions(), DMPlexGenerate()
87 @*/
88 PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts)
89 {
90   DM_Plex       *mesh = (DM_Plex*) dm->data;
91   PetscErrorCode ierr;
92 
93   PetscFunctionBegin;
94   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
95   PetscValidPointer(opts, 2);
96   ierr = PetscFree(mesh->triangleOpts);CHKERRQ(ierr);
97   ierr = PetscStrallocpy(opts, &mesh->triangleOpts);CHKERRQ(ierr);
98   PetscFunctionReturn(0);
99 }
100 
101 /*@C
102   DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator
103 
104   Not Collective
105 
106   Inputs Parameters:
107 + dm - The DMPlex object
108 - opts - The command line options
109 
110   Level: developer
111 
112 .keywords: mesh, points
113 .seealso: DMPlexTriangleSetOptions(), DMPlexGenerate()
114 @*/
115 PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts)
116 {
117   DM_Plex       *mesh = (DM_Plex*) dm->data;
118   PetscErrorCode ierr;
119 
120   PetscFunctionBegin;
121   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
122   PetscValidPointer(opts, 2);
123   ierr = PetscFree(mesh->tetgenOpts);CHKERRQ(ierr);
124   ierr = PetscStrallocpy(opts, &mesh->tetgenOpts);CHKERRQ(ierr);
125   PetscFunctionReturn(0);
126 }
127 
128 #if defined(PETSC_HAVE_TRIANGLE)
129 #include <triangle.h>
130 
131 static PetscErrorCode InitInput_Triangle(struct triangulateio *inputCtx)
132 {
133   PetscFunctionBegin;
134   inputCtx->numberofpoints             = 0;
135   inputCtx->numberofpointattributes    = 0;
136   inputCtx->pointlist                  = NULL;
137   inputCtx->pointattributelist         = NULL;
138   inputCtx->pointmarkerlist            = NULL;
139   inputCtx->numberofsegments           = 0;
140   inputCtx->segmentlist                = NULL;
141   inputCtx->segmentmarkerlist          = NULL;
142   inputCtx->numberoftriangleattributes = 0;
143   inputCtx->trianglelist               = NULL;
144   inputCtx->numberofholes              = 0;
145   inputCtx->holelist                   = NULL;
146   inputCtx->numberofregions            = 0;
147   inputCtx->regionlist                 = NULL;
148   PetscFunctionReturn(0);
149 }
150 
151 static PetscErrorCode InitOutput_Triangle(struct triangulateio *outputCtx)
152 {
153   PetscFunctionBegin;
154   outputCtx->numberofpoints        = 0;
155   outputCtx->pointlist             = NULL;
156   outputCtx->pointattributelist    = NULL;
157   outputCtx->pointmarkerlist       = NULL;
158   outputCtx->numberoftriangles     = 0;
159   outputCtx->trianglelist          = NULL;
160   outputCtx->triangleattributelist = NULL;
161   outputCtx->neighborlist          = NULL;
162   outputCtx->segmentlist           = NULL;
163   outputCtx->segmentmarkerlist     = NULL;
164   outputCtx->numberofedges         = 0;
165   outputCtx->edgelist              = NULL;
166   outputCtx->edgemarkerlist        = NULL;
167   PetscFunctionReturn(0);
168 }
169 
170 static PetscErrorCode FiniOutput_Triangle(struct triangulateio *outputCtx)
171 {
172   PetscFunctionBegin;
173   free(outputCtx->pointlist);
174   free(outputCtx->pointmarkerlist);
175   free(outputCtx->segmentlist);
176   free(outputCtx->segmentmarkerlist);
177   free(outputCtx->edgelist);
178   free(outputCtx->edgemarkerlist);
179   free(outputCtx->trianglelist);
180   free(outputCtx->neighborlist);
181   PetscFunctionReturn(0);
182 }
183 
184 static PetscErrorCode DMPlexGenerate_Triangle(DM boundary, PetscBool interpolate, DM *dm)
185 {
186   MPI_Comm             comm;
187   DM_Plex             *mesh             = (DM_Plex *) boundary->data;
188   PetscInt             dim              = 2;
189   const PetscBool      createConvexHull = PETSC_FALSE;
190   const PetscBool      constrained      = PETSC_FALSE;
191   const char          *labelName        = "marker";
192   struct triangulateio in;
193   struct triangulateio out;
194   DMLabel              label;
195   PetscInt             vStart, vEnd, v, eStart, eEnd, e;
196   PetscMPIInt          rank;
197   PetscErrorCode       ierr;
198 
199   PetscFunctionBegin;
200   ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr);
201   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
202   ierr = InitInput_Triangle(&in);CHKERRQ(ierr);
203   ierr = InitOutput_Triangle(&out);CHKERRQ(ierr);
204   ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr);
205   ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr);
206 
207   in.numberofpoints = vEnd - vStart;
208   if (in.numberofpoints > 0) {
209     PetscSection coordSection;
210     Vec          coordinates;
211     PetscScalar *array;
212 
213     ierr = PetscMalloc1(in.numberofpoints*dim, &in.pointlist);CHKERRQ(ierr);
214     ierr = PetscMalloc1(in.numberofpoints, &in.pointmarkerlist);CHKERRQ(ierr);
215     ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr);
216     ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr);
217     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
218     for (v = vStart; v < vEnd; ++v) {
219       const PetscInt idx = v - vStart;
220       PetscInt       off, d;
221 
222       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
223       for (d = 0; d < dim; ++d) {
224         in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
225       }
226       if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);}
227     }
228     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
229   }
230   ierr  = DMPlexGetHeightStratum(boundary, 0, &eStart, &eEnd);CHKERRQ(ierr);
231   in.numberofsegments = eEnd - eStart;
232   if (in.numberofsegments > 0) {
233     ierr = PetscMalloc1(in.numberofsegments*2, &in.segmentlist);CHKERRQ(ierr);
234     ierr = PetscMalloc1(in.numberofsegments, &in.segmentmarkerlist);CHKERRQ(ierr);
235     for (e = eStart; e < eEnd; ++e) {
236       const PetscInt  idx = e - eStart;
237       const PetscInt *cone;
238 
239       ierr = DMPlexGetCone(boundary, e, &cone);CHKERRQ(ierr);
240 
241       in.segmentlist[idx*2+0] = cone[0] - vStart;
242       in.segmentlist[idx*2+1] = cone[1] - vStart;
243 
244       if (label) {ierr = DMLabelGetValue(label, e, &in.segmentmarkerlist[idx]);CHKERRQ(ierr);}
245     }
246   }
247 #if 0 /* Do not currently support holes */
248   PetscReal *holeCoords;
249   PetscInt   h, d;
250 
251   ierr = DMPlexGetHoles(boundary, &in.numberofholes, &holeCords);CHKERRQ(ierr);
252   if (in.numberofholes > 0) {
253     ierr = PetscMalloc1(in.numberofholes*dim, &in.holelist);CHKERRQ(ierr);
254     for (h = 0; h < in.numberofholes; ++h) {
255       for (d = 0; d < dim; ++d) {
256         in.holelist[h*dim+d] = holeCoords[h*dim+d];
257       }
258     }
259   }
260 #endif
261   if (!rank) {
262     char args[32];
263 
264     /* Take away 'Q' for verbose output */
265     ierr = PetscStrcpy(args, "pqezQ");CHKERRQ(ierr);
266     if (createConvexHull)   {ierr = PetscStrcat(args, "c");CHKERRQ(ierr);}
267     if (constrained)        {ierr = PetscStrcpy(args, "zepDQ");CHKERRQ(ierr);}
268     if (mesh->triangleOpts) {triangulate(mesh->triangleOpts, &in, &out, NULL);}
269     else                    {triangulate(args, &in, &out, NULL);}
270   }
271   ierr = PetscFree(in.pointlist);CHKERRQ(ierr);
272   ierr = PetscFree(in.pointmarkerlist);CHKERRQ(ierr);
273   ierr = PetscFree(in.segmentlist);CHKERRQ(ierr);
274   ierr = PetscFree(in.segmentmarkerlist);CHKERRQ(ierr);
275   ierr = PetscFree(in.holelist);CHKERRQ(ierr);
276 
277   {
278     DMLabel        glabel      = NULL;
279     const PetscInt numCorners  = 3;
280     const PetscInt numCells    = out.numberoftriangles;
281     const PetscInt numVertices = out.numberofpoints;
282     const int     *cells      = out.trianglelist;
283     const double  *meshCoords = out.pointlist;
284 
285     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr);
286     if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);}
287     /* Set labels */
288     for (v = 0; v < numVertices; ++v) {
289       if (out.pointmarkerlist[v]) {
290         if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);}
291       }
292     }
293     if (interpolate) {
294       for (e = 0; e < out.numberofedges; e++) {
295         if (out.edgemarkerlist[e]) {
296           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
297           const PetscInt *edges;
298           PetscInt        numEdges;
299 
300           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
301           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
302           if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
303           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
304         }
305       }
306     }
307     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
308   }
309 #if 0 /* Do not currently support holes */
310   ierr = DMPlexCopyHoles(*dm, boundary);CHKERRQ(ierr);
311 #endif
312   ierr = FiniOutput_Triangle(&out);CHKERRQ(ierr);
313   PetscFunctionReturn(0);
314 }
315 
316 static PetscErrorCode DMPlexRefine_Triangle(DM dm, double *maxVolumes, DM *dmRefined)
317 {
318   MPI_Comm             comm;
319   PetscInt             dim       = 2;
320   const char          *labelName = "marker";
321   struct triangulateio in;
322   struct triangulateio out;
323   DMLabel              label;
324   PetscInt             vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal;
325   PetscMPIInt          rank;
326   PetscErrorCode       ierr;
327 
328   PetscFunctionBegin;
329   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
330   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
331   ierr = InitInput_Triangle(&in);CHKERRQ(ierr);
332   ierr = InitOutput_Triangle(&out);CHKERRQ(ierr);
333   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
334   ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
335   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
336   ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr);
337 
338   in.numberofpoints = vEnd - vStart;
339   if (in.numberofpoints > 0) {
340     PetscSection coordSection;
341     Vec          coordinates;
342     PetscScalar *array;
343 
344     ierr = PetscMalloc1(in.numberofpoints*dim, &in.pointlist);CHKERRQ(ierr);
345     ierr = PetscMalloc1(in.numberofpoints, &in.pointmarkerlist);CHKERRQ(ierr);
346     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
347     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
348     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
349     for (v = vStart; v < vEnd; ++v) {
350       const PetscInt idx = v - vStart;
351       PetscInt       off, d;
352 
353       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
354       for (d = 0; d < dim; ++d) {
355         in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
356       }
357       if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);}
358     }
359     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
360   }
361   ierr  = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
362 
363   in.numberofcorners   = 3;
364   in.numberoftriangles = cEnd - cStart;
365 
366   in.trianglearealist  = (double*) maxVolumes;
367   if (in.numberoftriangles > 0) {
368     ierr = PetscMalloc1(in.numberoftriangles*in.numberofcorners, &in.trianglelist);CHKERRQ(ierr);
369     for (c = cStart; c < cEnd; ++c) {
370       const PetscInt idx      = c - cStart;
371       PetscInt      *closure = NULL;
372       PetscInt       closureSize;
373 
374       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
375       if ((closureSize != 4) && (closureSize != 7)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a triangle, %D vertices in closure", closureSize);
376       for (v = 0; v < 3; ++v) {
377         in.trianglelist[idx*in.numberofcorners + v] = closure[(v+closureSize-3)*2] - vStart;
378       }
379       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
380     }
381   }
382   /* TODO: Segment markers are missing on input */
383 #if 0 /* Do not currently support holes */
384   PetscReal *holeCoords;
385   PetscInt   h, d;
386 
387   ierr = DMPlexGetHoles(boundary, &in.numberofholes, &holeCords);CHKERRQ(ierr);
388   if (in.numberofholes > 0) {
389     ierr = PetscMalloc1(in.numberofholes*dim, &in.holelist);CHKERRQ(ierr);
390     for (h = 0; h < in.numberofholes; ++h) {
391       for (d = 0; d < dim; ++d) {
392         in.holelist[h*dim+d] = holeCoords[h*dim+d];
393       }
394     }
395   }
396 #endif
397   if (!rank) {
398     char args[32];
399 
400     /* Take away 'Q' for verbose output */
401     ierr = PetscStrcpy(args, "pqezQra");CHKERRQ(ierr);
402     triangulate(args, &in, &out, NULL);
403   }
404   ierr = PetscFree(in.pointlist);CHKERRQ(ierr);
405   ierr = PetscFree(in.pointmarkerlist);CHKERRQ(ierr);
406   ierr = PetscFree(in.segmentlist);CHKERRQ(ierr);
407   ierr = PetscFree(in.segmentmarkerlist);CHKERRQ(ierr);
408   ierr = PetscFree(in.trianglelist);CHKERRQ(ierr);
409 
410   {
411     DMLabel        rlabel      = NULL;
412     const PetscInt numCorners  = 3;
413     const PetscInt numCells    = out.numberoftriangles;
414     const PetscInt numVertices = out.numberofpoints;
415     const int     *cells      = out.trianglelist;
416     const double  *meshCoords = out.pointlist;
417     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
418 
419     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
420     if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);}
421     /* Set labels */
422     for (v = 0; v < numVertices; ++v) {
423       if (out.pointmarkerlist[v]) {
424         if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);}
425       }
426     }
427     if (interpolate) {
428       PetscInt e;
429 
430       for (e = 0; e < out.numberofedges; e++) {
431         if (out.edgemarkerlist[e]) {
432           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
433           const PetscInt *edges;
434           PetscInt        numEdges;
435 
436           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
437           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
438           if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
439           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
440         }
441       }
442     }
443     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
444   }
445 #if 0 /* Do not currently support holes */
446   ierr = DMPlexCopyHoles(*dm, boundary);CHKERRQ(ierr);
447 #endif
448   ierr = FiniOutput_Triangle(&out);CHKERRQ(ierr);
449   PetscFunctionReturn(0);
450 }
451 #endif /* PETSC_HAVE_TRIANGLE */
452 
453 #if defined(PETSC_HAVE_TETGEN)
454 #include <tetgen.h>
455 static PetscErrorCode DMPlexGenerate_Tetgen(DM boundary, PetscBool interpolate, DM *dm)
456 {
457   MPI_Comm       comm;
458   DM_Plex       *mesh      = (DM_Plex *) boundary->data;
459   const PetscInt dim       = 3;
460   const char    *labelName = "marker";
461   ::tetgenio     in;
462   ::tetgenio     out;
463   DMLabel        label;
464   PetscInt       vStart, vEnd, v, fStart, fEnd, f;
465   PetscMPIInt    rank;
466   PetscErrorCode ierr;
467 
468   PetscFunctionBegin;
469   ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr);
470   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
471   ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr);
472   ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr);
473 
474   in.numberofpoints = vEnd - vStart;
475   if (in.numberofpoints > 0) {
476     PetscSection coordSection;
477     Vec          coordinates;
478     PetscScalar *array;
479 
480     in.pointlist       = new double[in.numberofpoints*dim];
481     in.pointmarkerlist = new int[in.numberofpoints];
482 
483     ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr);
484     ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr);
485     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
486     for (v = vStart; v < vEnd; ++v) {
487       const PetscInt idx = v - vStart;
488       PetscInt       off, d;
489 
490       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
491       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
492       if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);}
493     }
494     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
495   }
496   ierr  = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr);
497 
498   in.numberoffacets = fEnd - fStart;
499   if (in.numberoffacets > 0) {
500     in.facetlist       = new tetgenio::facet[in.numberoffacets];
501     in.facetmarkerlist = new int[in.numberoffacets];
502     for (f = fStart; f < fEnd; ++f) {
503       const PetscInt idx     = f - fStart;
504       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v;
505 
506       in.facetlist[idx].numberofpolygons = 1;
507       in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
508       in.facetlist[idx].numberofholes    = 0;
509       in.facetlist[idx].holelist         = NULL;
510 
511       ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
512       for (p = 0; p < numPoints*2; p += 2) {
513         const PetscInt point = points[p];
514         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
515       }
516 
517       tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
518       poly->numberofvertices = numVertices;
519       poly->vertexlist       = new int[poly->numberofvertices];
520       for (v = 0; v < numVertices; ++v) {
521         const PetscInt vIdx = points[v] - vStart;
522         poly->vertexlist[v] = vIdx;
523       }
524       if (label) {ierr = DMLabelGetValue(label, f, &in.facetmarkerlist[idx]);CHKERRQ(ierr);}
525       ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
526     }
527   }
528   if (!rank) {
529     char args[32];
530 
531     /* Take away 'Q' for verbose output */
532     ierr = PetscStrcpy(args, "pqezQ");CHKERRQ(ierr);
533     if (mesh->tetgenOpts) {::tetrahedralize(mesh->tetgenOpts, &in, &out);}
534     else                  {::tetrahedralize(args, &in, &out);}
535   }
536   {
537     DMLabel        glabel      = NULL;
538     const PetscInt numCorners  = 4;
539     const PetscInt numCells    = out.numberoftetrahedra;
540     const PetscInt numVertices = out.numberofpoints;
541     const double   *meshCoords = out.pointlist;
542     int            *cells      = out.tetrahedronlist;
543 
544     ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr);
545     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr);
546     if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);}
547     /* Set labels */
548     for (v = 0; v < numVertices; ++v) {
549       if (out.pointmarkerlist[v]) {
550         if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);}
551       }
552     }
553     if (interpolate) {
554       PetscInt e;
555 
556       for (e = 0; e < out.numberofedges; e++) {
557         if (out.edgemarkerlist[e]) {
558           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
559           const PetscInt *edges;
560           PetscInt        numEdges;
561 
562           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
563           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
564           if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
565           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
566         }
567       }
568       for (f = 0; f < out.numberoftrifaces; f++) {
569         if (out.trifacemarkerlist[f]) {
570           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
571           const PetscInt *faces;
572           PetscInt        numFaces;
573 
574           ierr = DMPlexGetJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
575           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
576           if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);}
577           ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
578         }
579       }
580     }
581     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
582   }
583   PetscFunctionReturn(0);
584 }
585 
586 static PetscErrorCode DMPlexRefine_Tetgen(DM dm, double *maxVolumes, DM *dmRefined)
587 {
588   MPI_Comm       comm;
589   const PetscInt dim       = 3;
590   const char    *labelName = "marker";
591   ::tetgenio     in;
592   ::tetgenio     out;
593   DMLabel        label;
594   PetscInt       vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal;
595   PetscMPIInt    rank;
596   PetscErrorCode ierr;
597 
598   PetscFunctionBegin;
599   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
600   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
601   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
602   ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
603   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
604   ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr);
605 
606   in.numberofpoints = vEnd - vStart;
607   if (in.numberofpoints > 0) {
608     PetscSection coordSection;
609     Vec          coordinates;
610     PetscScalar *array;
611 
612     in.pointlist       = new double[in.numberofpoints*dim];
613     in.pointmarkerlist = new int[in.numberofpoints];
614 
615     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
616     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
617     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
618     for (v = vStart; v < vEnd; ++v) {
619       const PetscInt idx = v - vStart;
620       PetscInt       off, d;
621 
622       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
623       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
624       if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);}
625     }
626     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
627   }
628   ierr  = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
629 
630   in.numberofcorners       = 4;
631   in.numberoftetrahedra    = cEnd - cStart;
632   in.tetrahedronvolumelist = (double*) maxVolumes;
633   if (in.numberoftetrahedra > 0) {
634     in.tetrahedronlist = new int[in.numberoftetrahedra*in.numberofcorners];
635     for (c = cStart; c < cEnd; ++c) {
636       const PetscInt idx      = c - cStart;
637       PetscInt      *closure = NULL;
638       PetscInt       closureSize;
639 
640       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
641       if ((closureSize != 5) && (closureSize != 15)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %D vertices in closure", closureSize);
642       for (v = 0; v < 4; ++v) {
643         in.tetrahedronlist[idx*in.numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
644       }
645       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
646     }
647   }
648   /* TODO: Put in boundary faces with markers */
649   if (!rank) {
650     char args[32];
651 
652     /* Take away 'Q' for verbose output */
653     /*ierr = PetscStrcpy(args, "qezQra");CHKERRQ(ierr); */
654     ierr = PetscStrcpy(args, "qezraVVVV");CHKERRQ(ierr);
655     ::tetrahedralize(args, &in, &out);
656   }
657   in.tetrahedronvolumelist = NULL;
658 
659   {
660     DMLabel        rlabel      = NULL;
661     const PetscInt numCorners  = 4;
662     const PetscInt numCells    = out.numberoftetrahedra;
663     const PetscInt numVertices = out.numberofpoints;
664     const double   *meshCoords = out.pointlist;
665     int            *cells      = out.tetrahedronlist;
666 
667     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
668 
669     ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr);
670     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
671     if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);}
672     /* Set labels */
673     for (v = 0; v < numVertices; ++v) {
674       if (out.pointmarkerlist[v]) {
675         if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);}
676       }
677     }
678     if (interpolate) {
679       PetscInt e, f;
680 
681       for (e = 0; e < out.numberofedges; e++) {
682         if (out.edgemarkerlist[e]) {
683           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
684           const PetscInt *edges;
685           PetscInt        numEdges;
686 
687           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
688           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
689           if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
690           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
691         }
692       }
693       for (f = 0; f < out.numberoftrifaces; f++) {
694         if (out.trifacemarkerlist[f]) {
695           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
696           const PetscInt *faces;
697           PetscInt        numFaces;
698 
699           ierr = DMPlexGetJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
700           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
701           if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);}
702           ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
703         }
704       }
705     }
706     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
707   }
708   PetscFunctionReturn(0);
709 }
710 #endif /* PETSC_HAVE_TETGEN */
711 
712 #if defined(PETSC_HAVE_CTETGEN)
713 #include <ctetgen.h>
714 
715 static PetscErrorCode DMPlexGenerate_CTetgen(DM boundary, PetscBool interpolate, DM *dm)
716 {
717   MPI_Comm       comm;
718   const PetscInt dim       = 3;
719   const char    *labelName = "marker";
720   PLC           *in, *out;
721   DMLabel        label;
722   PetscInt       verbose = 0, vStart, vEnd, v, fStart, fEnd, f;
723   PetscMPIInt    rank;
724   PetscErrorCode ierr;
725 
726   PetscFunctionBegin;
727   ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr);
728   ierr = PetscOptionsGetInt(NULL,((PetscObject) boundary)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr);
729   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
730   ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr);
731   ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr);
732   ierr = PLCCreate(&in);CHKERRQ(ierr);
733   ierr = PLCCreate(&out);CHKERRQ(ierr);
734 
735   in->numberofpoints = vEnd - vStart;
736   if (in->numberofpoints > 0) {
737     PetscSection coordSection;
738     Vec          coordinates;
739     PetscScalar *array;
740 
741     ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr);
742     ierr = PetscMalloc1(in->numberofpoints,       &in->pointmarkerlist);CHKERRQ(ierr);
743     ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr);
744     ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr);
745     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
746     for (v = vStart; v < vEnd; ++v) {
747       const PetscInt idx = v - vStart;
748       PetscInt       off, d, m = -1;
749 
750       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
751       for (d = 0; d < dim; ++d) {
752         in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
753       }
754       if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);}
755 
756       in->pointmarkerlist[idx] = (int) m;
757     }
758     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
759   }
760   ierr  = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr);
761 
762   in->numberoffacets = fEnd - fStart;
763   if (in->numberoffacets > 0) {
764     ierr = PetscMalloc1(in->numberoffacets, &in->facetlist);CHKERRQ(ierr);
765     ierr = PetscMalloc1(in->numberoffacets,   &in->facetmarkerlist);CHKERRQ(ierr);
766     for (f = fStart; f < fEnd; ++f) {
767       const PetscInt idx     = f - fStart;
768       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, m = -1;
769       polygon       *poly;
770 
771       in->facetlist[idx].numberofpolygons = 1;
772 
773       ierr = PetscMalloc1(in->facetlist[idx].numberofpolygons, &in->facetlist[idx].polygonlist);CHKERRQ(ierr);
774 
775       in->facetlist[idx].numberofholes    = 0;
776       in->facetlist[idx].holelist         = NULL;
777 
778       ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
779       for (p = 0; p < numPoints*2; p += 2) {
780         const PetscInt point = points[p];
781         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
782       }
783 
784       poly                   = in->facetlist[idx].polygonlist;
785       poly->numberofvertices = numVertices;
786       ierr                   = PetscMalloc1(poly->numberofvertices, &poly->vertexlist);CHKERRQ(ierr);
787       for (v = 0; v < numVertices; ++v) {
788         const PetscInt vIdx = points[v] - vStart;
789         poly->vertexlist[v] = vIdx;
790       }
791       if (label) {ierr = DMLabelGetValue(label, f, &m);CHKERRQ(ierr);}
792       in->facetmarkerlist[idx] = (int) m;
793       ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
794     }
795   }
796   if (!rank) {
797     TetGenOpts t;
798 
799     ierr        = TetGenOptsInitialize(&t);CHKERRQ(ierr);
800     t.in        = boundary; /* Should go away */
801     t.plc       = 1;
802     t.quality   = 1;
803     t.edgesout  = 1;
804     t.zeroindex = 1;
805     t.quiet     = 1;
806     t.verbose   = verbose;
807     ierr        = TetGenCheckOpts(&t);CHKERRQ(ierr);
808     ierr        = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr);
809   }
810   {
811     DMLabel        glabel      = NULL;
812     const PetscInt numCorners  = 4;
813     const PetscInt numCells    = out->numberoftetrahedra;
814     const PetscInt numVertices = out->numberofpoints;
815     double         *meshCoords;
816     int            *cells      = out->tetrahedronlist;
817 
818     if (sizeof (PetscReal) == sizeof (double)) {
819       meshCoords = (double *) out->pointlist;
820     }
821     else {
822       PetscInt i;
823 
824       ierr = PetscMalloc1(3 * numVertices,&meshCoords);CHKERRQ(ierr);
825       for (i = 0; i < 3 * numVertices; i++) {
826         meshCoords[i] = (PetscReal) out->pointlist[i];
827       }
828     }
829 
830     ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr);
831     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr);
832     if (sizeof (PetscReal) != sizeof (double)) {
833       ierr = PetscFree(meshCoords);CHKERRQ(ierr);
834     }
835     if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);}
836     /* Set labels */
837     for (v = 0; v < numVertices; ++v) {
838       if (out->pointmarkerlist[v]) {
839         if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);}
840       }
841     }
842     if (interpolate) {
843       PetscInt e;
844 
845       for (e = 0; e < out->numberofedges; e++) {
846         if (out->edgemarkerlist[e]) {
847           const PetscInt  vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells};
848           const PetscInt *edges;
849           PetscInt        numEdges;
850 
851           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
852           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
853           if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);}
854           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
855         }
856       }
857       for (f = 0; f < out->numberoftrifaces; f++) {
858         if (out->trifacemarkerlist[f]) {
859           const PetscInt  vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells};
860           const PetscInt *faces;
861           PetscInt        numFaces;
862 
863           ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
864           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
865           if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);}
866           ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
867         }
868       }
869     }
870     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
871   }
872 
873   ierr = PLCDestroy(&in);CHKERRQ(ierr);
874   ierr = PLCDestroy(&out);CHKERRQ(ierr);
875   PetscFunctionReturn(0);
876 }
877 
878 static PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined)
879 {
880   MPI_Comm       comm;
881   const PetscInt dim       = 3;
882   const char    *labelName = "marker";
883   PLC           *in, *out;
884   DMLabel        label;
885   PetscInt       verbose = 0, vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal;
886   PetscMPIInt    rank;
887   PetscErrorCode ierr;
888 
889   PetscFunctionBegin;
890   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
891   ierr = PetscOptionsGetInt(NULL,((PetscObject) dm)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr);
892   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
893   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
894   ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
895   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
896   ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr);
897   ierr = PLCCreate(&in);CHKERRQ(ierr);
898   ierr = PLCCreate(&out);CHKERRQ(ierr);
899 
900   in->numberofpoints = vEnd - vStart;
901   if (in->numberofpoints > 0) {
902     PetscSection coordSection;
903     Vec          coordinates;
904     PetscScalar *array;
905 
906     ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr);
907     ierr = PetscMalloc1(in->numberofpoints,       &in->pointmarkerlist);CHKERRQ(ierr);
908     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
909     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
910     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
911     for (v = vStart; v < vEnd; ++v) {
912       const PetscInt idx = v - vStart;
913       PetscInt       off, d, m = -1;
914 
915       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
916       for (d = 0; d < dim; ++d) {
917         in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
918       }
919       if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);}
920 
921       in->pointmarkerlist[idx] = (int) m;
922     }
923     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
924   }
925   ierr  = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
926 
927   in->numberofcorners       = 4;
928   in->numberoftetrahedra    = cEnd - cStart;
929   in->tetrahedronvolumelist = maxVolumes;
930   if (in->numberoftetrahedra > 0) {
931     ierr = PetscMalloc1(in->numberoftetrahedra*in->numberofcorners, &in->tetrahedronlist);CHKERRQ(ierr);
932     for (c = cStart; c < cEnd; ++c) {
933       const PetscInt idx      = c - cStart;
934       PetscInt      *closure = NULL;
935       PetscInt       closureSize;
936 
937       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
938       if ((closureSize != 5) && (closureSize != 15)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %D vertices in closure", closureSize);
939       for (v = 0; v < 4; ++v) {
940         in->tetrahedronlist[idx*in->numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
941       }
942       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
943     }
944   }
945   if (!rank) {
946     TetGenOpts t;
947 
948     ierr = TetGenOptsInitialize(&t);CHKERRQ(ierr);
949 
950     t.in        = dm; /* Should go away */
951     t.refine    = 1;
952     t.varvolume = 1;
953     t.quality   = 1;
954     t.edgesout  = 1;
955     t.zeroindex = 1;
956     t.quiet     = 1;
957     t.verbose   = verbose; /* Change this */
958 
959     ierr = TetGenCheckOpts(&t);CHKERRQ(ierr);
960     ierr = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr);
961   }
962   {
963     DMLabel        rlabel      = NULL;
964     const PetscInt numCorners  = 4;
965     const PetscInt numCells    = out->numberoftetrahedra;
966     const PetscInt numVertices = out->numberofpoints;
967     double         *meshCoords;
968     int            *cells      = out->tetrahedronlist;
969     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
970 
971     if (sizeof (PetscReal) == sizeof (double)) {
972       meshCoords = (double *) out->pointlist;
973     }
974     else {
975       PetscInt i;
976 
977       ierr = PetscMalloc1(3 * numVertices,&meshCoords);CHKERRQ(ierr);
978       for (i = 0; i < 3 * numVertices; i++) {
979         meshCoords[i] = (PetscReal) out->pointlist[i];
980       }
981     }
982 
983     ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr);
984     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
985     if (sizeof (PetscReal) != sizeof (double)) {
986       ierr = PetscFree(meshCoords);CHKERRQ(ierr);
987     }
988     if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);}
989     /* Set labels */
990     for (v = 0; v < numVertices; ++v) {
991       if (out->pointmarkerlist[v]) {
992         if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);}
993       }
994     }
995     if (interpolate) {
996       PetscInt e, f;
997 
998       for (e = 0; e < out->numberofedges; e++) {
999         if (out->edgemarkerlist[e]) {
1000           const PetscInt  vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells};
1001           const PetscInt *edges;
1002           PetscInt        numEdges;
1003 
1004           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1005           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
1006           if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);}
1007           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1008         }
1009       }
1010       for (f = 0; f < out->numberoftrifaces; f++) {
1011         if (out->trifacemarkerlist[f]) {
1012           const PetscInt  vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells};
1013           const PetscInt *faces;
1014           PetscInt        numFaces;
1015 
1016           ierr = DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1017           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
1018           if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);}
1019           ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1020         }
1021       }
1022     }
1023     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
1024   }
1025   ierr = PLCDestroy(&in);CHKERRQ(ierr);
1026   ierr = PLCDestroy(&out);CHKERRQ(ierr);
1027   PetscFunctionReturn(0);
1028 }
1029 #endif /* PETSC_HAVE_CTETGEN */
1030 
1031 /*@C
1032   DMPlexGenerate - Generates a mesh.
1033 
1034   Not Collective
1035 
1036   Input Parameters:
1037 + boundary - The DMPlex boundary object
1038 . name - The mesh generation package name
1039 - interpolate - Flag to create intermediate mesh elements
1040 
1041   Output Parameter:
1042 . mesh - The DMPlex object
1043 
1044   Level: intermediate
1045 
1046 .keywords: mesh, elements
1047 .seealso: DMPlexCreate(), DMRefine()
1048 @*/
1049 PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
1050 {
1051   PetscInt       dim;
1052   char           genname[1024];
1053   PetscBool      isTriangle = PETSC_FALSE, isTetgen = PETSC_FALSE, isCTetgen = PETSC_FALSE, flg;
1054   PetscErrorCode ierr;
1055 
1056   PetscFunctionBegin;
1057   PetscValidHeaderSpecific(boundary, DM_CLASSID, 1);
1058   PetscValidLogicalCollectiveBool(boundary, interpolate, 2);
1059   ierr = DMGetDimension(boundary, &dim);CHKERRQ(ierr);
1060   ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr);
1061   if (flg) name = genname;
1062   if (name) {
1063     ierr = PetscStrcmp(name, "triangle", &isTriangle);CHKERRQ(ierr);
1064     ierr = PetscStrcmp(name, "tetgen",   &isTetgen);CHKERRQ(ierr);
1065     ierr = PetscStrcmp(name, "ctetgen",  &isCTetgen);CHKERRQ(ierr);
1066   }
1067   switch (dim) {
1068   case 1:
1069     if (!name || isTriangle) {
1070 #if defined(PETSC_HAVE_TRIANGLE)
1071       ierr = DMPlexGenerate_Triangle(boundary, interpolate, mesh);CHKERRQ(ierr);
1072 #else
1073       SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Mesh generation needs external package support.\nPlease reconfigure with --download-triangle.");
1074 #endif
1075     } else SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Unknown 2D mesh generation package %s", name);
1076     break;
1077   case 2:
1078     if (!name || isCTetgen) {
1079 #if defined(PETSC_HAVE_CTETGEN)
1080       ierr = DMPlexGenerate_CTetgen(boundary, interpolate, mesh);CHKERRQ(ierr);
1081 #else
1082       SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "CTetgen needs external package support.\nPlease reconfigure with --download-ctetgen.");
1083 #endif
1084     } else if (isTetgen) {
1085 #if defined(PETSC_HAVE_TETGEN)
1086       ierr = DMPlexGenerate_Tetgen(boundary, interpolate, mesh);CHKERRQ(ierr);
1087 #else
1088       SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Tetgen needs external package support.\nPlease reconfigure with --download-tetgen.");
1089 #endif
1090     } else SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Unknown 3D mesh generation package %s", name);
1091     break;
1092   default:
1093     SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Mesh generation for a dimension %d boundary is not supported.", dim);
1094   }
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 #if defined(PETSC_HAVE_TRIANGLE) || defined(PETSC_HAVE_CTETGEN) || defined(PETSC_HAVE_TETGEN)
1099 static PetscErrorCode DMRefine_Plex_Label(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal maxVolumes[])
1100 {
1101   PetscInt       dim, c;
1102   PetscReal      refRatio;
1103   PetscErrorCode ierr;
1104 
1105   PetscFunctionBegin;
1106   ierr     = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1107   refRatio = (PetscReal) ((PetscInt) 1 << dim);
1108   for (c = cStart; c < cEnd; c++) {
1109     PetscReal vol;
1110     PetscInt  i, closureSize = 0;
1111     PetscInt  *closure = NULL;
1112     PetscBool anyRefine  = PETSC_FALSE;
1113     PetscBool anyCoarsen = PETSC_FALSE;
1114     PetscBool anyKeep    = PETSC_FALSE;
1115 
1116     ierr = DMPlexComputeCellGeometryFVM(dm,c,&vol,NULL,NULL);CHKERRQ(ierr);
1117     maxVolumes[c - cStart] = vol;
1118     ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1119     for (i = 0; i < closureSize; i++) {
1120       PetscInt point = closure[2 * i], refFlag;
1121 
1122       ierr = DMLabelGetValue(adaptLabel,point,&refFlag);CHKERRQ(ierr);
1123       switch (refFlag) {
1124       case DM_ADAPT_REFINE:
1125         anyRefine = PETSC_TRUE;
1126         break;
1127       case DM_ADAPT_COARSEN:
1128         anyCoarsen = PETSC_TRUE;
1129         break;
1130       case DM_ADAPT_KEEP:
1131         anyKeep = PETSC_TRUE;
1132         break;
1133       case DM_ADAPT_DETERMINE:
1134         break;
1135       default:
1136         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n",refFlag);
1137         break;
1138       }
1139       if (anyRefine) break;
1140     }
1141     ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1142     if (anyRefine) {
1143       maxVolumes[c - cStart] = vol / refRatio;
1144     } else if (anyKeep) {
1145       maxVolumes[c - cStart] = vol;
1146     } else if (anyCoarsen) {
1147       maxVolumes[c - cStart] = vol * refRatio;
1148     }
1149   }
1150   PetscFunctionReturn(0);
1151 }
1152 #endif
1153 
1154 static PetscErrorCode DMRefine_Plex_Internal(DM dm, MPI_Comm comm, DMLabel adaptLabel, DM *dmRefined)
1155 {
1156   PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *);
1157   PetscReal        refinementLimit;
1158   PetscInt         dim, cStart, cEnd;
1159   char             genname[1024], *name = NULL;
1160   PetscBool        isUniform, isTriangle = PETSC_FALSE, isTetgen = PETSC_FALSE, isCTetgen = PETSC_FALSE, flg, localized;
1161   PetscErrorCode   ierr;
1162 
1163   PetscFunctionBegin;
1164   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1165   ierr = DMPlexGetRefinementUniform(dm, &isUniform);CHKERRQ(ierr);
1166   if (isUniform) {
1167     CellRefiner cellRefiner;
1168 
1169     ierr = DMPlexGetCellRefiner_Internal(dm, &cellRefiner);CHKERRQ(ierr);
1170     ierr = DMPlexRefineUniform_Internal(dm, cellRefiner, dmRefined);CHKERRQ(ierr);
1171     ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr);
1172     if (localized) {
1173       ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);
1174     }
1175     PetscFunctionReturn(0);
1176   }
1177   ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr);
1178   ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr);
1179   if (refinementLimit == 0.0 && !refinementFunc) PetscFunctionReturn(0);
1180   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1181   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1182   ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr);
1183   if (flg) name = genname;
1184   if (name) {
1185     ierr = PetscStrcmp(name, "triangle", &isTriangle);CHKERRQ(ierr);
1186     ierr = PetscStrcmp(name, "tetgen",   &isTetgen);CHKERRQ(ierr);
1187     ierr = PetscStrcmp(name, "ctetgen",  &isCTetgen);CHKERRQ(ierr);
1188   }
1189   switch (dim) {
1190   case 2:
1191     if (!name || isTriangle) {
1192 #if defined(PETSC_HAVE_TRIANGLE)
1193       double  *maxVolumes;
1194       PetscInt c;
1195 
1196       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
1197       if (adaptLabel) {
1198         ierr = DMRefine_Plex_Label(dm,adaptLabel,cStart,cEnd,maxVolumes);CHKERRQ(ierr);
1199       } else if (refinementFunc) {
1200         for (c = cStart; c < cEnd; ++c) {
1201           PetscReal vol, centroid[3];
1202           PetscReal maxVol;
1203 
1204           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
1205           ierr = (*refinementFunc)(centroid, &maxVol);CHKERRQ(ierr);
1206           maxVolumes[c - cStart] = (double) maxVol;
1207         }
1208       } else {
1209         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
1210       }
1211       ierr = DMPlexRefine_Triangle(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
1212       ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
1213 #else
1214       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement needs external package support.\nPlease reconfigure with --download-triangle.");
1215 #endif
1216     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 2D mesh generation package %s", name);
1217     break;
1218   case 3:
1219     if (!name || isCTetgen) {
1220 #if defined(PETSC_HAVE_CTETGEN)
1221       PetscReal *maxVolumes;
1222       PetscInt   c;
1223 
1224       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
1225       if (adaptLabel) {
1226         ierr = DMRefine_Plex_Label(dm,adaptLabel,cStart,cEnd,maxVolumes);CHKERRQ(ierr);
1227       } else if (refinementFunc) {
1228         for (c = cStart; c < cEnd; ++c) {
1229           PetscReal vol, centroid[3];
1230 
1231           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
1232           ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr);
1233         }
1234       } else {
1235         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
1236       }
1237       ierr = DMPlexRefine_CTetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
1238 #else
1239       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "CTetgen needs external package support.\nPlease reconfigure with --download-ctetgen.");
1240 #endif
1241     } else if (isTetgen) {
1242 #if defined(PETSC_HAVE_TETGEN)
1243       double  *maxVolumes;
1244       PetscInt c;
1245 
1246       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
1247       if (adaptLabel) {
1248         ierr = DMRefine_Plex_Label(dm,adaptLabel,cStart,cEnd,maxVolumes);CHKERRQ(ierr);
1249       } else if (refinementFunc) {
1250         for (c = cStart; c < cEnd; ++c) {
1251           PetscReal vol, centroid[3];
1252 
1253           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
1254           ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr);
1255         }
1256       } else {
1257         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
1258       }
1259       ierr = DMPlexRefine_Tetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
1260 #else
1261       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tetgen needs external package support.\nPlease reconfigure with --with-c-language=cxx --download-tetgen.");
1262 #endif
1263     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 3D mesh generation package %s", name);
1264     break;
1265   default:
1266     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %d is not supported.", dim);
1267   }
1268   ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr);
1269   if (localized) {
1270     ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);
1271   }
1272   PetscFunctionReturn(0);
1273 }
1274 
1275 PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined)
1276 {
1277   PetscErrorCode ierr;
1278 
1279   PetscFunctionBegin;
1280   ierr = DMRefine_Plex_Internal(dm,comm,NULL,dmRefined);CHKERRQ(ierr);
1281   PetscFunctionReturn(0);
1282 }
1283 
1284 PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted)
1285 {
1286   PetscInt       defFlag, minFlag, maxFlag, numFlags, i;
1287   const PetscInt *flags;
1288   IS             flagIS;
1289   PetscErrorCode ierr;
1290 
1291   PetscFunctionBegin;
1292   ierr = DMLabelGetDefaultValue(adaptLabel,&defFlag);CHKERRQ(ierr);
1293   minFlag = defFlag;
1294   maxFlag = defFlag;
1295   ierr = DMLabelGetValueIS(adaptLabel,&flagIS);CHKERRQ(ierr);
1296   ierr = ISGetLocalSize(flagIS,&numFlags);CHKERRQ(ierr);
1297   ierr = ISGetIndices(flagIS,&flags);CHKERRQ(ierr);
1298   for (i = 0; i < numFlags; i++) {
1299     PetscInt flag = flags[i];
1300 
1301     minFlag = PetscMin(minFlag,flag);
1302     maxFlag = PetscMax(maxFlag,flag);
1303   }
1304   ierr = ISRestoreIndices(flagIS,&flags);CHKERRQ(ierr);
1305   ierr = ISDestroy(&flagIS);CHKERRQ(ierr);
1306   {
1307     PetscInt minMaxFlag[2], minMaxFlagGlobal[2];
1308 
1309     minMaxFlag[0] = minFlag;
1310     minMaxFlag[1] = -maxFlag;
1311     ierr = MPI_Allreduce(minMaxFlag,minMaxFlagGlobal,2,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
1312     minFlag = minMaxFlagGlobal[0];
1313     maxFlag = -minMaxFlagGlobal[1];
1314   }
1315   if (minFlag == maxFlag) {
1316     switch (minFlag) {
1317     case DM_ADAPT_DETERMINE:
1318       *dmAdapted = NULL;
1319       break;
1320     case DM_ADAPT_REFINE:
1321       ierr = DMPlexSetRefinementUniform(dm,PETSC_TRUE);CHKERRQ(ierr);
1322       ierr = DMRefine(dm,MPI_COMM_NULL,dmAdapted);CHKERRQ(ierr);
1323       break;
1324     case DM_ADAPT_COARSEN:
1325       ierr = DMCoarsen(dm,MPI_COMM_NULL,dmAdapted);CHKERRQ(ierr);
1326       break;
1327     default:
1328       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n",minFlag);
1329       break;
1330     }
1331   } else {
1332     ierr = DMPlexSetRefinementUniform(dm,PETSC_FALSE);CHKERRQ(ierr);
1333     ierr = DMRefine_Plex_Internal(dm,MPI_COMM_NULL,adaptLabel,dmAdapted);CHKERRQ(ierr);
1334   }
1335   PetscFunctionReturn(0);
1336 }
1337 
1338 PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[])
1339 {
1340   DM             cdm = dm;
1341   PetscInt       r;
1342   PetscBool      isUniform, localized;
1343   PetscErrorCode ierr;
1344 
1345   PetscFunctionBegin;
1346   ierr = DMPlexGetRefinementUniform(dm, &isUniform);CHKERRQ(ierr);
1347   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1348   if (!isUniform) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Non-uniform refinement is incompatible with the hierarchy");
1349   for (r = 0; r < nlevels; ++r) {
1350     CellRefiner cellRefiner;
1351 
1352     ierr = DMPlexGetCellRefiner_Internal(cdm, &cellRefiner);CHKERRQ(ierr);
1353     ierr = DMPlexRefineUniform_Internal(cdm, cellRefiner, &dmRefined[r]);CHKERRQ(ierr);
1354     ierr = DMCopyBoundary(cdm, dmRefined[r]);CHKERRQ(ierr);
1355     if (localized) {
1356       ierr = DMLocalizeCoordinates(dmRefined[r]);CHKERRQ(ierr);
1357     }
1358     ierr = DMSetCoarseDM(dmRefined[r], cdm);CHKERRQ(ierr);
1359     ierr = DMPlexSetRegularRefinement(dmRefined[r], PETSC_TRUE);CHKERRQ(ierr);
1360     cdm  = dmRefined[r];
1361   }
1362   PetscFunctionReturn(0);
1363 }
1364