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