xref: /petsc/src/dm/impls/plex/generators/tetgen/tetgenerate.cxx (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
13a074057SBarry Smith #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
23a074057SBarry Smith 
30fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
40fdc7489SMatthew Knepley #include <egads.h>
5c1cad2e7SMatthew G. Knepley /* Need to make EGADSLite header compatible */
6c1cad2e7SMatthew G. Knepley extern "C" int EGlite_getTopology(const ego, ego *, int *, int *, double *, int *, ego **, int **);
7c1cad2e7SMatthew G. Knepley extern "C" int EGlite_inTopology(const ego, const double *);
80fdc7489SMatthew Knepley #endif
90fdc7489SMatthew Knepley 
10b83bf3d7SVaclav Hapla #if defined(PETSC_HAVE_TETGEN_TETLIBRARY_NEEDED)
11b83bf3d7SVaclav Hapla #define TETLIBRARY
12b83bf3d7SVaclav Hapla #endif
133a074057SBarry Smith #include <tetgen.h>
143a074057SBarry Smith 
153a074057SBarry Smith /* This is to fix the tetrahedron orientation from TetGen */
169371c9d4SSatish Balay static PetscErrorCode DMPlexInvertCells_Tetgen(PetscInt numCells, PetscInt numCorners, PetscInt cells[]) {
173a074057SBarry Smith   PetscInt bound = numCells * numCorners, coff;
183a074057SBarry Smith 
193a074057SBarry Smith   PetscFunctionBegin;
209371c9d4SSatish Balay #define SWAP(a, b) \
219371c9d4SSatish Balay   do { \
229371c9d4SSatish Balay     PetscInt tmp = (a); \
239371c9d4SSatish Balay     (a)          = (b); \
249371c9d4SSatish Balay     (b)          = tmp; \
259371c9d4SSatish Balay   } while (0)
2696ca5757SLisandro Dalcin   for (coff = 0; coff < bound; coff += numCorners) SWAP(cells[coff], cells[coff + 1]);
2796ca5757SLisandro Dalcin #undef SWAP
283a074057SBarry Smith   PetscFunctionReturn(0);
293a074057SBarry Smith }
303a074057SBarry Smith 
319371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexGenerate_Tetgen(DM boundary, PetscBool interpolate, DM *dm) {
323a074057SBarry Smith   MPI_Comm               comm;
333a074057SBarry Smith   const PetscInt         dim = 3;
343a074057SBarry Smith   ::tetgenio             in;
353a074057SBarry Smith   ::tetgenio             out;
369318fe57SMatthew G. Knepley   PetscContainer         modelObj;
370fdc7489SMatthew Knepley   DMUniversalLabel       universal;
38af226901SMatthew G. Knepley   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, defVal;
390fdc7489SMatthew Knepley   DMPlexInterpolatedFlag isInterpolated;
403a074057SBarry Smith   PetscMPIInt            rank;
413a074057SBarry Smith 
423a074057SBarry Smith   PetscFunctionBegin;
439566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)boundary, &comm));
449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
459566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolatedCollective(boundary, &isInterpolated));
469566063dSJacob Faibussowitsch   PetscCall(DMUniversalLabelCreate(boundary, &universal));
47af226901SMatthew G. Knepley   PetscCall(DMLabelGetDefaultValue(universal->label, &defVal));
483a074057SBarry Smith 
499566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd));
503a074057SBarry Smith   in.numberofpoints = vEnd - vStart;
513a074057SBarry Smith   if (in.numberofpoints > 0) {
523a074057SBarry Smith     PetscSection       coordSection;
533a074057SBarry Smith     Vec                coordinates;
540fdc7489SMatthew Knepley     const PetscScalar *array;
553a074057SBarry Smith 
563a074057SBarry Smith     in.pointlist       = new double[in.numberofpoints * dim];
573a074057SBarry Smith     in.pointmarkerlist = new int[in.numberofpoints];
583a074057SBarry Smith 
5941e9d8b5SMatthew G. Knepley     PetscCall(PetscArrayzero(in.pointmarkerlist, (size_t)in.numberofpoints));
609566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(boundary, &coordinates));
619566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(boundary, &coordSection));
629566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coordinates, &array));
633a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
643a074057SBarry Smith       const PetscInt idx = v - vStart;
650fdc7489SMatthew Knepley       PetscInt       off, d, val;
663a074057SBarry Smith 
679566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
683a074057SBarry Smith       for (d = 0; d < dim; ++d) in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]);
699566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, v, &val));
70af226901SMatthew G. Knepley       if (val != defVal) in.pointmarkerlist[idx] = (int)val;
713a074057SBarry Smith     }
729566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coordinates, &array));
733a074057SBarry Smith   }
743a074057SBarry Smith 
759566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(boundary, 1, &eStart, &eEnd));
760fdc7489SMatthew Knepley   in.numberofedges = eEnd - eStart;
770fdc7489SMatthew Knepley   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberofedges > 0) {
780fdc7489SMatthew Knepley     in.edgelist       = new int[in.numberofedges * 2];
790fdc7489SMatthew Knepley     in.edgemarkerlist = new int[in.numberofedges];
800fdc7489SMatthew Knepley     for (e = eStart; e < eEnd; ++e) {
810fdc7489SMatthew Knepley       const PetscInt  idx = e - eStart;
820fdc7489SMatthew Knepley       const PetscInt *cone;
830fdc7489SMatthew Knepley       PetscInt        coneSize, val;
840fdc7489SMatthew Knepley 
859566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(boundary, e, &coneSize));
869566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(boundary, e, &cone));
870fdc7489SMatthew Knepley       in.edgelist[idx * 2]     = cone[0] - vStart;
880fdc7489SMatthew Knepley       in.edgelist[idx * 2 + 1] = cone[1] - vStart;
890fdc7489SMatthew Knepley 
909566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, e, &val));
91af226901SMatthew G. Knepley       if (val != defVal) in.edgemarkerlist[idx] = (int)val;
920fdc7489SMatthew Knepley     }
930fdc7489SMatthew Knepley   }
940fdc7489SMatthew Knepley 
959566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd));
963a074057SBarry Smith   in.numberoffacets = fEnd - fStart;
973a074057SBarry Smith   if (in.numberoffacets > 0) {
983a074057SBarry Smith     in.facetlist       = new tetgenio::facet[in.numberoffacets];
993a074057SBarry Smith     in.facetmarkerlist = new int[in.numberoffacets];
1003a074057SBarry Smith     for (f = fStart; f < fEnd; ++f) {
1013a074057SBarry Smith       const PetscInt idx    = f - fStart;
1020fdc7489SMatthew Knepley       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, val = -1;
1033a074057SBarry Smith 
1043a074057SBarry Smith       in.facetlist[idx].numberofpolygons = 1;
1053a074057SBarry Smith       in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
1063a074057SBarry Smith       in.facetlist[idx].numberofholes    = 0;
1073a074057SBarry Smith       in.facetlist[idx].holelist         = NULL;
1083a074057SBarry Smith 
1099566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points));
1103a074057SBarry Smith       for (p = 0; p < numPoints * 2; p += 2) {
1113a074057SBarry Smith         const PetscInt point = points[p];
1123a074057SBarry Smith         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
1133a074057SBarry Smith       }
1143a074057SBarry Smith 
1153a074057SBarry Smith       tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
1163a074057SBarry Smith       poly->numberofvertices  = numVertices;
1173a074057SBarry Smith       poly->vertexlist        = new int[poly->numberofvertices];
1183a074057SBarry Smith       for (v = 0; v < numVertices; ++v) {
1193a074057SBarry Smith         const PetscInt vIdx = points[v] - vStart;
1203a074057SBarry Smith         poly->vertexlist[v] = vIdx;
1213a074057SBarry Smith       }
1229566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, f, &val));
123af226901SMatthew G. Knepley       if (val != defVal) in.facetmarkerlist[idx] = (int)val;
1249566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points));
1253a074057SBarry Smith     }
1263a074057SBarry Smith   }
127dd400576SPatrick Sanan   if (rank == 0) {
1280fdc7489SMatthew Knepley     DM_Plex *mesh = (DM_Plex *)boundary->data;
1293a074057SBarry Smith     char     args[32];
1303a074057SBarry Smith 
1313a074057SBarry Smith     /* Take away 'Q' for verbose output */
1320fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
1339566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(args, "pqezQY"));
1340fdc7489SMatthew Knepley #else
1359566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(args, "pqezQ"));
1360fdc7489SMatthew Knepley #endif
1379371c9d4SSatish Balay     if (mesh->tetgenOpts) {
1389371c9d4SSatish Balay       ::tetrahedralize(mesh->tetgenOpts, &in, &out);
1399371c9d4SSatish Balay     } else {
1409371c9d4SSatish Balay       ::tetrahedralize(args, &in, &out);
1419371c9d4SSatish Balay     }
1423a074057SBarry Smith   }
1433a074057SBarry Smith   {
1443a074057SBarry Smith     const PetscInt numCorners  = 4;
1453a074057SBarry Smith     const PetscInt numCells    = out.numberoftetrahedra;
1463a074057SBarry Smith     const PetscInt numVertices = out.numberofpoints;
147a4a685f2SJacob Faibussowitsch     PetscReal     *meshCoords  = NULL;
148a4a685f2SJacob Faibussowitsch     PetscInt      *cells       = NULL;
149a4a685f2SJacob Faibussowitsch 
150a4a685f2SJacob Faibussowitsch     if (sizeof(PetscReal) == sizeof(out.pointlist[0])) {
151a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *)out.pointlist;
152a4a685f2SJacob Faibussowitsch     } else {
153a4a685f2SJacob Faibussowitsch       PetscInt i;
154a4a685f2SJacob Faibussowitsch 
155a4a685f2SJacob Faibussowitsch       meshCoords = new PetscReal[dim * numVertices];
1560fdc7489SMatthew Knepley       for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal)out.pointlist[i];
157a4a685f2SJacob Faibussowitsch     }
158a4a685f2SJacob Faibussowitsch     if (sizeof(PetscInt) == sizeof(out.tetrahedronlist[0])) {
159a4a685f2SJacob Faibussowitsch       cells = (PetscInt *)out.tetrahedronlist;
160a4a685f2SJacob Faibussowitsch     } else {
161a4a685f2SJacob Faibussowitsch       PetscInt i;
162a4a685f2SJacob Faibussowitsch 
163a4a685f2SJacob Faibussowitsch       cells = new PetscInt[numCells * numCorners];
1640fdc7489SMatthew Knepley       for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out.tetrahedronlist[i];
165a4a685f2SJacob Faibussowitsch     }
1663a074057SBarry Smith 
1679566063dSJacob Faibussowitsch     PetscCall(DMPlexInvertCells_Tetgen(numCells, numCorners, cells));
1689566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm));
1690fdc7489SMatthew Knepley 
1703a074057SBarry Smith     /* Set labels */
1719566063dSJacob Faibussowitsch     PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dm));
1723a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
173*48a46eb9SPierre Jolivet       if (out.pointmarkerlist[v]) PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, v + numCells, out.pointmarkerlist[v]));
1743a074057SBarry Smith     }
1753a074057SBarry Smith     if (interpolate) {
1763a074057SBarry Smith       PetscInt e;
1773a074057SBarry Smith 
1783a074057SBarry Smith       for (e = 0; e < out.numberofedges; e++) {
1793a074057SBarry Smith         if (out.edgemarkerlist[e]) {
1803a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells};
1813a074057SBarry Smith           const PetscInt *edges;
1823a074057SBarry Smith           PetscInt        numEdges;
1833a074057SBarry Smith 
1849566063dSJacob Faibussowitsch           PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges));
18563a3b9bcSJacob Faibussowitsch           PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
1869566063dSJacob Faibussowitsch           PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, edges[0], out.edgemarkerlist[e]));
1879566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges));
1883a074057SBarry Smith         }
1893a074057SBarry Smith       }
1903a074057SBarry Smith       for (f = 0; f < out.numberoftrifaces; f++) {
1913a074057SBarry Smith         if (out.trifacemarkerlist[f]) {
1923a074057SBarry Smith           const PetscInt  vertices[3] = {out.trifacelist[f * 3 + 0] + numCells, out.trifacelist[f * 3 + 1] + numCells, out.trifacelist[f * 3 + 2] + numCells};
1933a074057SBarry Smith           const PetscInt *faces;
1943a074057SBarry Smith           PetscInt        numFaces;
1953a074057SBarry Smith 
1969566063dSJacob Faibussowitsch           PetscCall(DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces));
19763a3b9bcSJacob Faibussowitsch           PetscCheck(numFaces == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces);
1989566063dSJacob Faibussowitsch           PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, faces[0], out.trifacemarkerlist[f]));
1999566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces));
2003a074057SBarry Smith         }
2013a074057SBarry Smith       }
2023a074057SBarry Smith     }
2030fdc7489SMatthew Knepley 
2049566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)boundary, "EGADS Model", (PetscObject *)&modelObj));
2059318fe57SMatthew G. Knepley     if (modelObj) {
2060fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
2070fdc7489SMatthew Knepley       DMLabel   bodyLabel;
2080fdc7489SMatthew Knepley       PetscInt  cStart, cEnd, c, eStart, eEnd, fStart, fEnd;
209c1cad2e7SMatthew G. Knepley       PetscBool islite = PETSC_FALSE;
2100fdc7489SMatthew Knepley       ego      *bodies;
2110fdc7489SMatthew Knepley       ego       model, geom;
2120fdc7489SMatthew Knepley       int       Nb, oclass, mtype, *senses;
2130fdc7489SMatthew Knepley 
2140fdc7489SMatthew Knepley       /* Get Attached EGADS Model from Original DMPlex */
2159566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery((PetscObject)boundary, "EGADS Model", (PetscObject *)&modelObj));
216c1cad2e7SMatthew G. Knepley       if (modelObj) {
2179566063dSJacob Faibussowitsch         PetscCall(PetscContainerGetPointer(modelObj, (void **)&model));
2189566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
2190fdc7489SMatthew Knepley         /* Transfer EGADS Model to Volumetric Mesh */
2209566063dSJacob Faibussowitsch         PetscCall(PetscObjectCompose((PetscObject)*dm, "EGADS Model", (PetscObject)modelObj));
221c1cad2e7SMatthew G. Knepley       } else {
2229566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)boundary, "EGADSLite Model", (PetscObject *)&modelObj));
223c1cad2e7SMatthew G. Knepley         if (modelObj) {
2249566063dSJacob Faibussowitsch           PetscCall(PetscContainerGetPointer(modelObj, (void **)&model));
2259566063dSJacob Faibussowitsch           PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
226c1cad2e7SMatthew G. Knepley           /* Transfer EGADS Model to Volumetric Mesh */
2279566063dSJacob Faibussowitsch           PetscCall(PetscObjectCompose((PetscObject)*dm, "EGADSLite Model", (PetscObject)modelObj));
228c1cad2e7SMatthew G. Knepley           islite = PETSC_TRUE;
229c1cad2e7SMatthew G. Knepley         }
230c1cad2e7SMatthew G. Knepley       }
231c1cad2e7SMatthew G. Knepley       if (!modelObj) goto skip_egads;
2320fdc7489SMatthew Knepley 
2330fdc7489SMatthew Knepley       /* Set Cell Labels */
2349566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dm, "EGADS Body ID", &bodyLabel));
2359566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd));
2369566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd));
2379566063dSJacob Faibussowitsch       PetscCall(DMPlexGetDepthStratum(*dm, 1, &eStart, &eEnd));
2380fdc7489SMatthew Knepley 
2390fdc7489SMatthew Knepley       for (c = cStart; c < cEnd; ++c) {
2400fdc7489SMatthew Knepley         PetscReal centroid[3] = {0., 0., 0.};
2410fdc7489SMatthew Knepley         PetscInt  b;
2420fdc7489SMatthew Knepley 
2430fdc7489SMatthew Knepley         /* Deterimine what body the cell's centroid is located in */
2440fdc7489SMatthew Knepley         if (!interpolate) {
2450fdc7489SMatthew Knepley           PetscSection coordSection;
2460fdc7489SMatthew Knepley           Vec          coordinates;
2470fdc7489SMatthew Knepley           PetscScalar *coords = NULL;
2480fdc7489SMatthew Knepley           PetscInt     coordSize, s, d;
2490fdc7489SMatthew Knepley 
2509566063dSJacob Faibussowitsch           PetscCall(DMGetCoordinatesLocal(*dm, &coordinates));
2519566063dSJacob Faibussowitsch           PetscCall(DMGetCoordinateSection(*dm, &coordSection));
2529566063dSJacob Faibussowitsch           PetscCall(DMPlexVecGetClosure(*dm, coordSection, coordinates, c, &coordSize, &coords));
2539371c9d4SSatish Balay           for (s = 0; s < coordSize; ++s)
2549371c9d4SSatish Balay             for (d = 0; d < dim; ++d) centroid[d] += coords[s * dim + d];
2559566063dSJacob Faibussowitsch           PetscCall(DMPlexVecRestoreClosure(*dm, coordSection, coordinates, c, &coordSize, &coords));
2561baa6e33SBarry Smith         } else PetscCall(DMPlexComputeCellGeometryFVM(*dm, c, NULL, centroid, NULL));
2570fdc7489SMatthew Knepley         for (b = 0; b < Nb; ++b) {
2589371c9d4SSatish Balay           if (islite) {
2599371c9d4SSatish Balay             if (EGlite_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;
2609371c9d4SSatish Balay           } else {
2619371c9d4SSatish Balay             if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;
2629371c9d4SSatish Balay           }
2630fdc7489SMatthew Knepley         }
2640fdc7489SMatthew Knepley         if (b < Nb) {
2650fdc7489SMatthew Knepley           PetscInt  cval    = b, eVal, fVal;
2660fdc7489SMatthew Knepley           PetscInt *closure = NULL, Ncl, cl;
2670fdc7489SMatthew Knepley 
2689566063dSJacob Faibussowitsch           PetscCall(DMLabelSetValue(bodyLabel, c, cval));
2699566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure));
2700fdc7489SMatthew Knepley           for (cl = 0; cl < Ncl; cl += 2) {
2710fdc7489SMatthew Knepley             const PetscInt p = closure[cl];
2720fdc7489SMatthew Knepley 
2730fdc7489SMatthew Knepley             if (p >= eStart && p < eEnd) {
2749566063dSJacob Faibussowitsch               PetscCall(DMLabelGetValue(bodyLabel, p, &eVal));
2759566063dSJacob Faibussowitsch               if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
2760fdc7489SMatthew Knepley             }
2770fdc7489SMatthew Knepley             if (p >= fStart && p < fEnd) {
2789566063dSJacob Faibussowitsch               PetscCall(DMLabelGetValue(bodyLabel, p, &fVal));
2799566063dSJacob Faibussowitsch               if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
2800fdc7489SMatthew Knepley             }
2810fdc7489SMatthew Knepley           }
2829566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure));
2830fdc7489SMatthew Knepley         }
2840fdc7489SMatthew Knepley       }
285c1cad2e7SMatthew G. Knepley     skip_egads:;
2860fdc7489SMatthew Knepley #endif
2879318fe57SMatthew G. Knepley     }
2889566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
2893a074057SBarry Smith   }
2909566063dSJacob Faibussowitsch   PetscCall(DMUniversalLabelDestroy(&universal));
2913a074057SBarry Smith   PetscFunctionReturn(0);
2923a074057SBarry Smith }
2933a074057SBarry Smith 
2949371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexRefine_Tetgen(DM dm, double *maxVolumes, DM *dmRefined) {
2953a074057SBarry Smith   MPI_Comm               comm;
2963a074057SBarry Smith   const PetscInt         dim = 3;
2973a074057SBarry Smith   ::tetgenio             in;
2983a074057SBarry Smith   ::tetgenio             out;
2999318fe57SMatthew G. Knepley   PetscContainer         modelObj;
3000fdc7489SMatthew Knepley   DMUniversalLabel       universal;
301af226901SMatthew G. Knepley   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, cStart, cEnd, c, defVal;
3020fdc7489SMatthew Knepley   DMPlexInterpolatedFlag isInterpolated;
3033a074057SBarry Smith   PetscMPIInt            rank;
3043a074057SBarry Smith 
3053a074057SBarry Smith   PetscFunctionBegin;
3069566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3089566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolatedCollective(dm, &isInterpolated));
3099566063dSJacob Faibussowitsch   PetscCall(DMUniversalLabelCreate(dm, &universal));
310af226901SMatthew G. Knepley   PetscCall(DMLabelGetDefaultValue(universal->label, &defVal));
3113a074057SBarry Smith 
3129566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
3133a074057SBarry Smith   in.numberofpoints = vEnd - vStart;
3143a074057SBarry Smith   if (in.numberofpoints > 0) {
3153a074057SBarry Smith     PetscSection coordSection;
3163a074057SBarry Smith     Vec          coordinates;
3173a074057SBarry Smith     PetscScalar *array;
3183a074057SBarry Smith 
3193a074057SBarry Smith     in.pointlist       = new double[in.numberofpoints * dim];
3203a074057SBarry Smith     in.pointmarkerlist = new int[in.numberofpoints];
3213a074057SBarry Smith 
32241e9d8b5SMatthew G. Knepley     PetscCall(PetscArrayzero(in.pointmarkerlist, (size_t)in.numberofpoints));
3239566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
3249566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dm, &coordSection));
3259566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &array));
3263a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
3273a074057SBarry Smith       const PetscInt idx = v - vStart;
3280fdc7489SMatthew Knepley       PetscInt       off, d, val;
3293a074057SBarry Smith 
3309566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
3313a074057SBarry Smith       for (d = 0; d < dim; ++d) in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]);
3329566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, v, &val));
333af226901SMatthew G. Knepley       if (val != defVal) in.pointmarkerlist[idx] = (int)val;
3343a074057SBarry Smith     }
3359566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &array));
3363a074057SBarry Smith   }
3373a074057SBarry Smith 
3389566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
3390fdc7489SMatthew Knepley   in.numberofedges = eEnd - eStart;
3400fdc7489SMatthew Knepley   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberofedges > 0) {
3410fdc7489SMatthew Knepley     in.edgelist       = new int[in.numberofedges * 2];
3420fdc7489SMatthew Knepley     in.edgemarkerlist = new int[in.numberofedges];
3430fdc7489SMatthew Knepley     for (e = eStart; e < eEnd; ++e) {
3440fdc7489SMatthew Knepley       const PetscInt  idx = e - eStart;
3450fdc7489SMatthew Knepley       const PetscInt *cone;
3460fdc7489SMatthew Knepley       PetscInt        coneSize, val;
3470fdc7489SMatthew Knepley 
3489566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, e, &coneSize));
3499566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, e, &cone));
3500fdc7489SMatthew Knepley       in.edgelist[idx * 2]     = cone[0] - vStart;
3510fdc7489SMatthew Knepley       in.edgelist[idx * 2 + 1] = cone[1] - vStart;
3520fdc7489SMatthew Knepley 
3539566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, e, &val));
354af226901SMatthew G. Knepley       if (val != defVal) in.edgemarkerlist[idx] = (int)val;
3550fdc7489SMatthew Knepley     }
3560fdc7489SMatthew Knepley   }
3570fdc7489SMatthew Knepley 
3589566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
3590fdc7489SMatthew Knepley   in.numberoffacets = fEnd - fStart;
3600fdc7489SMatthew Knepley   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberoffacets > 0) {
3610fdc7489SMatthew Knepley     in.facetlist       = new tetgenio::facet[in.numberoffacets];
3620fdc7489SMatthew Knepley     in.facetmarkerlist = new int[in.numberoffacets];
3630fdc7489SMatthew Knepley     for (f = fStart; f < fEnd; ++f) {
3640fdc7489SMatthew Knepley       const PetscInt idx    = f - fStart;
3650fdc7489SMatthew Knepley       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, val;
3660fdc7489SMatthew Knepley 
3670fdc7489SMatthew Knepley       in.facetlist[idx].numberofpolygons = 1;
3680fdc7489SMatthew Knepley       in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
3690fdc7489SMatthew Knepley       in.facetlist[idx].numberofholes    = 0;
3700fdc7489SMatthew Knepley       in.facetlist[idx].holelist         = NULL;
3710fdc7489SMatthew Knepley 
3729566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points));
3730fdc7489SMatthew Knepley       for (p = 0; p < numPoints * 2; p += 2) {
3740fdc7489SMatthew Knepley         const PetscInt point = points[p];
3750fdc7489SMatthew Knepley         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
3760fdc7489SMatthew Knepley       }
3770fdc7489SMatthew Knepley 
3780fdc7489SMatthew Knepley       tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
3790fdc7489SMatthew Knepley       poly->numberofvertices  = numVertices;
3800fdc7489SMatthew Knepley       poly->vertexlist        = new int[poly->numberofvertices];
3810fdc7489SMatthew Knepley       for (v = 0; v < numVertices; ++v) {
3820fdc7489SMatthew Knepley         const PetscInt vIdx = points[v] - vStart;
3830fdc7489SMatthew Knepley         poly->vertexlist[v] = vIdx;
3840fdc7489SMatthew Knepley       }
3850fdc7489SMatthew Knepley 
3869566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, f, &val));
387af226901SMatthew G. Knepley       if (val != defVal) in.facetmarkerlist[idx] = (int)val;
3880fdc7489SMatthew Knepley 
3899566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points));
3900fdc7489SMatthew Knepley     }
3910fdc7489SMatthew Knepley   }
3920fdc7489SMatthew Knepley 
3939566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
3943a074057SBarry Smith   in.numberofcorners       = 4;
3953a074057SBarry Smith   in.numberoftetrahedra    = cEnd - cStart;
3963a074057SBarry Smith   in.tetrahedronvolumelist = (double *)maxVolumes;
3973a074057SBarry Smith   if (in.numberoftetrahedra > 0) {
3983a074057SBarry Smith     in.tetrahedronlist = new int[in.numberoftetrahedra * in.numberofcorners];
3993a074057SBarry Smith     for (c = cStart; c < cEnd; ++c) {
4003a074057SBarry Smith       const PetscInt idx     = c - cStart;
4013a074057SBarry Smith       PetscInt      *closure = NULL;
4023a074057SBarry Smith       PetscInt       closureSize;
4033a074057SBarry Smith 
4049566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
40563a3b9bcSJacob Faibussowitsch       PetscCheck(!(closureSize != 5) || !(closureSize != 15), comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %" PetscInt_FMT " vertices in closure", closureSize);
4060fdc7489SMatthew Knepley       for (v = 0; v < 4; ++v) in.tetrahedronlist[idx * in.numberofcorners + v] = closure[(v + closureSize - 4) * 2] - vStart;
4079566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
4083a074057SBarry Smith     }
4093a074057SBarry Smith   }
4100fdc7489SMatthew Knepley 
411dd400576SPatrick Sanan   if (rank == 0) {
4123a074057SBarry Smith     char args[32];
4133a074057SBarry Smith 
4143a074057SBarry Smith     /* Take away 'Q' for verbose output */
4159566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(args, "qezQra"));
4163a074057SBarry Smith     ::tetrahedralize(args, &in, &out);
4173a074057SBarry Smith   }
4183a074057SBarry Smith 
4190fdc7489SMatthew Knepley   in.tetrahedronvolumelist = NULL;
4203a074057SBarry Smith   {
4213a074057SBarry Smith     const PetscInt numCorners  = 4;
4223a074057SBarry Smith     const PetscInt numCells    = out.numberoftetrahedra;
4233a074057SBarry Smith     const PetscInt numVertices = out.numberofpoints;
424a4a685f2SJacob Faibussowitsch     PetscReal     *meshCoords  = NULL;
425a4a685f2SJacob Faibussowitsch     PetscInt      *cells       = NULL;
4260fdc7489SMatthew Knepley     PetscBool      interpolate = isInterpolated == DMPLEX_INTERPOLATED_FULL ? PETSC_TRUE : PETSC_FALSE;
4273a074057SBarry Smith 
428a4a685f2SJacob Faibussowitsch     if (sizeof(PetscReal) == sizeof(out.pointlist[0])) {
429a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *)out.pointlist;
430a4a685f2SJacob Faibussowitsch     } else {
431a4a685f2SJacob Faibussowitsch       PetscInt i;
432a4a685f2SJacob Faibussowitsch 
433a4a685f2SJacob Faibussowitsch       meshCoords = new PetscReal[dim * numVertices];
4340fdc7489SMatthew Knepley       for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal)out.pointlist[i];
435a4a685f2SJacob Faibussowitsch     }
436a4a685f2SJacob Faibussowitsch     if (sizeof(PetscInt) == sizeof(out.tetrahedronlist[0])) {
437a4a685f2SJacob Faibussowitsch       cells = (PetscInt *)out.tetrahedronlist;
438a4a685f2SJacob Faibussowitsch     } else {
439a4a685f2SJacob Faibussowitsch       PetscInt i;
440a4a685f2SJacob Faibussowitsch 
441a4a685f2SJacob Faibussowitsch       cells = new PetscInt[numCells * numCorners];
4420fdc7489SMatthew Knepley       for (i = 0; i < numCells * numCorners; ++i) cells[i] = (PetscInt)out.tetrahedronlist[i];
443a4a685f2SJacob Faibussowitsch     }
444a4a685f2SJacob Faibussowitsch 
4459566063dSJacob Faibussowitsch     PetscCall(DMPlexInvertCells_Tetgen(numCells, numCorners, cells));
4469566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined));
4470fdc7489SMatthew Knepley     if (sizeof(PetscReal) != sizeof(out.pointlist[0])) { delete[] meshCoords; }
4480fdc7489SMatthew Knepley     if (sizeof(PetscInt) != sizeof(out.tetrahedronlist[0])) { delete[] cells; }
4490fdc7489SMatthew Knepley 
4503a074057SBarry Smith     /* Set labels */
4519566063dSJacob Faibussowitsch     PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dmRefined));
4523a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
453*48a46eb9SPierre Jolivet       if (out.pointmarkerlist[v]) PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, v + numCells, out.pointmarkerlist[v]));
4543a074057SBarry Smith     }
4553a074057SBarry Smith     if (interpolate) {
4560fdc7489SMatthew Knepley       PetscInt e, f;
4573a074057SBarry Smith 
4580fdc7489SMatthew Knepley       for (e = 0; e < out.numberofedges; ++e) {
4593a074057SBarry Smith         if (out.edgemarkerlist[e]) {
4603a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells};
4613a074057SBarry Smith           const PetscInt *edges;
4623a074057SBarry Smith           PetscInt        numEdges;
4633a074057SBarry Smith 
4649566063dSJacob Faibussowitsch           PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges));
46563a3b9bcSJacob Faibussowitsch           PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
4669566063dSJacob Faibussowitsch           PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, edges[0], out.edgemarkerlist[e]));
4679566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges));
4683a074057SBarry Smith         }
4693a074057SBarry Smith       }
4700fdc7489SMatthew Knepley       for (f = 0; f < out.numberoftrifaces; ++f) {
4713a074057SBarry Smith         if (out.trifacemarkerlist[f]) {
4723a074057SBarry Smith           const PetscInt  vertices[3] = {out.trifacelist[f * 3 + 0] + numCells, out.trifacelist[f * 3 + 1] + numCells, out.trifacelist[f * 3 + 2] + numCells};
4733a074057SBarry Smith           const PetscInt *faces;
4743a074057SBarry Smith           PetscInt        numFaces;
4753a074057SBarry Smith 
4769566063dSJacob Faibussowitsch           PetscCall(DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces));
47763a3b9bcSJacob Faibussowitsch           PetscCheck(numFaces == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces);
4789566063dSJacob Faibussowitsch           PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, faces[0], out.trifacemarkerlist[f]));
4799566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces));
4803a074057SBarry Smith         }
4813a074057SBarry Smith       }
4823a074057SBarry Smith     }
4830fdc7489SMatthew Knepley 
4849566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
4859318fe57SMatthew G. Knepley     if (modelObj) {
4860fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
4870fdc7489SMatthew Knepley       DMLabel   bodyLabel;
4880fdc7489SMatthew Knepley       PetscInt  cStart, cEnd, c, eStart, eEnd, fStart, fEnd;
489c1cad2e7SMatthew G. Knepley       PetscBool islite = PETSC_FALSE;
4900fdc7489SMatthew Knepley       ego      *bodies;
4910fdc7489SMatthew Knepley       ego       model, geom;
4920fdc7489SMatthew Knepley       int       Nb, oclass, mtype, *senses;
4930fdc7489SMatthew Knepley 
4940fdc7489SMatthew Knepley       /* Get Attached EGADS Model from Original DMPlex */
4959566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
496c1cad2e7SMatthew G. Knepley       if (modelObj) {
4979566063dSJacob Faibussowitsch         PetscCall(PetscContainerGetPointer(modelObj, (void **)&model));
4989566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
4990fdc7489SMatthew Knepley         /* Transfer EGADS Model to Volumetric Mesh */
5009566063dSJacob Faibussowitsch         PetscCall(PetscObjectCompose((PetscObject)*dmRefined, "EGADS Model", (PetscObject)modelObj));
501c1cad2e7SMatthew G. Knepley       } else {
5029566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSLite Model", (PetscObject *)&modelObj));
503c1cad2e7SMatthew G. Knepley         if (modelObj) {
5049566063dSJacob Faibussowitsch           PetscCall(PetscContainerGetPointer(modelObj, (void **)&model));
5059566063dSJacob Faibussowitsch           PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
506c1cad2e7SMatthew G. Knepley           /* Transfer EGADS Model to Volumetric Mesh */
5079566063dSJacob Faibussowitsch           PetscCall(PetscObjectCompose((PetscObject)*dmRefined, "EGADSLite Model", (PetscObject)modelObj));
508c1cad2e7SMatthew G. Knepley           islite = PETSC_TRUE;
509c1cad2e7SMatthew G. Knepley         }
510c1cad2e7SMatthew G. Knepley       }
511c1cad2e7SMatthew G. Knepley       if (!modelObj) goto skip_egads;
5120fdc7489SMatthew Knepley 
5130fdc7489SMatthew Knepley       /* Set Cell Labels */
5149566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dmRefined, "EGADS Body ID", &bodyLabel));
5159566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(*dmRefined, 0, &cStart, &cEnd));
5169566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(*dmRefined, 1, &fStart, &fEnd));
5179566063dSJacob Faibussowitsch       PetscCall(DMPlexGetDepthStratum(*dmRefined, 1, &eStart, &eEnd));
5180fdc7489SMatthew Knepley 
5190fdc7489SMatthew Knepley       for (c = cStart; c < cEnd; ++c) {
5200fdc7489SMatthew Knepley         PetscReal centroid[3] = {0., 0., 0.};
5210fdc7489SMatthew Knepley         PetscInt  b;
5220fdc7489SMatthew Knepley 
5230fdc7489SMatthew Knepley         /* Deterimine what body the cell's centroid is located in */
5240fdc7489SMatthew Knepley         if (!interpolate) {
5250fdc7489SMatthew Knepley           PetscSection coordSection;
5260fdc7489SMatthew Knepley           Vec          coordinates;
5270fdc7489SMatthew Knepley           PetscScalar *coords = NULL;
5280fdc7489SMatthew Knepley           PetscInt     coordSize, s, d;
5290fdc7489SMatthew Knepley 
5309566063dSJacob Faibussowitsch           PetscCall(DMGetCoordinatesLocal(*dmRefined, &coordinates));
5319566063dSJacob Faibussowitsch           PetscCall(DMGetCoordinateSection(*dmRefined, &coordSection));
5329566063dSJacob Faibussowitsch           PetscCall(DMPlexVecGetClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords));
5339371c9d4SSatish Balay           for (s = 0; s < coordSize; ++s)
5349371c9d4SSatish Balay             for (d = 0; d < dim; ++d) centroid[d] += coords[s * dim + d];
5359566063dSJacob Faibussowitsch           PetscCall(DMPlexVecRestoreClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords));
5361baa6e33SBarry Smith         } else PetscCall(DMPlexComputeCellGeometryFVM(*dmRefined, c, NULL, centroid, NULL));
5370fdc7489SMatthew Knepley         for (b = 0; b < Nb; ++b) {
5389371c9d4SSatish Balay           if (islite) {
5399371c9d4SSatish Balay             if (EGlite_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;
5409371c9d4SSatish Balay           } else {
5419371c9d4SSatish Balay             if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;
5429371c9d4SSatish Balay           }
5430fdc7489SMatthew Knepley         }
5440fdc7489SMatthew Knepley         if (b < Nb) {
5450fdc7489SMatthew Knepley           PetscInt  cval    = b, eVal, fVal;
5460fdc7489SMatthew Knepley           PetscInt *closure = NULL, Ncl, cl;
5470fdc7489SMatthew Knepley 
5489566063dSJacob Faibussowitsch           PetscCall(DMLabelSetValue(bodyLabel, c, cval));
5499566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure));
5500fdc7489SMatthew Knepley           for (cl = 0; cl < Ncl; cl += 2) {
5510fdc7489SMatthew Knepley             const PetscInt p = closure[cl];
5520fdc7489SMatthew Knepley 
5530fdc7489SMatthew Knepley             if (p >= eStart && p < eEnd) {
5549566063dSJacob Faibussowitsch               PetscCall(DMLabelGetValue(bodyLabel, p, &eVal));
5559566063dSJacob Faibussowitsch               if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
5560fdc7489SMatthew Knepley             }
5570fdc7489SMatthew Knepley             if (p >= fStart && p < fEnd) {
5589566063dSJacob Faibussowitsch               PetscCall(DMLabelGetValue(bodyLabel, p, &fVal));
5599566063dSJacob Faibussowitsch               if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
5600fdc7489SMatthew Knepley             }
5610fdc7489SMatthew Knepley           }
5629566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure));
5630fdc7489SMatthew Knepley         }
5640fdc7489SMatthew Knepley       }
565c1cad2e7SMatthew G. Knepley     skip_egads:;
5660fdc7489SMatthew Knepley #endif
5679318fe57SMatthew G. Knepley     }
5689566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE));
5693a074057SBarry Smith   }
5703a074057SBarry Smith   PetscFunctionReturn(0);
5713a074057SBarry Smith }
572