xref: /petsc/src/dm/impls/plex/generators/tetgen/tetgenerate.cxx (revision 6fcee9a8fefc5a8db0df9ee45d6f015910f1fa48)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 
3 #ifdef PETSC_HAVE_EGADS
4 #include <egads.h>
5 /* Need to make EGADSLite header compatible */
6 extern "C" int EGlite_getTopology(const ego, ego *, int *, int *, double *, int *, ego **, int **);
7 extern "C" int EGlite_inTopology(const ego, const double *);
8 #endif
9 
10 #if defined(PETSC_HAVE_TETGEN_TETLIBRARY_NEEDED)
11 #define TETLIBRARY
12 #endif
13 #include <tetgen.h>
14 
15 /* This is to fix the tetrahedron orientation from TetGen */
16 static PetscErrorCode DMPlexInvertCells_Tetgen(PetscInt numCells, PetscInt numCorners, PetscInt cells[])
17 {
18   PetscInt bound = numCells*numCorners, coff;
19 
20   PetscFunctionBegin;
21 #define SWAP(a,b) do { PetscInt tmp = (a); (a) = (b); (b) = tmp; } while (0)
22   for (coff = 0; coff < bound; coff += numCorners) SWAP(cells[coff],cells[coff+1]);
23 #undef SWAP
24   PetscFunctionReturn(0);
25 }
26 
27 PETSC_EXTERN PetscErrorCode DMPlexGenerate_Tetgen(DM boundary, PetscBool interpolate, DM *dm)
28 {
29   MPI_Comm               comm;
30   const PetscInt         dim = 3;
31   ::tetgenio             in;
32   ::tetgenio             out;
33   PetscContainer         modelObj;
34   DMUniversalLabel       universal;
35   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, defVal;
36   DMPlexInterpolatedFlag isInterpolated;
37   PetscMPIInt            rank;
38 
39   PetscFunctionBegin;
40   PetscCall(PetscObjectGetComm((PetscObject)boundary,&comm));
41   PetscCallMPI(MPI_Comm_rank(comm, &rank));
42   PetscCall(DMPlexIsInterpolatedCollective(boundary, &isInterpolated));
43   PetscCall(DMUniversalLabelCreate(boundary, &universal));
44   PetscCall(DMLabelGetDefaultValue(universal->label, &defVal));
45 
46   PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd));
47   in.numberofpoints = vEnd - vStart;
48   if (in.numberofpoints > 0) {
49     PetscSection       coordSection;
50     Vec                coordinates;
51     const PetscScalar *array;
52 
53     in.pointlist       = new double[in.numberofpoints*dim];
54     in.pointmarkerlist = new int[in.numberofpoints];
55 
56     PetscCall(PetscArrayzero(in.pointmarkerlist, (size_t) in.numberofpoints));
57     PetscCall(DMGetCoordinatesLocal(boundary, &coordinates));
58     PetscCall(DMGetCoordinateSection(boundary, &coordSection));
59     PetscCall(VecGetArrayRead(coordinates, &array));
60     for (v = vStart; v < vEnd; ++v) {
61       const PetscInt idx = v - vStart;
62       PetscInt       off, d, val;
63 
64       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
65       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
66       PetscCall(DMLabelGetValue(universal->label, v, &val));
67       if (val != defVal) in.pointmarkerlist[idx] = (int) val;
68     }
69     PetscCall(VecRestoreArrayRead(coordinates, &array));
70   }
71 
72   PetscCall(DMPlexGetHeightStratum(boundary, 1, &eStart, &eEnd));
73   in.numberofedges = eEnd - eStart;
74   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberofedges > 0) {
75     in.edgelist       = new int[in.numberofedges * 2];
76     in.edgemarkerlist = new int[in.numberofedges];
77     for (e = eStart; e < eEnd; ++e) {
78       const PetscInt  idx = e - eStart;
79       const PetscInt *cone;
80       PetscInt        coneSize, val;
81 
82       PetscCall(DMPlexGetConeSize(boundary, e, &coneSize));
83       PetscCall(DMPlexGetCone(boundary, e, &cone));
84       in.edgelist[idx*2]     = cone[0] - vStart;
85       in.edgelist[idx*2 + 1] = cone[1] - vStart;
86 
87       PetscCall(DMLabelGetValue(universal->label, e, &val));
88       if (val != defVal) in.edgemarkerlist[idx] = (int) val;
89     }
90   }
91 
92   PetscCall(DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd));
93   in.numberoffacets = fEnd - fStart;
94   if (in.numberoffacets > 0) {
95     in.facetlist       = new tetgenio::facet[in.numberoffacets];
96     in.facetmarkerlist = new int[in.numberoffacets];
97     for (f = fStart; f < fEnd; ++f) {
98       const PetscInt idx    = f - fStart;
99       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, val = -1;
100 
101       in.facetlist[idx].numberofpolygons = 1;
102       in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
103       in.facetlist[idx].numberofholes    = 0;
104       in.facetlist[idx].holelist         = NULL;
105 
106       PetscCall(DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points));
107       for (p = 0; p < numPoints*2; p += 2) {
108         const PetscInt point = points[p];
109         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
110       }
111 
112       tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
113       poly->numberofvertices = numVertices;
114       poly->vertexlist       = new int[poly->numberofvertices];
115       for (v = 0; v < numVertices; ++v) {
116         const PetscInt vIdx = points[v] - vStart;
117         poly->vertexlist[v] = vIdx;
118       }
119       PetscCall(DMLabelGetValue(universal->label, f, &val));
120       if (val != defVal) in.facetmarkerlist[idx] = (int) val;
121       PetscCall(DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points));
122     }
123   }
124   if (rank == 0) {
125     DM_Plex *mesh = (DM_Plex *) boundary->data;
126     char     args[32];
127 
128     /* Take away 'Q' for verbose output */
129 #ifdef PETSC_HAVE_EGADS
130     PetscCall(PetscStrcpy(args, "pqezQY"));
131 #else
132     PetscCall(PetscStrcpy(args, "pqezQ"));
133 #endif
134     if (mesh->tetgenOpts) {::tetrahedralize(mesh->tetgenOpts, &in, &out);}
135     else                  {::tetrahedralize(args, &in, &out);}
136   }
137   {
138     const PetscInt   numCorners  = 4;
139     const PetscInt   numCells    = out.numberoftetrahedra;
140     const PetscInt   numVertices = out.numberofpoints;
141     PetscReal        *meshCoords = NULL;
142     PetscInt         *cells      = NULL;
143 
144     if (sizeof (PetscReal) == sizeof (out.pointlist[0])) {
145       meshCoords = (PetscReal *) out.pointlist;
146     } else {
147       PetscInt i;
148 
149       meshCoords = new PetscReal[dim * numVertices];
150       for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out.pointlist[i];
151     }
152     if (sizeof (PetscInt) == sizeof (out.tetrahedronlist[0])) {
153       cells = (PetscInt *) out.tetrahedronlist;
154     } else {
155       PetscInt i;
156 
157       cells = new PetscInt[numCells * numCorners];
158       for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt) out.tetrahedronlist[i];
159     }
160 
161     PetscCall(DMPlexInvertCells_Tetgen(numCells, numCorners, cells));
162     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm));
163 
164     /* Set labels */
165     PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dm));
166     for (v = 0; v < numVertices; ++v) {
167       if (out.pointmarkerlist[v]) {
168         PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, v+numCells, out.pointmarkerlist[v]));
169       }
170     }
171     if (interpolate) {
172       PetscInt e;
173 
174       for (e = 0; e < out.numberofedges; e++) {
175         if (out.edgemarkerlist[e]) {
176           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
177           const PetscInt *edges;
178           PetscInt        numEdges;
179 
180           PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges));
181           PetscCheck(numEdges == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
182           PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, edges[0], out.edgemarkerlist[e]));
183           PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges));
184         }
185       }
186       for (f = 0; f < out.numberoftrifaces; f++) {
187         if (out.trifacemarkerlist[f]) {
188           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
189           const PetscInt *faces;
190           PetscInt        numFaces;
191 
192           PetscCall(DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces));
193           PetscCheck(numFaces == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces);
194           PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, faces[0], out.trifacemarkerlist[f]));
195           PetscCall(DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces));
196         }
197       }
198     }
199 
200     PetscCall(PetscObjectQuery((PetscObject) boundary, "EGADS Model", (PetscObject *) &modelObj));
201     if (modelObj) {
202 #ifdef PETSC_HAVE_EGADS
203       DMLabel        bodyLabel;
204       PetscInt       cStart, cEnd, c, eStart, eEnd, fStart, fEnd;
205       PetscBool      islite = PETSC_FALSE;
206       ego           *bodies;
207       ego            model, geom;
208       int            Nb, oclass, mtype, *senses;
209 
210       /* Get Attached EGADS Model from Original DMPlex */
211       PetscCall(PetscObjectQuery((PetscObject) boundary, "EGADS Model", (PetscObject *) &modelObj));
212       if (modelObj) {
213         PetscCall(PetscContainerGetPointer(modelObj, (void **) &model));
214         PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
215         /* Transfer EGADS Model to Volumetric Mesh */
216         PetscCall(PetscObjectCompose((PetscObject) *dm, "EGADS Model", (PetscObject) modelObj));
217       } else {
218         PetscCall(PetscObjectQuery((PetscObject) boundary, "EGADSLite Model", (PetscObject *) &modelObj));
219         if (modelObj) {
220           PetscCall(PetscContainerGetPointer(modelObj, (void **) &model));
221           PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
222           /* Transfer EGADS Model to Volumetric Mesh */
223           PetscCall(PetscObjectCompose((PetscObject) *dm, "EGADSLite Model", (PetscObject) modelObj));
224           islite = PETSC_TRUE;
225         }
226       }
227       if (!modelObj) goto skip_egads;
228 
229       /* Set Cell Labels */
230       PetscCall(DMGetLabel(*dm, "EGADS Body ID", &bodyLabel));
231       PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd));
232       PetscCall(DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd));
233       PetscCall(DMPlexGetDepthStratum(*dm, 1, &eStart, &eEnd));
234 
235       for (c = cStart; c < cEnd; ++c) {
236         PetscReal centroid[3] = {0., 0., 0.};
237         PetscInt  b;
238 
239         /* Deterimine what body the cell's centroid is located in */
240         if (!interpolate) {
241           PetscSection   coordSection;
242           Vec            coordinates;
243           PetscScalar   *coords = NULL;
244           PetscInt       coordSize, s, d;
245 
246           PetscCall(DMGetCoordinatesLocal(*dm, &coordinates));
247           PetscCall(DMGetCoordinateSection(*dm, &coordSection));
248           PetscCall(DMPlexVecGetClosure(*dm, coordSection, coordinates, c, &coordSize, &coords));
249           for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d];
250           PetscCall(DMPlexVecRestoreClosure(*dm, coordSection, coordinates, c, &coordSize, &coords));
251         } else PetscCall(DMPlexComputeCellGeometryFVM(*dm, c, NULL, centroid, NULL));
252         for (b = 0; b < Nb; ++b) {
253           if (islite) {if (EGlite_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
254           else        {if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
255         }
256         if (b < Nb) {
257           PetscInt   cval = b, eVal, fVal;
258           PetscInt *closure = NULL, Ncl, cl;
259 
260           PetscCall(DMLabelSetValue(bodyLabel, c, cval));
261           PetscCall(DMPlexGetTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure));
262           for (cl = 0; cl < Ncl; cl += 2) {
263             const PetscInt p = closure[cl];
264 
265             if (p >= eStart && p < eEnd) {
266               PetscCall(DMLabelGetValue(bodyLabel, p, &eVal));
267               if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
268             }
269             if (p >= fStart && p < fEnd) {
270               PetscCall(DMLabelGetValue(bodyLabel, p, &fVal));
271               if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
272             }
273           }
274           PetscCall(DMPlexRestoreTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure));
275         }
276       }
277 skip_egads: ;
278 #endif
279     }
280     PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
281   }
282   PetscCall(DMUniversalLabelDestroy(&universal));
283   PetscFunctionReturn(0);
284 }
285 
286 PETSC_EXTERN PetscErrorCode DMPlexRefine_Tetgen(DM dm, double *maxVolumes, DM *dmRefined)
287 {
288   MPI_Comm               comm;
289   const PetscInt         dim = 3;
290   ::tetgenio             in;
291   ::tetgenio             out;
292   PetscContainer         modelObj;
293   DMUniversalLabel       universal;
294   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, cStart, cEnd, c, defVal;
295   DMPlexInterpolatedFlag isInterpolated;
296   PetscMPIInt            rank;
297 
298   PetscFunctionBegin;
299   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
300   PetscCallMPI(MPI_Comm_rank(comm, &rank));
301   PetscCall(DMPlexIsInterpolatedCollective(dm, &isInterpolated));
302   PetscCall(DMUniversalLabelCreate(dm, &universal));
303   PetscCall(DMLabelGetDefaultValue(universal->label, &defVal));
304 
305   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
306   in.numberofpoints = vEnd - vStart;
307   if (in.numberofpoints > 0) {
308     PetscSection coordSection;
309     Vec          coordinates;
310     PetscScalar *array;
311 
312     in.pointlist       = new double[in.numberofpoints*dim];
313     in.pointmarkerlist = new int[in.numberofpoints];
314 
315     PetscCall(PetscArrayzero(in.pointmarkerlist, (size_t) in.numberofpoints));
316     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
317     PetscCall(DMGetCoordinateSection(dm, &coordSection));
318     PetscCall(VecGetArray(coordinates, &array));
319     for (v = vStart; v < vEnd; ++v) {
320       const PetscInt idx = v - vStart;
321       PetscInt       off, d, val;
322 
323       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
324       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
325       PetscCall(DMLabelGetValue(universal->label, v, &val));
326       if (val != defVal) in.pointmarkerlist[idx] = (int) val;
327     }
328     PetscCall(VecRestoreArray(coordinates, &array));
329   }
330 
331   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
332   in.numberofedges = eEnd - eStart;
333   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberofedges > 0) {
334     in.edgelist       = new int[in.numberofedges * 2];
335     in.edgemarkerlist = new int[in.numberofedges];
336     for (e = eStart; e < eEnd; ++e) {
337       const PetscInt  idx = e - eStart;
338       const PetscInt *cone;
339       PetscInt        coneSize, val;
340 
341       PetscCall(DMPlexGetConeSize(dm, e, &coneSize));
342       PetscCall(DMPlexGetCone(dm, e, &cone));
343       in.edgelist[idx*2]     = cone[0] - vStart;
344       in.edgelist[idx*2 + 1] = cone[1] - vStart;
345 
346       PetscCall(DMLabelGetValue(universal->label, e, &val));
347       if (val != defVal) in.edgemarkerlist[idx] = (int) val;
348     }
349   }
350 
351   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
352   in.numberoffacets = fEnd - fStart;
353   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberoffacets > 0) {
354     in.facetlist       = new tetgenio::facet[in.numberoffacets];
355     in.facetmarkerlist = new int[in.numberoffacets];
356     for (f = fStart; f < fEnd; ++f) {
357       const PetscInt idx    = f - fStart;
358       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, val;
359 
360       in.facetlist[idx].numberofpolygons = 1;
361       in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
362       in.facetlist[idx].numberofholes    = 0;
363       in.facetlist[idx].holelist         = NULL;
364 
365       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points));
366       for (p = 0; p < numPoints*2; p += 2) {
367         const PetscInt point = points[p];
368         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
369       }
370 
371       tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
372       poly->numberofvertices = numVertices;
373       poly->vertexlist       = new int[poly->numberofvertices];
374       for (v = 0; v < numVertices; ++v) {
375         const PetscInt vIdx = points[v] - vStart;
376         poly->vertexlist[v] = vIdx;
377       }
378 
379       PetscCall(DMLabelGetValue(universal->label, f, &val));
380       if (val != defVal) in.facetmarkerlist[idx] = (int) val;
381 
382       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points));
383     }
384   }
385 
386   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
387   in.numberofcorners       = 4;
388   in.numberoftetrahedra    = cEnd - cStart;
389   in.tetrahedronvolumelist = (double *) maxVolumes;
390   if (in.numberoftetrahedra > 0) {
391     in.tetrahedronlist = new int[in.numberoftetrahedra*in.numberofcorners];
392     for (c = cStart; c < cEnd; ++c) {
393       const PetscInt idx     = c - cStart;
394       PetscInt      *closure = NULL;
395       PetscInt       closureSize;
396 
397       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
398       PetscCheck(!(closureSize != 5) || !(closureSize != 15),comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %" PetscInt_FMT " vertices in closure", closureSize);
399       for (v = 0; v < 4; ++v) in.tetrahedronlist[idx*in.numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
400       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
401     }
402   }
403 
404   if (rank == 0) {
405     char args[32];
406 
407     /* Take away 'Q' for verbose output */
408     PetscCall(PetscStrcpy(args, "qezQra"));
409     ::tetrahedralize(args, &in, &out);
410   }
411 
412   in.tetrahedronvolumelist = NULL;
413   {
414     const PetscInt   numCorners  = 4;
415     const PetscInt   numCells    = out.numberoftetrahedra;
416     const PetscInt   numVertices = out.numberofpoints;
417     PetscReal        *meshCoords = NULL;
418     PetscInt         *cells      = NULL;
419     PetscBool        interpolate = isInterpolated == DMPLEX_INTERPOLATED_FULL ? PETSC_TRUE : PETSC_FALSE;
420 
421     if (sizeof (PetscReal) == sizeof (out.pointlist[0])) {
422       meshCoords = (PetscReal *) out.pointlist;
423     } else {
424       PetscInt i;
425 
426       meshCoords = new PetscReal[dim * numVertices];
427       for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out.pointlist[i];
428     }
429     if (sizeof (PetscInt) == sizeof (out.tetrahedronlist[0])) {
430       cells = (PetscInt *) out.tetrahedronlist;
431     } else {
432       PetscInt i;
433 
434       cells = new PetscInt[numCells * numCorners];
435       for (i = 0; i < numCells * numCorners; ++i)cells[i] = (PetscInt) out.tetrahedronlist[i];
436     }
437 
438     PetscCall(DMPlexInvertCells_Tetgen(numCells, numCorners, cells));
439     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined));
440     if (sizeof (PetscReal) != sizeof (out.pointlist[0])) {delete [] meshCoords;}
441     if (sizeof (PetscInt) != sizeof (out.tetrahedronlist[0])) {delete [] cells;}
442 
443     /* Set labels */
444     PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dmRefined));
445     for (v = 0; v < numVertices; ++v) {
446       if (out.pointmarkerlist[v]) {
447         PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, v+numCells, out.pointmarkerlist[v]));
448       }
449     }
450     if (interpolate) {
451       PetscInt e, f;
452 
453       for (e = 0; e < out.numberofedges; ++e) {
454         if (out.edgemarkerlist[e]) {
455           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
456           const PetscInt *edges;
457           PetscInt        numEdges;
458 
459           PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges));
460           PetscCheck(numEdges == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
461           PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, edges[0], out.edgemarkerlist[e]));
462           PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges));
463         }
464       }
465       for (f = 0; f < out.numberoftrifaces; ++f) {
466         if (out.trifacemarkerlist[f]) {
467           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
468           const PetscInt *faces;
469           PetscInt        numFaces;
470 
471           PetscCall(DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces));
472           PetscCheck(numFaces == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces);
473           PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, faces[0], out.trifacemarkerlist[f]));
474           PetscCall(DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces));
475         }
476       }
477     }
478 
479     PetscCall(PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj));
480     if (modelObj) {
481 #ifdef PETSC_HAVE_EGADS
482       DMLabel        bodyLabel;
483       PetscInt       cStart, cEnd, c, eStart, eEnd, fStart, fEnd;
484       PetscBool      islite = PETSC_FALSE;
485       ego           *bodies;
486       ego            model, geom;
487       int            Nb, oclass, mtype, *senses;
488 
489       /* Get Attached EGADS Model from Original DMPlex */
490       PetscCall(PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj));
491       if (modelObj) {
492         PetscCall(PetscContainerGetPointer(modelObj, (void **) &model));
493         PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
494         /* Transfer EGADS Model to Volumetric Mesh */
495         PetscCall(PetscObjectCompose((PetscObject) *dmRefined, "EGADS Model", (PetscObject) modelObj));
496       } else {
497         PetscCall(PetscObjectQuery((PetscObject) dm, "EGADSLite Model", (PetscObject *) &modelObj));
498         if (modelObj) {
499           PetscCall(PetscContainerGetPointer(modelObj, (void **) &model));
500           PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
501           /* Transfer EGADS Model to Volumetric Mesh */
502           PetscCall(PetscObjectCompose((PetscObject) *dmRefined, "EGADSLite Model", (PetscObject) modelObj));
503           islite = PETSC_TRUE;
504         }
505       }
506       if (!modelObj) goto skip_egads;
507 
508       /* Set Cell Labels */
509       PetscCall(DMGetLabel(*dmRefined, "EGADS Body ID", &bodyLabel));
510       PetscCall(DMPlexGetHeightStratum(*dmRefined, 0, &cStart, &cEnd));
511       PetscCall(DMPlexGetHeightStratum(*dmRefined, 1, &fStart, &fEnd));
512       PetscCall(DMPlexGetDepthStratum(*dmRefined, 1, &eStart, &eEnd));
513 
514       for (c = cStart; c < cEnd; ++c) {
515         PetscReal centroid[3] = {0., 0., 0.};
516         PetscInt  b;
517 
518         /* Deterimine what body the cell's centroid is located in */
519         if (!interpolate) {
520           PetscSection   coordSection;
521           Vec            coordinates;
522           PetscScalar   *coords = NULL;
523           PetscInt       coordSize, s, d;
524 
525           PetscCall(DMGetCoordinatesLocal(*dmRefined, &coordinates));
526           PetscCall(DMGetCoordinateSection(*dmRefined, &coordSection));
527           PetscCall(DMPlexVecGetClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords));
528           for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d];
529           PetscCall(DMPlexVecRestoreClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords));
530         } else PetscCall(DMPlexComputeCellGeometryFVM(*dmRefined, c, NULL, centroid, NULL));
531         for (b = 0; b < Nb; ++b) {
532           if (islite) {if (EGlite_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
533           else        {if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
534         }
535         if (b < Nb) {
536           PetscInt   cval = b, eVal, fVal;
537           PetscInt *closure = NULL, Ncl, cl;
538 
539           PetscCall(DMLabelSetValue(bodyLabel, c, cval));
540           PetscCall(DMPlexGetTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure));
541           for (cl = 0; cl < Ncl; cl += 2) {
542             const PetscInt p = closure[cl];
543 
544             if (p >= eStart && p < eEnd) {
545               PetscCall(DMLabelGetValue(bodyLabel, p, &eVal));
546               if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
547             }
548             if (p >= fStart && p < fEnd) {
549               PetscCall(DMLabelGetValue(bodyLabel, p, &fVal));
550               if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
551             }
552           }
553           PetscCall(DMPlexRestoreTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure));
554         }
555       }
556 skip_egads: ;
557 #endif
558     }
559     PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE));
560   }
561   PetscFunctionReturn(0);
562 }
563