xref: /petsc/src/dm/impls/plex/plexgenerate.c (revision 84572febd70de60755b1f0d1e6de8e999bfb2bb0)
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 /*@
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 /*@
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 #if 0
555       PetscInt e;
556 
557       /* This check is never actually executed for ctetgen (which never returns edgemarkers) and seems to be broken for
558        * tetgen */
559       for (e = 0; e < out.numberofedges; e++) {
560         if (out.edgemarkerlist[e]) {
561           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
562           const PetscInt *edges;
563           PetscInt        numEdges;
564 
565           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
566           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
567           if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
568           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
569         }
570       }
571 #endif
572       for (f = 0; f < out.numberoftrifaces; f++) {
573         if (out.trifacemarkerlist[f]) {
574           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
575           const PetscInt *faces;
576           PetscInt        numFaces;
577 
578           ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
579           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
580           if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);}
581           ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
582         }
583       }
584     }
585     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
586   }
587   PetscFunctionReturn(0);
588 }
589 
590 static PetscErrorCode DMPlexRefine_Tetgen(DM dm, double *maxVolumes, DM *dmRefined)
591 {
592   MPI_Comm       comm;
593   const PetscInt dim       = 3;
594   const char    *labelName = "marker";
595   ::tetgenio     in;
596   ::tetgenio     out;
597   DMLabel        label;
598   PetscInt       vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal;
599   PetscMPIInt    rank;
600   PetscErrorCode ierr;
601 
602   PetscFunctionBegin;
603   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
604   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
605   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
606   ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
607   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
608   ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr);
609 
610   in.numberofpoints = vEnd - vStart;
611   if (in.numberofpoints > 0) {
612     PetscSection coordSection;
613     Vec          coordinates;
614     PetscScalar *array;
615 
616     in.pointlist       = new double[in.numberofpoints*dim];
617     in.pointmarkerlist = new int[in.numberofpoints];
618 
619     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
620     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
621     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
622     for (v = vStart; v < vEnd; ++v) {
623       const PetscInt idx = v - vStart;
624       PetscInt       off, d;
625 
626       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
627       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
628       if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);}
629     }
630     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
631   }
632   ierr  = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
633 
634   in.numberofcorners       = 4;
635   in.numberoftetrahedra    = cEnd - cStart;
636   in.tetrahedronvolumelist = (double*) maxVolumes;
637   if (in.numberoftetrahedra > 0) {
638     in.tetrahedronlist = new int[in.numberoftetrahedra*in.numberofcorners];
639     for (c = cStart; c < cEnd; ++c) {
640       const PetscInt idx      = c - cStart;
641       PetscInt      *closure = NULL;
642       PetscInt       closureSize;
643 
644       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
645       if ((closureSize != 5) && (closureSize != 15)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %D vertices in closure", closureSize);
646       for (v = 0; v < 4; ++v) {
647         in.tetrahedronlist[idx*in.numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
648       }
649       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
650     }
651   }
652   /* TODO: Put in boundary faces with markers */
653   if (!rank) {
654     char args[32];
655 
656 #if 1
657     /* Take away 'Q' for verbose output */
658     ierr = PetscStrcpy(args, "qezQra");CHKERRQ(ierr);
659 #else
660     ierr = PetscStrcpy(args, "qezraVVVV");CHKERRQ(ierr);
661 #endif
662     ::tetrahedralize(args, &in, &out);
663   }
664   in.tetrahedronvolumelist = NULL;
665 
666   {
667     DMLabel        rlabel      = NULL;
668     const PetscInt numCorners  = 4;
669     const PetscInt numCells    = out.numberoftetrahedra;
670     const PetscInt numVertices = out.numberofpoints;
671     const double   *meshCoords = out.pointlist;
672     int            *cells      = out.tetrahedronlist;
673 
674     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
675 
676     ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr);
677     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
678     if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);}
679     /* Set labels */
680     for (v = 0; v < numVertices; ++v) {
681       if (out.pointmarkerlist[v]) {
682         if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);}
683       }
684     }
685     if (interpolate) {
686       PetscInt f;
687 #if 0
688       PetscInt e;
689 
690       for (e = 0; e < out.numberofedges; e++) {
691         if (out.edgemarkerlist[e]) {
692           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
693           const PetscInt *edges;
694           PetscInt        numEdges;
695 
696           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
697           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
698           if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
699           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
700         }
701       }
702 #endif
703       for (f = 0; f < out.numberoftrifaces; f++) {
704         if (out.trifacemarkerlist[f]) {
705           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
706           const PetscInt *faces;
707           PetscInt        numFaces;
708 
709           ierr = DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
710           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
711           if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);}
712           ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
713         }
714       }
715     }
716     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
717   }
718   PetscFunctionReturn(0);
719 }
720 #endif /* PETSC_HAVE_TETGEN */
721 
722 #if defined(PETSC_HAVE_CTETGEN)
723 #include <ctetgen.h>
724 
725 static PetscErrorCode DMPlexGenerate_CTetgen(DM boundary, PetscBool interpolate, DM *dm)
726 {
727   MPI_Comm       comm;
728   const PetscInt dim       = 3;
729   const char    *labelName = "marker";
730   PLC           *in, *out;
731   DMLabel        label;
732   PetscInt       verbose = 0, vStart, vEnd, v, fStart, fEnd, f;
733   PetscMPIInt    rank;
734   PetscErrorCode ierr;
735 
736   PetscFunctionBegin;
737   ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr);
738   ierr = PetscOptionsGetInt(NULL,((PetscObject) boundary)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr);
739   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
740   ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr);
741   ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr);
742   ierr = PLCCreate(&in);CHKERRQ(ierr);
743   ierr = PLCCreate(&out);CHKERRQ(ierr);
744 
745   in->numberofpoints = vEnd - vStart;
746   if (in->numberofpoints > 0) {
747     PetscSection coordSection;
748     Vec          coordinates;
749     PetscScalar *array;
750 
751     ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr);
752     ierr = PetscMalloc1(in->numberofpoints,       &in->pointmarkerlist);CHKERRQ(ierr);
753     ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr);
754     ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr);
755     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
756     for (v = vStart; v < vEnd; ++v) {
757       const PetscInt idx = v - vStart;
758       PetscInt       off, d, m = -1;
759 
760       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
761       for (d = 0; d < dim; ++d) {
762         in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
763       }
764       if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);}
765 
766       in->pointmarkerlist[idx] = (int) m;
767     }
768     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
769   }
770   ierr  = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr);
771 
772   in->numberoffacets = fEnd - fStart;
773   if (in->numberoffacets > 0) {
774     ierr = PetscMalloc1(in->numberoffacets, &in->facetlist);CHKERRQ(ierr);
775     ierr = PetscMalloc1(in->numberoffacets,   &in->facetmarkerlist);CHKERRQ(ierr);
776     for (f = fStart; f < fEnd; ++f) {
777       const PetscInt idx     = f - fStart;
778       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, m = -1;
779       polygon       *poly;
780 
781       in->facetlist[idx].numberofpolygons = 1;
782 
783       ierr = PetscMalloc1(in->facetlist[idx].numberofpolygons, &in->facetlist[idx].polygonlist);CHKERRQ(ierr);
784 
785       in->facetlist[idx].numberofholes    = 0;
786       in->facetlist[idx].holelist         = NULL;
787 
788       ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
789       for (p = 0; p < numPoints*2; p += 2) {
790         const PetscInt point = points[p];
791         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
792       }
793 
794       poly                   = in->facetlist[idx].polygonlist;
795       poly->numberofvertices = numVertices;
796       ierr                   = PetscMalloc1(poly->numberofvertices, &poly->vertexlist);CHKERRQ(ierr);
797       for (v = 0; v < numVertices; ++v) {
798         const PetscInt vIdx = points[v] - vStart;
799         poly->vertexlist[v] = vIdx;
800       }
801       if (label) {ierr = DMLabelGetValue(label, f, &m);CHKERRQ(ierr);}
802       in->facetmarkerlist[idx] = (int) m;
803       ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
804     }
805   }
806   if (!rank) {
807     TetGenOpts t;
808 
809     ierr        = TetGenOptsInitialize(&t);CHKERRQ(ierr);
810     t.in        = boundary; /* Should go away */
811     t.plc       = 1;
812     t.quality   = 1;
813     t.edgesout  = 1;
814     t.zeroindex = 1;
815     t.quiet     = 1;
816     t.verbose   = verbose;
817     ierr        = TetGenCheckOpts(&t);CHKERRQ(ierr);
818     ierr        = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr);
819   }
820   {
821     DMLabel        glabel      = NULL;
822     const PetscInt numCorners  = 4;
823     const PetscInt numCells    = out->numberoftetrahedra;
824     const PetscInt numVertices = out->numberofpoints;
825     double         *meshCoords;
826     int            *cells      = out->tetrahedronlist;
827 
828     if (sizeof (PetscReal) == sizeof (double)) {
829       meshCoords = (double *) out->pointlist;
830     }
831     else {
832       PetscInt i;
833 
834       ierr = PetscMalloc1(3 * numVertices,&meshCoords);CHKERRQ(ierr);
835       for (i = 0; i < 3 * numVertices; i++) {
836         meshCoords[i] = (PetscReal) out->pointlist[i];
837       }
838     }
839 
840     ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr);
841     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr);
842     if (sizeof (PetscReal) != sizeof (double)) {
843       ierr = PetscFree(meshCoords);CHKERRQ(ierr);
844     }
845     if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);}
846     /* Set labels */
847     for (v = 0; v < numVertices; ++v) {
848       if (out->pointmarkerlist[v]) {
849         if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);}
850       }
851     }
852     if (interpolate) {
853       PetscInt e;
854 
855       for (e = 0; e < out->numberofedges; e++) {
856         if (out->edgemarkerlist[e]) {
857           const PetscInt  vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells};
858           const PetscInt *edges;
859           PetscInt        numEdges;
860 
861           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
862           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
863           if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);}
864           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
865         }
866       }
867       for (f = 0; f < out->numberoftrifaces; f++) {
868         if (out->trifacemarkerlist[f]) {
869           const PetscInt  vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells};
870           const PetscInt *faces;
871           PetscInt        numFaces;
872 
873           ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
874           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
875           if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);}
876           ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
877         }
878       }
879     }
880     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
881   }
882 
883   ierr = PLCDestroy(&in);CHKERRQ(ierr);
884   ierr = PLCDestroy(&out);CHKERRQ(ierr);
885   PetscFunctionReturn(0);
886 }
887 
888 static PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined)
889 {
890   MPI_Comm       comm;
891   const PetscInt dim       = 3;
892   const char    *labelName = "marker";
893   PLC           *in, *out;
894   DMLabel        label;
895   PetscInt       verbose = 0, vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal;
896   PetscMPIInt    rank;
897   PetscErrorCode ierr;
898 
899   PetscFunctionBegin;
900   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
901   ierr = PetscOptionsGetInt(NULL,((PetscObject) dm)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr);
902   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
903   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
904   ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
905   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
906   ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr);
907   ierr = PLCCreate(&in);CHKERRQ(ierr);
908   ierr = PLCCreate(&out);CHKERRQ(ierr);
909 
910   in->numberofpoints = vEnd - vStart;
911   if (in->numberofpoints > 0) {
912     PetscSection coordSection;
913     Vec          coordinates;
914     PetscScalar *array;
915 
916     ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr);
917     ierr = PetscMalloc1(in->numberofpoints,       &in->pointmarkerlist);CHKERRQ(ierr);
918     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
919     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
920     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
921     for (v = vStart; v < vEnd; ++v) {
922       const PetscInt idx = v - vStart;
923       PetscInt       off, d, m = -1;
924 
925       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
926       for (d = 0; d < dim; ++d) {
927         in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
928       }
929       if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);}
930 
931       in->pointmarkerlist[idx] = (int) m;
932     }
933     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
934   }
935   ierr  = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
936 
937   in->numberofcorners       = 4;
938   in->numberoftetrahedra    = cEnd - cStart;
939   in->tetrahedronvolumelist = maxVolumes;
940   if (in->numberoftetrahedra > 0) {
941     ierr = PetscMalloc1(in->numberoftetrahedra*in->numberofcorners, &in->tetrahedronlist);CHKERRQ(ierr);
942     for (c = cStart; c < cEnd; ++c) {
943       const PetscInt idx      = c - cStart;
944       PetscInt      *closure = NULL;
945       PetscInt       closureSize;
946 
947       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
948       if ((closureSize != 5) && (closureSize != 15)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %D vertices in closure", closureSize);
949       for (v = 0; v < 4; ++v) {
950         in->tetrahedronlist[idx*in->numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
951       }
952       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
953     }
954   }
955   if (!rank) {
956     TetGenOpts t;
957 
958     ierr = TetGenOptsInitialize(&t);CHKERRQ(ierr);
959 
960     t.in        = dm; /* Should go away */
961     t.refine    = 1;
962     t.varvolume = 1;
963     t.quality   = 1;
964     t.edgesout  = 1;
965     t.zeroindex = 1;
966     t.quiet     = 1;
967     t.verbose   = verbose; /* Change this */
968 
969     ierr = TetGenCheckOpts(&t);CHKERRQ(ierr);
970     ierr = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr);
971   }
972   in->tetrahedronvolumelist = NULL;
973   {
974     DMLabel        rlabel      = NULL;
975     const PetscInt numCorners  = 4;
976     const PetscInt numCells    = out->numberoftetrahedra;
977     const PetscInt numVertices = out->numberofpoints;
978     double         *meshCoords;
979     int            *cells      = out->tetrahedronlist;
980     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
981 
982     if (sizeof (PetscReal) == sizeof (double)) {
983       meshCoords = (double *) out->pointlist;
984     }
985     else {
986       PetscInt i;
987 
988       ierr = PetscMalloc1(3 * numVertices,&meshCoords);CHKERRQ(ierr);
989       for (i = 0; i < 3 * numVertices; i++) {
990         meshCoords[i] = (PetscReal) out->pointlist[i];
991       }
992     }
993 
994     ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr);
995     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
996     if (sizeof (PetscReal) != sizeof (double)) {
997       ierr = PetscFree(meshCoords);CHKERRQ(ierr);
998     }
999     if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);}
1000     /* Set labels */
1001     for (v = 0; v < numVertices; ++v) {
1002       if (out->pointmarkerlist[v]) {
1003         if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);}
1004       }
1005     }
1006     if (interpolate) {
1007       PetscInt e, f;
1008 
1009       for (e = 0; e < out->numberofedges; e++) {
1010         if (out->edgemarkerlist[e]) {
1011           const PetscInt  vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells};
1012           const PetscInt *edges;
1013           PetscInt        numEdges;
1014 
1015           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1016           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
1017           if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);}
1018           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1019         }
1020       }
1021       for (f = 0; f < out->numberoftrifaces; f++) {
1022         if (out->trifacemarkerlist[f]) {
1023           const PetscInt  vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells};
1024           const PetscInt *faces;
1025           PetscInt        numFaces;
1026 
1027           ierr = DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1028           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
1029           if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);}
1030           ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1031         }
1032       }
1033     }
1034     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
1035   }
1036   ierr = PLCDestroy(&in);CHKERRQ(ierr);
1037   ierr = PLCDestroy(&out);CHKERRQ(ierr);
1038   PetscFunctionReturn(0);
1039 }
1040 #endif /* PETSC_HAVE_CTETGEN */
1041 
1042 /*@C
1043   DMPlexGenerate - Generates a mesh.
1044 
1045   Not Collective
1046 
1047   Input Parameters:
1048 + boundary - The DMPlex boundary object
1049 . name - The mesh generation package name
1050 - interpolate - Flag to create intermediate mesh elements
1051 
1052   Output Parameter:
1053 . mesh - The DMPlex object
1054 
1055   Level: intermediate
1056 
1057 .keywords: mesh, elements
1058 .seealso: DMPlexCreate(), DMRefine()
1059 @*/
1060 PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
1061 {
1062   PetscInt       dim;
1063   char           genname[1024];
1064   PetscBool      isTriangle = PETSC_FALSE, isTetgen = PETSC_FALSE, isCTetgen = PETSC_FALSE, flg;
1065   PetscErrorCode ierr;
1066 
1067   PetscFunctionBegin;
1068   PetscValidHeaderSpecific(boundary, DM_CLASSID, 1);
1069   PetscValidLogicalCollectiveBool(boundary, interpolate, 2);
1070   ierr = DMGetDimension(boundary, &dim);CHKERRQ(ierr);
1071   ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr);
1072   if (flg) name = genname;
1073   if (name) {
1074     ierr = PetscStrcmp(name, "triangle", &isTriangle);CHKERRQ(ierr);
1075     ierr = PetscStrcmp(name, "tetgen",   &isTetgen);CHKERRQ(ierr);
1076     ierr = PetscStrcmp(name, "ctetgen",  &isCTetgen);CHKERRQ(ierr);
1077   }
1078   switch (dim) {
1079   case 1:
1080     if (!name || isTriangle) {
1081 #if defined(PETSC_HAVE_TRIANGLE)
1082       ierr = DMPlexGenerate_Triangle(boundary, interpolate, mesh);CHKERRQ(ierr);
1083 #else
1084       SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Mesh generation needs external package support.\nPlease reconfigure with --download-triangle.");
1085 #endif
1086     } else SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Unknown 2D mesh generation package %s", name);
1087     break;
1088   case 2:
1089     if (!name) {
1090 #if defined(PETSC_HAVE_CTETGEN)
1091       ierr = DMPlexGenerate_CTetgen(boundary, interpolate, mesh);CHKERRQ(ierr);
1092 #elif defined(PETSC_HAVE_TETGEN)
1093       ierr = DMPlexGenerate_Tetgen(boundary, interpolate, mesh);CHKERRQ(ierr);
1094 #else
1095       SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "External package CTetgen or Tetgen needed.\nPlease reconfigure with '--download-ctetgen' or '--with-clanguage=cxx --download-tetgen'.");
1096 #endif
1097     } else if (isCTetgen) {
1098 #if defined(PETSC_HAVE_CTETGEN)
1099       ierr = DMPlexGenerate_CTetgen(boundary, interpolate, mesh);CHKERRQ(ierr);
1100 #else
1101       SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "CTetgen needs external package support.\nPlease reconfigure with --download-ctetgen.");
1102 #endif
1103     } else if (isTetgen) {
1104 #if defined(PETSC_HAVE_TETGEN)
1105       ierr = DMPlexGenerate_Tetgen(boundary, interpolate, mesh);CHKERRQ(ierr);
1106 #else
1107       SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Tetgen needs external package support.\nPlease reconfigure with --download-tetgen.");
1108 #endif
1109     } else SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Unknown 3D mesh generation package %s", name);
1110     break;
1111   default:
1112     SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Mesh generation for a dimension %d boundary is not supported.", dim);
1113   }
1114   PetscFunctionReturn(0);
1115 }
1116 
1117 #if defined(PETSC_HAVE_TRIANGLE) || defined(PETSC_HAVE_CTETGEN) || defined(PETSC_HAVE_TETGEN)
1118 static PetscErrorCode DMRefine_Plex_Label(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal maxVolumes[])
1119 {
1120   PetscInt       dim, c;
1121   PetscReal      refRatio;
1122   PetscErrorCode ierr;
1123 
1124   PetscFunctionBegin;
1125   ierr     = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1126   refRatio = (PetscReal) ((PetscInt) 1 << dim);
1127   for (c = cStart; c < cEnd; c++) {
1128     PetscReal vol;
1129     PetscInt  i, closureSize = 0;
1130     PetscInt  *closure = NULL;
1131     PetscBool anyRefine  = PETSC_FALSE;
1132     PetscBool anyCoarsen = PETSC_FALSE;
1133     PetscBool anyKeep    = PETSC_FALSE;
1134 
1135     ierr = DMPlexComputeCellGeometryFVM(dm,c,&vol,NULL,NULL);CHKERRQ(ierr);
1136     maxVolumes[c - cStart] = vol;
1137     ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1138     for (i = 0; i < closureSize; i++) {
1139       PetscInt point = closure[2 * i], refFlag;
1140 
1141       ierr = DMLabelGetValue(adaptLabel,point,&refFlag);CHKERRQ(ierr);
1142       switch (refFlag) {
1143       case DM_ADAPT_REFINE:
1144         anyRefine = PETSC_TRUE;
1145         break;
1146       case DM_ADAPT_COARSEN:
1147         anyCoarsen = PETSC_TRUE;
1148         break;
1149       case DM_ADAPT_KEEP:
1150         anyKeep = PETSC_TRUE;
1151         break;
1152       case DM_ADAPT_DETERMINE:
1153         break;
1154       default:
1155         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n",refFlag);
1156         break;
1157       }
1158       if (anyRefine) break;
1159     }
1160     ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1161     if (anyRefine) {
1162       maxVolumes[c - cStart] = vol / refRatio;
1163     } else if (anyKeep) {
1164       maxVolumes[c - cStart] = vol;
1165     } else if (anyCoarsen) {
1166       maxVolumes[c - cStart] = vol * refRatio;
1167     }
1168   }
1169   PetscFunctionReturn(0);
1170 }
1171 #endif
1172 
1173 static PetscErrorCode DMRefine_Plex_Internal(DM dm, MPI_Comm comm, DMLabel adaptLabel, DM *dmRefined)
1174 {
1175   PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *);
1176   PetscReal        refinementLimit;
1177   PetscInt         dim, cStart, cEnd;
1178   char             genname[1024], *name = NULL;
1179   PetscBool        isUniform, isTriangle = PETSC_FALSE, isTetgen = PETSC_FALSE, isCTetgen = PETSC_FALSE, flg, localized;
1180   PetscErrorCode   ierr;
1181 
1182   PetscFunctionBegin;
1183   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1184   ierr = DMPlexGetRefinementUniform(dm, &isUniform);CHKERRQ(ierr);
1185   if (isUniform) {
1186     CellRefiner cellRefiner;
1187 
1188     ierr = DMPlexGetCellRefiner_Internal(dm, &cellRefiner);CHKERRQ(ierr);
1189     ierr = DMPlexRefineUniform_Internal(dm, cellRefiner, dmRefined);CHKERRQ(ierr);
1190     ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr);
1191     if (localized) {
1192       ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);
1193     }
1194     PetscFunctionReturn(0);
1195   }
1196   ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr);
1197   ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr);
1198   if (refinementLimit == 0.0 && !refinementFunc) PetscFunctionReturn(0);
1199   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1200   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1201   ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr);
1202   if (flg) name = genname;
1203   if (name) {
1204     ierr = PetscStrcmp(name, "triangle", &isTriangle);CHKERRQ(ierr);
1205     ierr = PetscStrcmp(name, "tetgen",   &isTetgen);CHKERRQ(ierr);
1206     ierr = PetscStrcmp(name, "ctetgen",  &isCTetgen);CHKERRQ(ierr);
1207   }
1208   switch (dim) {
1209   case 2:
1210     if (!name || isTriangle) {
1211 #if defined(PETSC_HAVE_TRIANGLE)
1212       PetscReal *maxVolumes;
1213       PetscInt  c;
1214 
1215       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
1216       if (adaptLabel) {
1217         ierr = DMRefine_Plex_Label(dm,adaptLabel,cStart,cEnd,maxVolumes);CHKERRQ(ierr);
1218       } else if (refinementFunc) {
1219         for (c = cStart; c < cEnd; ++c) {
1220           PetscReal vol, centroid[3];
1221           PetscReal maxVol;
1222 
1223           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
1224           ierr = (*refinementFunc)(centroid, &maxVol);CHKERRQ(ierr);
1225           maxVolumes[c - cStart] = (double) maxVol;
1226         }
1227       } else {
1228         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
1229       }
1230 #if !defined(PETSC_USE_REAL_DOUBLE)
1231       {
1232         double *mvols;
1233         ierr = PetscMalloc1(cEnd - cStart,&mvols);CHKERRQ(ierr);
1234         for (c = 0; c < cEnd-cStart; ++c) mvols[c] = (double)maxVolumes[c];
1235         ierr = DMPlexRefine_Triangle(dm, mvols, dmRefined);CHKERRQ(ierr);
1236         ierr = PetscFree(mvols);CHKERRQ(ierr);
1237       }
1238 #else
1239       ierr = DMPlexRefine_Triangle(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
1240 #endif
1241       ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
1242 #else
1243       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement needs external package support.\nPlease reconfigure with --download-triangle.");
1244 #endif
1245     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 2D mesh generation package %s", name);
1246     break;
1247   case 3:
1248     if (!name || isCTetgen || isTetgen) {
1249       PetscReal *maxVolumes;
1250       PetscInt   c;
1251 
1252       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
1253       if (adaptLabel) {
1254         ierr = DMRefine_Plex_Label(dm,adaptLabel,cStart,cEnd,maxVolumes);CHKERRQ(ierr);
1255       } else if (refinementFunc) {
1256         for (c = cStart; c < cEnd; ++c) {
1257           PetscReal vol, centroid[3];
1258 
1259           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
1260           ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr);
1261         }
1262       } else {
1263         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
1264       }
1265       if (!name) {
1266 #if defined(PETSC_HAVE_CTETGEN)
1267         ierr = DMPlexRefine_CTetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
1268 #elif defined(PETSC_HAVE_TETGEN)
1269         ierr = DMPlexRefine_Tetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
1270 #else
1271         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "External package CTetgen or Tetgen needed.\nPlease reconfigure with '--download-ctetgen' or '--with-clanguage=cxx --download-tetgen'.");
1272 #endif
1273       } else if (isCTetgen) {
1274 #if defined(PETSC_HAVE_CTETGEN)
1275         ierr = DMPlexRefine_CTetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
1276 #else
1277         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "CTetgen needs external package support.\nPlease reconfigure with --download-ctetgen.");
1278 #endif
1279       } else {
1280 #if defined(PETSC_HAVE_TETGEN)
1281         ierr = DMPlexRefine_Tetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
1282 #else
1283         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tetgen needs external package support.\nPlease reconfigure with --download-tetgen.");
1284 #endif
1285       }
1286       ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
1287     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 3D mesh generation package %s", name);
1288     break;
1289   default:
1290     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %d is not supported.", dim);
1291   }
1292   ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr);
1293   if (localized) {
1294     ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);
1295   }
1296   PetscFunctionReturn(0);
1297 }
1298 
1299 PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined)
1300 {
1301   PetscErrorCode ierr;
1302 
1303   PetscFunctionBegin;
1304   ierr = DMRefine_Plex_Internal(dm,comm,NULL,dmRefined);CHKERRQ(ierr);
1305   PetscFunctionReturn(0);
1306 }
1307 
1308 PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted)
1309 {
1310   PetscInt       defFlag, minFlag, maxFlag, numFlags, i;
1311   const PetscInt *flags;
1312   IS             flagIS;
1313   PetscErrorCode ierr;
1314 
1315   PetscFunctionBegin;
1316   ierr = DMLabelGetDefaultValue(adaptLabel,&defFlag);CHKERRQ(ierr);
1317   minFlag = defFlag;
1318   maxFlag = defFlag;
1319   ierr = DMLabelGetValueIS(adaptLabel,&flagIS);CHKERRQ(ierr);
1320   ierr = ISGetLocalSize(flagIS,&numFlags);CHKERRQ(ierr);
1321   ierr = ISGetIndices(flagIS,&flags);CHKERRQ(ierr);
1322   for (i = 0; i < numFlags; i++) {
1323     PetscInt flag = flags[i];
1324 
1325     minFlag = PetscMin(minFlag,flag);
1326     maxFlag = PetscMax(maxFlag,flag);
1327   }
1328   ierr = ISRestoreIndices(flagIS,&flags);CHKERRQ(ierr);
1329   ierr = ISDestroy(&flagIS);CHKERRQ(ierr);
1330   {
1331     PetscInt minMaxFlag[2], minMaxFlagGlobal[2];
1332 
1333     minMaxFlag[0] = minFlag;
1334     minMaxFlag[1] = -maxFlag;
1335     ierr = MPI_Allreduce(minMaxFlag,minMaxFlagGlobal,2,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
1336     minFlag = minMaxFlagGlobal[0];
1337     maxFlag = -minMaxFlagGlobal[1];
1338   }
1339   if (minFlag == maxFlag) {
1340     switch (minFlag) {
1341     case DM_ADAPT_DETERMINE:
1342       *dmAdapted = NULL;
1343       break;
1344     case DM_ADAPT_REFINE:
1345       ierr = DMPlexSetRefinementUniform(dm,PETSC_TRUE);CHKERRQ(ierr);
1346       ierr = DMRefine(dm,MPI_COMM_NULL,dmAdapted);CHKERRQ(ierr);
1347       break;
1348     case DM_ADAPT_COARSEN:
1349       ierr = DMCoarsen(dm,MPI_COMM_NULL,dmAdapted);CHKERRQ(ierr);
1350       break;
1351     default:
1352       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n",minFlag);
1353       break;
1354     }
1355   } else {
1356     ierr = DMPlexSetRefinementUniform(dm,PETSC_FALSE);CHKERRQ(ierr);
1357     ierr = DMRefine_Plex_Internal(dm,MPI_COMM_NULL,adaptLabel,dmAdapted);CHKERRQ(ierr);
1358   }
1359   PetscFunctionReturn(0);
1360 }
1361 
1362 PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[])
1363 {
1364   DM             cdm = dm;
1365   PetscInt       r;
1366   PetscBool      isUniform, localized;
1367   PetscErrorCode ierr;
1368 
1369   PetscFunctionBegin;
1370   ierr = DMPlexGetRefinementUniform(dm, &isUniform);CHKERRQ(ierr);
1371   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1372   if (!isUniform) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Non-uniform refinement is incompatible with the hierarchy");
1373   for (r = 0; r < nlevels; ++r) {
1374     CellRefiner cellRefiner;
1375 
1376     ierr = DMPlexGetCellRefiner_Internal(cdm, &cellRefiner);CHKERRQ(ierr);
1377     ierr = DMPlexRefineUniform_Internal(cdm, cellRefiner, &dmRefined[r]);CHKERRQ(ierr);
1378     ierr = DMCopyBoundary(cdm, dmRefined[r]);CHKERRQ(ierr);
1379     if (localized) {
1380       ierr = DMLocalizeCoordinates(dmRefined[r]);CHKERRQ(ierr);
1381     }
1382     ierr = DMSetCoarseDM(dmRefined[r], cdm);CHKERRQ(ierr);
1383     ierr = DMPlexSetRegularRefinement(dmRefined[r], PETSC_TRUE);CHKERRQ(ierr);
1384     cdm  = dmRefined[r];
1385   }
1386   PetscFunctionReturn(0);
1387 }
1388