Actual source code: ctetgenerate.c

  1: #include <petsc/private/dmpleximpl.h>

  3: #ifdef PETSC_HAVE_EGADS
  4:   #include <egads.h>
  5: #endif

  7: #include <ctetgen.h>

  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;

 14:   PetscFunctionBegin;
 15: #define SWAP(a, b) \
 16:   do { \
 17:     PetscInt tmp = (a); \
 18:     (a)          = (b); \
 19:     (b)          = tmp; \
 20:   } while (0)
 21:   for (coff = 0; coff < bound; coff += numCorners) SWAP(cells[coff], cells[coff + 1]);
 22: #undef SWAP
 23:   PetscFunctionReturn(PETSC_SUCCESS);
 24: }

 26: PETSC_EXTERN PetscErrorCode DMPlexGenerate_CTetgen(DM boundary, PetscBool interpolate, DM *dm)
 27: {
 28:   MPI_Comm               comm;
 29:   const PetscInt         dim = 3;
 30:   PLC                   *in, *out;
 31:   DMUniversalLabel       universal;
 32:   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, defVal, verbose = 0;
 33:   DMPlexInterpolatedFlag isInterpolated;
 34:   PetscMPIInt            rank;

 36:   PetscFunctionBegin;
 37:   PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)boundary)->prefix, "-ctetgen_verbose", &verbose, NULL));
 38:   PetscCall(PetscObjectGetComm((PetscObject)boundary, &comm));
 39:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
 40:   PetscCall(DMPlexIsInterpolatedCollective(boundary, &isInterpolated));
 41:   PetscCall(DMUniversalLabelCreate(boundary, &universal));
 42:   PetscCall(DMLabelGetDefaultValue(universal->label, &defVal));

 44:   PetscCall(PLCCreate(&in));
 45:   PetscCall(PLCCreate(&out));

 47:   PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd));
 48:   in->numberofpoints = vEnd - vStart;
 49:   if (in->numberofpoints > 0) {
 50:     PetscSection       coordSection;
 51:     Vec                coordinates;
 52:     const PetscScalar *array;

 54:     PetscCall(PetscMalloc1(in->numberofpoints * dim, &in->pointlist));
 55:     PetscCall(PetscCalloc1(in->numberofpoints, &in->pointmarkerlist));
 56:     PetscCall(DMGetCoordinatesLocal(boundary, &coordinates));
 57:     PetscCall(DMGetCoordinateSection(boundary, &coordSection));
 58:     PetscCall(VecGetArrayRead(coordinates, &array));
 59:     for (v = vStart; v < vEnd; ++v) {
 60:       const PetscInt idx = v - vStart;
 61:       PetscInt       off, d, m;

 63:       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
 64:       for (d = 0; d < dim; ++d) in->pointlist[idx * dim + d] = PetscRealPart(array[off + d]);
 65:       PetscCall(DMLabelGetValue(universal->label, v, &m));
 66:       if (m != defVal) in->pointmarkerlist[idx] = (int)m;
 67:     }
 68:     PetscCall(VecRestoreArrayRead(coordinates, &array));
 69:   }

 71:   PetscCall(DMPlexGetHeightStratum(boundary, 1, &eStart, &eEnd));
 72:   in->numberofedges = eEnd - eStart;
 73:   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberofedges > 0) {
 74:     PetscCall(PetscMalloc1(in->numberofedges * 2, &in->edgelist));
 75:     PetscCall(PetscMalloc1(in->numberofedges, &in->edgemarkerlist));
 76:     for (e = eStart; e < eEnd; ++e) {
 77:       const PetscInt  idx = e - eStart;
 78:       const PetscInt *cone;
 79:       PetscInt        coneSize, val;

 81:       PetscCall(DMPlexGetConeSize(boundary, e, &coneSize));
 82:       PetscCall(DMPlexGetCone(boundary, e, &cone));
 83:       in->edgelist[idx * 2]     = cone[0] - vStart;
 84:       in->edgelist[idx * 2 + 1] = cone[1] - vStart;

 86:       PetscCall(DMLabelGetValue(universal->label, e, &val));
 87:       if (val != defVal) in->edgemarkerlist[idx] = (int)val;
 88:     }
 89:   }

 91:   PetscCall(DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd));
 92:   in->numberoffacets = fEnd - fStart;
 93:   if (in->numberoffacets > 0) {
 94:     PetscCall(PetscMalloc1(in->numberoffacets, &in->facetlist));
 95:     PetscCall(PetscMalloc1(in->numberoffacets, &in->facetmarkerlist));
 96:     for (f = fStart; f < fEnd; ++f) {
 97:       const PetscInt idx    = f - fStart;
 98:       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, m = -1;
 99:       polygon       *poly;

101:       in->facetlist[idx].numberofpolygons = 1;
102:       PetscCall(PetscMalloc1(in->facetlist[idx].numberofpolygons, &in->facetlist[idx].polygonlist));
103:       in->facetlist[idx].numberofholes = 0;
104:       in->facetlist[idx].holelist      = NULL;

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:       }

112:       poly                   = in->facetlist[idx].polygonlist;
113:       poly->numberofvertices = numVertices;
114:       PetscCall(PetscMalloc1(poly->numberofvertices, &poly->vertexlist));
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, &m));
120:       if (m != defVal) in->facetmarkerlist[idx] = (int)m;
121:       PetscCall(DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points));
122:     }
123:   }
124:   if (rank == 0) {
125:     TetGenOpts t;

127:     PetscCall(TetGenOptsInitialize(&t));
128:     t.in        = boundary; /* Should go away */
129:     t.plc       = 1;
130:     t.quality   = 1;
131:     t.edgesout  = 1;
132:     t.zeroindex = 1;
133:     t.quiet     = 1;
134:     t.verbose   = verbose;
135: #if 0
136:   #ifdef PETSC_HAVE_EGADS
137:     /* Need to add in more TetGen code */
138:     t.nobisect  = 1; /* Add Y to preserve Surface Mesh for EGADS */
139:   #endif
140: #endif

142:     PetscCall(TetGenCheckOpts(&t));
143:     PetscCall(TetGenTetrahedralize(&t, in, out));
144:   }
145:   {
146:     const PetscInt numCorners  = 4;
147:     const PetscInt numCells    = out->numberoftetrahedra;
148:     const PetscInt numVertices = out->numberofpoints;
149:     PetscReal     *meshCoords  = NULL;
150:     PetscInt      *cells       = NULL;

152:     if (sizeof(PetscReal) == sizeof(out->pointlist[0])) {
153:       meshCoords = (PetscReal *)out->pointlist;
154:     } else {
155:       PetscCall(PetscMalloc1(dim * numVertices, &meshCoords));
156:       for (PetscInt i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal)out->pointlist[i];
157:     }
158:     if (sizeof(PetscInt) == sizeof(out->tetrahedronlist[0])) {
159:       cells = (PetscInt *)out->tetrahedronlist;
160:     } else {
161:       PetscCall(PetscMalloc1(numCells * numCorners, &cells));
162:       for (PetscInt i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out->tetrahedronlist[i];
163:     }

165:     PetscCall(DMPlexInvertCells_CTetgen(numCells, numCorners, cells));
166:     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm));
167:     if (sizeof(PetscReal) != sizeof(out->pointlist[0])) PetscCall(PetscFree(meshCoords));
168:     if (sizeof(PetscInt) != sizeof(out->tetrahedronlist[0])) PetscCall(PetscFree(cells));

170:     /* Set labels */
171:     PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dm));
172:     for (v = 0; v < numVertices; ++v) {
173:       if (out->pointmarkerlist[v]) PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, v + numCells, out->pointmarkerlist[v]));
174:     }
175:     if (interpolate) {
176:       for (PetscInt e = 0; e < out->numberofedges; e++) {
177:         if (out->edgemarkerlist[e]) {
178:           const PetscInt  vertices[2] = {out->edgelist[e * 2 + 0] + numCells, out->edgelist[e * 2 + 1] + numCells};
179:           const PetscInt *edges;
180:           PetscInt        numEdges;

182:           PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges));
183:           PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
184:           PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, edges[0], out->edgemarkerlist[e]));
185:           PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges));
186:         }
187:       }
188:       for (f = 0; f < out->numberoftrifaces; f++) {
189:         if (out->trifacemarkerlist[f]) {
190:           const PetscInt  vertices[3] = {out->trifacelist[f * 3 + 0] + numCells, out->trifacelist[f * 3 + 1] + numCells, out->trifacelist[f * 3 + 2] + numCells};
191:           const PetscInt *faces;
192:           PetscInt        numFaces;

194:           PetscCall(DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces));
195:           PetscCheck(numFaces == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces);
196:           PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, faces[0], out->trifacemarkerlist[f]));
197:           PetscCall(DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces));
198:         }
199:       }
200:     }

202: #ifdef PETSC_HAVE_EGADS
203:     {
204:       DMLabel        bodyLabel;
205:       PetscContainer modelObj;
206:       PetscInt       cStart, cEnd, c, eStart, eEnd, fStart, fEnd;
207:       ego           *bodies;
208:       ego            model, geom;
209:       int            Nb, oclass, mtype, *senses;

211:       /* Get Attached EGADS Model from Original DMPlex */
212:       PetscCall(PetscObjectQuery((PetscObject)boundary, "EGADS Model", (PetscObject *)&modelObj));
213:       if (modelObj) {
214:         PetscCall(PetscContainerGetPointer(modelObj, &model));
215:         PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
216:         /* Transfer EGADS Model to Volumetric Mesh */
217:         PetscCall(PetscObjectCompose((PetscObject)*dm, "EGADS Model", (PetscObject)modelObj));

219:         /* Set Cell Labels */
220:         PetscCall(DMGetLabel(*dm, "EGADS Body ID", &bodyLabel));
221:         PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd));
222:         PetscCall(DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd));
223:         PetscCall(DMPlexGetDepthStratum(*dm, 1, &eStart, &eEnd));

225:         for (c = cStart; c < cEnd; ++c) {
226:           PetscReal centroid[3] = {0., 0., 0.};
227:           PetscInt  b;

229:           /* Determine what body the cell's centroid is located in */
230:           if (!interpolate) {
231:             PetscSection coordSection;
232:             Vec          coordinates;
233:             PetscScalar *coords = NULL;
234:             PetscInt     coordSize, s, d;

236:             PetscCall(DMGetCoordinatesLocal(*dm, &coordinates));
237:             PetscCall(DMGetCoordinateSection(*dm, &coordSection));
238:             PetscCall(DMPlexVecGetClosure(*dm, coordSection, coordinates, c, &coordSize, &coords));
239:             for (s = 0; s < coordSize; ++s)
240:               for (d = 0; d < dim; ++d) centroid[d] += coords[s * dim + d];
241:             PetscCall(DMPlexVecRestoreClosure(*dm, coordSection, coordinates, c, &coordSize, &coords));
242:           } else PetscCall(DMPlexComputeCellGeometryFVM(*dm, c, NULL, centroid, NULL));
243:           for (b = 0; b < Nb; ++b) {
244:             if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;
245:           }
246:           if (b < Nb) {
247:             PetscInt  cval    = b, eVal, fVal;
248:             PetscInt *closure = NULL, Ncl, cl;

250:             PetscCall(DMLabelSetValue(bodyLabel, c, cval));
251:             PetscCall(DMPlexGetTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure));
252:             for (cl = 0; cl < Ncl; ++cl) {
253:               const PetscInt p = closure[cl * 2];

255:               if (p >= eStart && p < eEnd) {
256:                 PetscCall(DMLabelGetValue(bodyLabel, p, &eVal));
257:                 if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
258:               }
259:               if (p >= fStart && p < fEnd) {
260:                 PetscCall(DMLabelGetValue(bodyLabel, p, &fVal));
261:                 if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
262:               }
263:             }
264:             PetscCall(DMPlexRestoreTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure));
265:           }
266:         }
267:       }
268:     }
269: #endif
270:     PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
271:   }

273:   PetscCall(DMUniversalLabelDestroy(&universal));
274:   PetscCall(PLCDestroy(&in));
275:   PetscCall(PLCDestroy(&out));
276:   PetscFunctionReturn(PETSC_SUCCESS);
277: }

279: PETSC_EXTERN PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined)
280: {
281:   MPI_Comm               comm;
282:   const PetscInt         dim = 3;
283:   PLC                   *in, *out;
284:   DMUniversalLabel       universal;
285:   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, cStart, cEnd, c, defVal, verbose = 0;
286:   DMPlexInterpolatedFlag isInterpolated;
287:   PetscMPIInt            rank;

289:   PetscFunctionBegin;
290:   PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)dm)->prefix, "-ctetgen_verbose", &verbose, NULL));
291:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
292:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
293:   PetscCall(DMPlexIsInterpolatedCollective(dm, &isInterpolated));
294:   PetscCall(DMUniversalLabelCreate(dm, &universal));
295:   PetscCall(DMLabelGetDefaultValue(universal->label, &defVal));

297:   PetscCall(PLCCreate(&in));
298:   PetscCall(PLCCreate(&out));

300:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
301:   in->numberofpoints = vEnd - vStart;
302:   if (in->numberofpoints > 0) {
303:     PetscSection coordSection;
304:     Vec          coordinates;
305:     PetscScalar *array;

307:     PetscCall(PetscMalloc1(in->numberofpoints * dim, &in->pointlist));
308:     PetscCall(PetscCalloc1(in->numberofpoints, &in->pointmarkerlist));
309:     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
310:     PetscCall(DMGetCoordinateSection(dm, &coordSection));
311:     PetscCall(VecGetArray(coordinates, &array));
312:     for (v = vStart; v < vEnd; ++v) {
313:       const PetscInt idx = v - vStart;
314:       PetscInt       off, d, m;

316:       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
317:       for (d = 0; d < dim; ++d) in->pointlist[idx * dim + d] = PetscRealPart(array[off + d]);
318:       PetscCall(DMLabelGetValue(universal->label, v, &m));
319:       if (m != defVal) in->pointmarkerlist[idx] = (int)m;
320:     }
321:     PetscCall(VecRestoreArray(coordinates, &array));
322:   }

324:   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
325:   in->numberofedges = eEnd - eStart;
326:   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberofedges > 0) {
327:     PetscCall(PetscMalloc1(in->numberofedges * 2, &in->edgelist));
328:     PetscCall(PetscMalloc1(in->numberofedges, &in->edgemarkerlist));
329:     for (e = eStart; e < eEnd; ++e) {
330:       const PetscInt  idx = e - eStart;
331:       const PetscInt *cone;
332:       PetscInt        coneSize, val;

334:       PetscCall(DMPlexGetConeSize(dm, e, &coneSize));
335:       PetscCall(DMPlexGetCone(dm, e, &cone));
336:       in->edgelist[idx * 2]     = cone[0] - vStart;
337:       in->edgelist[idx * 2 + 1] = cone[1] - vStart;

339:       PetscCall(DMLabelGetValue(universal->label, e, &val));
340:       if (val != defVal) in->edgemarkerlist[idx] = (int)val;
341:     }
342:   }

344:   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
345:   in->numberoftrifaces = 0;
346:   for (f = fStart; f < fEnd; ++f) {
347:     PetscInt supportSize;

349:     PetscCall(DMPlexGetSupportSize(dm, f, &supportSize));
350:     if (supportSize == 1) ++in->numberoftrifaces;
351:   }
352:   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberoftrifaces > 0) {
353:     PetscInt tf = 0;

355:     PetscCall(PetscMalloc1(in->numberoftrifaces * 3, &in->trifacelist));
356:     PetscCall(PetscMalloc1(in->numberoftrifaces, &in->trifacemarkerlist));
357:     for (f = fStart; f < fEnd; ++f) {
358:       PetscInt *points = NULL;
359:       PetscInt  supportSize, numPoints, p, Nv = 0, val;

361:       PetscCall(DMPlexGetSupportSize(dm, f, &supportSize));
362:       if (supportSize != 1) continue;
363:       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points));
364:       for (p = 0; p < numPoints * 2; p += 2) {
365:         const PetscInt point = points[p];
366:         if ((point >= vStart) && (point < vEnd)) in->trifacelist[tf * 3 + Nv++] = point - vStart;
367:       }
368:       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points));
369:       PetscCheck(Nv == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " has %" PetscInt_FMT " vertices, not 3", f, Nv);
370:       PetscCall(DMLabelGetValue(universal->label, f, &val));
371:       if (val != defVal) in->trifacemarkerlist[tf] = (int)val;
372:       ++tf;
373:     }
374:   }

376:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
377:   in->numberofcorners       = 4;
378:   in->numberoftetrahedra    = cEnd - cStart;
379:   in->tetrahedronvolumelist = maxVolumes;
380:   if (in->numberoftetrahedra > 0) {
381:     PetscCall(PetscMalloc1(in->numberoftetrahedra * in->numberofcorners, &in->tetrahedronlist));
382:     for (c = cStart; c < cEnd; ++c) {
383:       const PetscInt idx     = c - cStart;
384:       PetscInt      *closure = NULL;
385:       PetscInt       closureSize;

387:       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
388:       PetscCheck((closureSize == 5) || (closureSize == 15), comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %" PetscInt_FMT " vertices in closure", closureSize);
389:       for (v = 0; v < 4; ++v) in->tetrahedronlist[idx * in->numberofcorners + v] = closure[(v + closureSize - 4) * 2] - vStart;
390:       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
391:     }
392:   }

394:   if (rank == 0) {
395:     TetGenOpts t;

397:     PetscCall(TetGenOptsInitialize(&t));

399:     t.in        = dm; /* Should go away */
400:     t.refine    = 1;
401:     t.varvolume = 1;
402:     t.quality   = 1;
403:     t.edgesout  = 1;
404:     t.zeroindex = 1;
405:     t.quiet     = 1;
406:     t.verbose   = verbose; /* Change this */

408:     PetscCall(TetGenCheckOpts(&t));
409:     PetscCall(TetGenTetrahedralize(&t, in, out));
410:   }

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;

421:     if (sizeof(PetscReal) == sizeof(out->pointlist[0])) {
422:       meshCoords = (PetscReal *)out->pointlist;
423:     } else {
424:       PetscCall(PetscMalloc1(dim * numVertices, &meshCoords));
425:       for (PetscInt i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal)out->pointlist[i];
426:     }
427:     if (sizeof(PetscInt) == sizeof(out->tetrahedronlist[0])) {
428:       cells = (PetscInt *)out->tetrahedronlist;
429:     } else {
430:       PetscCall(PetscMalloc1(numCells * numCorners, &cells));
431:       for (PetscInt i = 0; i < numCells * numCorners; ++i) cells[i] = (PetscInt)out->tetrahedronlist[i];
432:     }

434:     PetscCall(DMPlexInvertCells_CTetgen(numCells, numCorners, cells));
435:     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined));
436:     if (sizeof(PetscReal) != sizeof(out->pointlist[0])) PetscCall(PetscFree(meshCoords));
437:     if (sizeof(PetscInt) != sizeof(out->tetrahedronlist[0])) PetscCall(PetscFree(cells));

439:     /* Set labels */
440:     PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dmRefined));
441:     for (v = 0; v < numVertices; ++v) {
442:       if (out->pointmarkerlist[v]) PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, v + numCells, out->pointmarkerlist[v]));
443:     }
444:     if (interpolate) {
445:       for (PetscInt e = 0; e < out->numberofedges; e++) {
446:         if (out->edgemarkerlist[e]) {
447:           const PetscInt  vertices[2] = {out->edgelist[e * 2 + 0] + numCells, out->edgelist[e * 2 + 1] + numCells};
448:           const PetscInt *edges;
449:           PetscInt        numEdges;

451:           PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges));
452:           PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
453:           PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, edges[0], out->edgemarkerlist[e]));
454:           PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges));
455:         }
456:       }
457:       for (PetscInt f = 0; f < out->numberoftrifaces; f++) {
458:         if (out->trifacemarkerlist[f]) {
459:           const PetscInt  vertices[3] = {out->trifacelist[f * 3 + 0] + numCells, out->trifacelist[f * 3 + 1] + numCells, out->trifacelist[f * 3 + 2] + numCells};
460:           const PetscInt *faces;
461:           PetscInt        numFaces;

463:           PetscCall(DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces));
464:           PetscCheck(numFaces == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces);
465:           PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, faces[0], out->trifacemarkerlist[f]));
466:           PetscCall(DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces));
467:         }
468:       }
469:     }

471: #ifdef PETSC_HAVE_EGADS
472:     {
473:       DMLabel        bodyLabel;
474:       PetscContainer modelObj;
475:       PetscInt       cStart, cEnd, c, eStart, eEnd, fStart, fEnd;
476:       ego           *bodies;
477:       ego            model, geom;
478:       int            Nb, oclass, mtype, *senses;

480:       /* Get Attached EGADS Model from Original DMPlex */
481:       PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
482:       if (modelObj) {
483:         PetscCall(PetscContainerGetPointer(modelObj, &model));
484:         PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
485:         /* Transfer EGADS Model to Volumetric Mesh */
486:         PetscCall(PetscObjectCompose((PetscObject)*dmRefined, "EGADS Model", (PetscObject)modelObj));

488:         /* Set Cell Labels */
489:         PetscCall(DMGetLabel(*dmRefined, "EGADS Body ID", &bodyLabel));
490:         PetscCall(DMPlexGetHeightStratum(*dmRefined, 0, &cStart, &cEnd));
491:         PetscCall(DMPlexGetHeightStratum(*dmRefined, 1, &fStart, &fEnd));
492:         PetscCall(DMPlexGetDepthStratum(*dmRefined, 1, &eStart, &eEnd));

494:         for (c = cStart; c < cEnd; ++c) {
495:           PetscReal centroid[3] = {0., 0., 0.};
496:           PetscInt  b;

498:           /* Determine what body the cell's centroid is located in */
499:           if (!interpolate) {
500:             PetscSection coordSection;
501:             Vec          coordinates;
502:             PetscScalar *coords = NULL;
503:             PetscInt     coordSize, s, d;

505:             PetscCall(DMGetCoordinatesLocal(*dmRefined, &coordinates));
506:             PetscCall(DMGetCoordinateSection(*dmRefined, &coordSection));
507:             PetscCall(DMPlexVecGetClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords));
508:             for (s = 0; s < coordSize; ++s)
509:               for (d = 0; d < dim; ++d) centroid[d] += coords[s * dim + d];
510:             PetscCall(DMPlexVecRestoreClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords));
511:           } else PetscCall(DMPlexComputeCellGeometryFVM(*dmRefined, c, NULL, centroid, NULL));
512:           for (b = 0; b < Nb; ++b) {
513:             if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;
514:           }
515:           if (b < Nb) {
516:             PetscInt  cval    = b, eVal, fVal;
517:             PetscInt *closure = NULL, Ncl, cl;

519:             PetscCall(DMLabelSetValue(bodyLabel, c, cval));
520:             PetscCall(DMPlexGetTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure));
521:             for (cl = 0; cl < Ncl; cl += 2) {
522:               const PetscInt p = closure[cl];

524:               if (p >= eStart && p < eEnd) {
525:                 PetscCall(DMLabelGetValue(bodyLabel, p, &eVal));
526:                 if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
527:               }
528:               if (p >= fStart && p < fEnd) {
529:                 PetscCall(DMLabelGetValue(bodyLabel, p, &fVal));
530:                 if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
531:               }
532:             }
533:             PetscCall(DMPlexRestoreTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure));
534:           }
535:         }
536:       }
537:     }
538: #endif
539:     PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE));
540:   }
541:   PetscCall(DMUniversalLabelDestroy(&universal));
542:   PetscCall(PLCDestroy(&in));
543:   PetscCall(PLCDestroy(&out));
544:   PetscFunctionReturn(PETSC_SUCCESS);
545: }