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