xref: /petsc/src/dm/impls/plex/generators/tetgen/tetgenerate.cxx (revision 41e9d8b5bae7cddec70f4add541c9e66dc5b0939)
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 */
16a4a685f2SJacob Faibussowitsch static PetscErrorCode DMPlexInvertCells_Tetgen(PetscInt numCells, PetscInt numCorners, PetscInt cells[])
173a074057SBarry Smith {
183a074057SBarry Smith   PetscInt bound = numCells*numCorners, coff;
193a074057SBarry Smith 
203a074057SBarry Smith   PetscFunctionBegin;
21a4a685f2SJacob Faibussowitsch #define SWAP(a,b) do { PetscInt tmp = (a); (a) = (b); (b) = tmp; } while (0)
2296ca5757SLisandro Dalcin   for (coff = 0; coff < bound; coff += numCorners) SWAP(cells[coff],cells[coff+1]);
2396ca5757SLisandro Dalcin #undef SWAP
243a074057SBarry Smith   PetscFunctionReturn(0);
253a074057SBarry Smith }
263a074057SBarry Smith 
273a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerate_Tetgen(DM boundary, PetscBool interpolate, DM *dm)
283a074057SBarry Smith {
293a074057SBarry Smith   MPI_Comm               comm;
303a074057SBarry Smith   const PetscInt         dim = 3;
313a074057SBarry Smith   ::tetgenio             in;
323a074057SBarry Smith   ::tetgenio             out;
339318fe57SMatthew G. Knepley   PetscContainer         modelObj;
340fdc7489SMatthew Knepley   DMUniversalLabel       universal;
35af226901SMatthew G. Knepley   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, defVal;
360fdc7489SMatthew Knepley   DMPlexInterpolatedFlag isInterpolated;
373a074057SBarry Smith   PetscMPIInt            rank;
383a074057SBarry Smith 
393a074057SBarry Smith   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)boundary,&comm));
419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
429566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolatedCollective(boundary, &isInterpolated));
439566063dSJacob Faibussowitsch   PetscCall(DMUniversalLabelCreate(boundary, &universal));
44af226901SMatthew G. Knepley   PetscCall(DMLabelGetDefaultValue(universal->label, &defVal));
453a074057SBarry Smith 
469566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd));
473a074057SBarry Smith   in.numberofpoints = vEnd - vStart;
483a074057SBarry Smith   if (in.numberofpoints > 0) {
493a074057SBarry Smith     PetscSection       coordSection;
503a074057SBarry Smith     Vec                coordinates;
510fdc7489SMatthew Knepley     const PetscScalar *array;
523a074057SBarry Smith 
533a074057SBarry Smith     in.pointlist       = new double[in.numberofpoints*dim];
543a074057SBarry Smith     in.pointmarkerlist = new int[in.numberofpoints];
553a074057SBarry Smith 
56*41e9d8b5SMatthew G. Knepley     PetscCall(PetscArrayzero(in.pointmarkerlist, (size_t) in.numberofpoints));
579566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(boundary, &coordinates));
589566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(boundary, &coordSection));
599566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coordinates, &array));
603a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
613a074057SBarry Smith       const PetscInt idx = v - vStart;
620fdc7489SMatthew Knepley       PetscInt       off, d, val;
633a074057SBarry Smith 
649566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
653a074057SBarry Smith       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
669566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, v, &val));
67af226901SMatthew G. Knepley       if (val != defVal) in.pointmarkerlist[idx] = (int) val;
683a074057SBarry Smith     }
699566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coordinates, &array));
703a074057SBarry Smith   }
713a074057SBarry Smith 
729566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(boundary, 1, &eStart, &eEnd));
730fdc7489SMatthew Knepley   in.numberofedges = eEnd - eStart;
740fdc7489SMatthew Knepley   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberofedges > 0) {
750fdc7489SMatthew Knepley     in.edgelist       = new int[in.numberofedges * 2];
760fdc7489SMatthew Knepley     in.edgemarkerlist = new int[in.numberofedges];
770fdc7489SMatthew Knepley     for (e = eStart; e < eEnd; ++e) {
780fdc7489SMatthew Knepley       const PetscInt  idx = e - eStart;
790fdc7489SMatthew Knepley       const PetscInt *cone;
800fdc7489SMatthew Knepley       PetscInt        coneSize, val;
810fdc7489SMatthew Knepley 
829566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(boundary, e, &coneSize));
839566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(boundary, e, &cone));
840fdc7489SMatthew Knepley       in.edgelist[idx*2]     = cone[0] - vStart;
850fdc7489SMatthew Knepley       in.edgelist[idx*2 + 1] = cone[1] - vStart;
860fdc7489SMatthew Knepley 
879566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, e, &val));
88af226901SMatthew G. Knepley       if (val != defVal) in.edgemarkerlist[idx] = (int) val;
890fdc7489SMatthew Knepley     }
900fdc7489SMatthew Knepley   }
910fdc7489SMatthew Knepley 
929566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd));
933a074057SBarry Smith   in.numberoffacets = fEnd - fStart;
943a074057SBarry Smith   if (in.numberoffacets > 0) {
953a074057SBarry Smith     in.facetlist       = new tetgenio::facet[in.numberoffacets];
963a074057SBarry Smith     in.facetmarkerlist = new int[in.numberoffacets];
973a074057SBarry Smith     for (f = fStart; f < fEnd; ++f) {
983a074057SBarry Smith       const PetscInt idx    = f - fStart;
990fdc7489SMatthew Knepley       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, val = -1;
1003a074057SBarry Smith 
1013a074057SBarry Smith       in.facetlist[idx].numberofpolygons = 1;
1023a074057SBarry Smith       in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
1033a074057SBarry Smith       in.facetlist[idx].numberofholes    = 0;
1043a074057SBarry Smith       in.facetlist[idx].holelist         = NULL;
1053a074057SBarry Smith 
1069566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points));
1073a074057SBarry Smith       for (p = 0; p < numPoints*2; p += 2) {
1083a074057SBarry Smith         const PetscInt point = points[p];
1093a074057SBarry Smith         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
1103a074057SBarry Smith       }
1113a074057SBarry Smith 
1123a074057SBarry Smith       tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
1133a074057SBarry Smith       poly->numberofvertices = numVertices;
1143a074057SBarry Smith       poly->vertexlist       = new int[poly->numberofvertices];
1153a074057SBarry Smith       for (v = 0; v < numVertices; ++v) {
1163a074057SBarry Smith         const PetscInt vIdx = points[v] - vStart;
1173a074057SBarry Smith         poly->vertexlist[v] = vIdx;
1183a074057SBarry Smith       }
1199566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, f, &val));
120af226901SMatthew G. Knepley       if (val != defVal) in.facetmarkerlist[idx] = (int) val;
1219566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points));
1223a074057SBarry Smith     }
1233a074057SBarry Smith   }
124dd400576SPatrick Sanan   if (rank == 0) {
1250fdc7489SMatthew Knepley     DM_Plex *mesh = (DM_Plex *) boundary->data;
1263a074057SBarry Smith     char     args[32];
1273a074057SBarry Smith 
1283a074057SBarry Smith     /* Take away 'Q' for verbose output */
1290fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
1309566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(args, "pqezQY"));
1310fdc7489SMatthew Knepley #else
1329566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(args, "pqezQ"));
1330fdc7489SMatthew Knepley #endif
1343a074057SBarry Smith     if (mesh->tetgenOpts) {::tetrahedralize(mesh->tetgenOpts, &in, &out);}
1353a074057SBarry Smith     else                  {::tetrahedralize(args, &in, &out);}
1363a074057SBarry Smith   }
1373a074057SBarry Smith   {
1383a074057SBarry Smith     const PetscInt   numCorners  = 4;
1393a074057SBarry Smith     const PetscInt   numCells    = out.numberoftetrahedra;
1403a074057SBarry Smith     const PetscInt   numVertices = out.numberofpoints;
141a4a685f2SJacob Faibussowitsch     PetscReal        *meshCoords = NULL;
142a4a685f2SJacob Faibussowitsch     PetscInt         *cells      = NULL;
143a4a685f2SJacob Faibussowitsch 
144a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) == sizeof (out.pointlist[0])) {
145a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *) out.pointlist;
146a4a685f2SJacob Faibussowitsch     } else {
147a4a685f2SJacob Faibussowitsch       PetscInt i;
148a4a685f2SJacob Faibussowitsch 
149a4a685f2SJacob Faibussowitsch       meshCoords = new PetscReal[dim * numVertices];
1500fdc7489SMatthew Knepley       for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out.pointlist[i];
151a4a685f2SJacob Faibussowitsch     }
152a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) == sizeof (out.tetrahedronlist[0])) {
153a4a685f2SJacob Faibussowitsch       cells = (PetscInt *) out.tetrahedronlist;
154a4a685f2SJacob Faibussowitsch     } else {
155a4a685f2SJacob Faibussowitsch       PetscInt i;
156a4a685f2SJacob Faibussowitsch 
157a4a685f2SJacob Faibussowitsch       cells = new PetscInt[numCells * numCorners];
1580fdc7489SMatthew Knepley       for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt) out.tetrahedronlist[i];
159a4a685f2SJacob Faibussowitsch     }
1603a074057SBarry Smith 
1619566063dSJacob Faibussowitsch     PetscCall(DMPlexInvertCells_Tetgen(numCells, numCorners, cells));
1629566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm));
1630fdc7489SMatthew Knepley 
1643a074057SBarry Smith     /* Set labels */
1659566063dSJacob Faibussowitsch     PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dm));
1663a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
1673a074057SBarry Smith       if (out.pointmarkerlist[v]) {
1689566063dSJacob Faibussowitsch         PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, v+numCells, out.pointmarkerlist[v]));
1693a074057SBarry Smith       }
1703a074057SBarry Smith     }
1713a074057SBarry Smith     if (interpolate) {
1723a074057SBarry Smith       PetscInt e;
1733a074057SBarry Smith 
1743a074057SBarry Smith       for (e = 0; e < out.numberofedges; e++) {
1753a074057SBarry Smith         if (out.edgemarkerlist[e]) {
1763a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
1773a074057SBarry Smith           const PetscInt *edges;
1783a074057SBarry Smith           PetscInt        numEdges;
1793a074057SBarry Smith 
1809566063dSJacob Faibussowitsch           PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges));
18163a3b9bcSJacob Faibussowitsch           PetscCheck(numEdges == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
1829566063dSJacob Faibussowitsch           PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, edges[0], out.edgemarkerlist[e]));
1839566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges));
1843a074057SBarry Smith         }
1853a074057SBarry Smith       }
1863a074057SBarry Smith       for (f = 0; f < out.numberoftrifaces; f++) {
1873a074057SBarry Smith         if (out.trifacemarkerlist[f]) {
1883a074057SBarry Smith           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
1893a074057SBarry Smith           const PetscInt *faces;
1903a074057SBarry Smith           PetscInt        numFaces;
1913a074057SBarry Smith 
1929566063dSJacob Faibussowitsch           PetscCall(DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces));
19363a3b9bcSJacob Faibussowitsch           PetscCheck(numFaces == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces);
1949566063dSJacob Faibussowitsch           PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, faces[0], out.trifacemarkerlist[f]));
1959566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces));
1963a074057SBarry Smith         }
1973a074057SBarry Smith       }
1983a074057SBarry Smith     }
1990fdc7489SMatthew Knepley 
2009566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject) boundary, "EGADS Model", (PetscObject *) &modelObj));
2019318fe57SMatthew G. Knepley     if (modelObj) {
2020fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
2030fdc7489SMatthew Knepley       DMLabel        bodyLabel;
2040fdc7489SMatthew Knepley       PetscInt       cStart, cEnd, c, eStart, eEnd, fStart, fEnd;
205c1cad2e7SMatthew G. Knepley       PetscBool      islite = PETSC_FALSE;
2060fdc7489SMatthew Knepley       ego           *bodies;
2070fdc7489SMatthew Knepley       ego            model, geom;
2080fdc7489SMatthew Knepley       int            Nb, oclass, mtype, *senses;
2090fdc7489SMatthew Knepley 
2100fdc7489SMatthew Knepley       /* Get Attached EGADS Model from Original DMPlex */
2119566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery((PetscObject) boundary, "EGADS Model", (PetscObject *) &modelObj));
212c1cad2e7SMatthew G. Knepley       if (modelObj) {
2139566063dSJacob Faibussowitsch         PetscCall(PetscContainerGetPointer(modelObj, (void **) &model));
2149566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
2150fdc7489SMatthew Knepley         /* Transfer EGADS Model to Volumetric Mesh */
2169566063dSJacob Faibussowitsch         PetscCall(PetscObjectCompose((PetscObject) *dm, "EGADS Model", (PetscObject) modelObj));
217c1cad2e7SMatthew G. Knepley       } else {
2189566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject) boundary, "EGADSLite Model", (PetscObject *) &modelObj));
219c1cad2e7SMatthew G. Knepley         if (modelObj) {
2209566063dSJacob Faibussowitsch           PetscCall(PetscContainerGetPointer(modelObj, (void **) &model));
2219566063dSJacob Faibussowitsch           PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
222c1cad2e7SMatthew G. Knepley           /* Transfer EGADS Model to Volumetric Mesh */
2239566063dSJacob Faibussowitsch           PetscCall(PetscObjectCompose((PetscObject) *dm, "EGADSLite Model", (PetscObject) modelObj));
224c1cad2e7SMatthew G. Knepley           islite = PETSC_TRUE;
225c1cad2e7SMatthew G. Knepley         }
226c1cad2e7SMatthew G. Knepley       }
227c1cad2e7SMatthew G. Knepley       if (!modelObj) goto skip_egads;
2280fdc7489SMatthew Knepley 
2290fdc7489SMatthew Knepley       /* Set Cell Labels */
2309566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dm, "EGADS Body ID", &bodyLabel));
2319566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd));
2329566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd));
2339566063dSJacob Faibussowitsch       PetscCall(DMPlexGetDepthStratum(*dm, 1, &eStart, &eEnd));
2340fdc7489SMatthew Knepley 
2350fdc7489SMatthew Knepley       for (c = cStart; c < cEnd; ++c) {
2360fdc7489SMatthew Knepley         PetscReal centroid[3] = {0., 0., 0.};
2370fdc7489SMatthew Knepley         PetscInt  b;
2380fdc7489SMatthew Knepley 
2390fdc7489SMatthew Knepley         /* Deterimine what body the cell's centroid is located in */
2400fdc7489SMatthew Knepley         if (!interpolate) {
2410fdc7489SMatthew Knepley           PetscSection   coordSection;
2420fdc7489SMatthew Knepley           Vec            coordinates;
2430fdc7489SMatthew Knepley           PetscScalar   *coords = NULL;
2440fdc7489SMatthew Knepley           PetscInt       coordSize, s, d;
2450fdc7489SMatthew Knepley 
2469566063dSJacob Faibussowitsch           PetscCall(DMGetCoordinatesLocal(*dm, &coordinates));
2479566063dSJacob Faibussowitsch           PetscCall(DMGetCoordinateSection(*dm, &coordSection));
2489566063dSJacob Faibussowitsch           PetscCall(DMPlexVecGetClosure(*dm, coordSection, coordinates, c, &coordSize, &coords));
2490fdc7489SMatthew Knepley           for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d];
2509566063dSJacob Faibussowitsch           PetscCall(DMPlexVecRestoreClosure(*dm, coordSection, coordinates, c, &coordSize, &coords));
2511baa6e33SBarry Smith         } else PetscCall(DMPlexComputeCellGeometryFVM(*dm, c, NULL, centroid, NULL));
2520fdc7489SMatthew Knepley         for (b = 0; b < Nb; ++b) {
253c1cad2e7SMatthew G. Knepley           if (islite) {if (EGlite_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
254c1cad2e7SMatthew G. Knepley           else        {if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
2550fdc7489SMatthew Knepley         }
2560fdc7489SMatthew Knepley         if (b < Nb) {
2570fdc7489SMatthew Knepley           PetscInt   cval = b, eVal, fVal;
2580fdc7489SMatthew Knepley           PetscInt *closure = NULL, Ncl, cl;
2590fdc7489SMatthew Knepley 
2609566063dSJacob Faibussowitsch           PetscCall(DMLabelSetValue(bodyLabel, c, cval));
2619566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure));
2620fdc7489SMatthew Knepley           for (cl = 0; cl < Ncl; cl += 2) {
2630fdc7489SMatthew Knepley             const PetscInt p = closure[cl];
2640fdc7489SMatthew Knepley 
2650fdc7489SMatthew Knepley             if (p >= eStart && p < eEnd) {
2669566063dSJacob Faibussowitsch               PetscCall(DMLabelGetValue(bodyLabel, p, &eVal));
2679566063dSJacob Faibussowitsch               if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
2680fdc7489SMatthew Knepley             }
2690fdc7489SMatthew Knepley             if (p >= fStart && p < fEnd) {
2709566063dSJacob Faibussowitsch               PetscCall(DMLabelGetValue(bodyLabel, p, &fVal));
2719566063dSJacob Faibussowitsch               if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
2720fdc7489SMatthew Knepley             }
2730fdc7489SMatthew Knepley           }
2749566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure));
2750fdc7489SMatthew Knepley         }
2760fdc7489SMatthew Knepley       }
277c1cad2e7SMatthew G. Knepley skip_egads: ;
2780fdc7489SMatthew Knepley #endif
2799318fe57SMatthew G. Knepley     }
2809566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
2813a074057SBarry Smith   }
2829566063dSJacob Faibussowitsch   PetscCall(DMUniversalLabelDestroy(&universal));
2833a074057SBarry Smith   PetscFunctionReturn(0);
2843a074057SBarry Smith }
2853a074057SBarry Smith 
2863a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexRefine_Tetgen(DM dm, double *maxVolumes, DM *dmRefined)
2873a074057SBarry Smith {
2883a074057SBarry Smith   MPI_Comm               comm;
2893a074057SBarry Smith   const PetscInt         dim = 3;
2903a074057SBarry Smith   ::tetgenio             in;
2913a074057SBarry Smith   ::tetgenio             out;
2929318fe57SMatthew G. Knepley   PetscContainer         modelObj;
2930fdc7489SMatthew Knepley   DMUniversalLabel       universal;
294af226901SMatthew G. Knepley   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, cStart, cEnd, c, defVal;
2950fdc7489SMatthew Knepley   DMPlexInterpolatedFlag isInterpolated;
2963a074057SBarry Smith   PetscMPIInt            rank;
2973a074057SBarry Smith 
2983a074057SBarry Smith   PetscFunctionBegin;
2999566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
3009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3019566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolatedCollective(dm, &isInterpolated));
3029566063dSJacob Faibussowitsch   PetscCall(DMUniversalLabelCreate(dm, &universal));
303af226901SMatthew G. Knepley   PetscCall(DMLabelGetDefaultValue(universal->label, &defVal));
3043a074057SBarry Smith 
3059566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
3063a074057SBarry Smith   in.numberofpoints = vEnd - vStart;
3073a074057SBarry Smith   if (in.numberofpoints > 0) {
3083a074057SBarry Smith     PetscSection coordSection;
3093a074057SBarry Smith     Vec          coordinates;
3103a074057SBarry Smith     PetscScalar *array;
3113a074057SBarry Smith 
3123a074057SBarry Smith     in.pointlist       = new double[in.numberofpoints*dim];
3133a074057SBarry Smith     in.pointmarkerlist = new int[in.numberofpoints];
3143a074057SBarry Smith 
315*41e9d8b5SMatthew G. Knepley     PetscCall(PetscArrayzero(in.pointmarkerlist, (size_t) in.numberofpoints));
3169566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
3179566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dm, &coordSection));
3189566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &array));
3193a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
3203a074057SBarry Smith       const PetscInt idx = v - vStart;
3210fdc7489SMatthew Knepley       PetscInt       off, d, val;
3223a074057SBarry Smith 
3239566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
3243a074057SBarry Smith       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
3259566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, v, &val));
326af226901SMatthew G. Knepley       if (val != defVal) in.pointmarkerlist[idx] = (int) val;
3273a074057SBarry Smith     }
3289566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &array));
3293a074057SBarry Smith   }
3303a074057SBarry Smith 
3319566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
3320fdc7489SMatthew Knepley   in.numberofedges = eEnd - eStart;
3330fdc7489SMatthew Knepley   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberofedges > 0) {
3340fdc7489SMatthew Knepley     in.edgelist       = new int[in.numberofedges * 2];
3350fdc7489SMatthew Knepley     in.edgemarkerlist = new int[in.numberofedges];
3360fdc7489SMatthew Knepley     for (e = eStart; e < eEnd; ++e) {
3370fdc7489SMatthew Knepley       const PetscInt  idx = e - eStart;
3380fdc7489SMatthew Knepley       const PetscInt *cone;
3390fdc7489SMatthew Knepley       PetscInt        coneSize, val;
3400fdc7489SMatthew Knepley 
3419566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, e, &coneSize));
3429566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, e, &cone));
3430fdc7489SMatthew Knepley       in.edgelist[idx*2]     = cone[0] - vStart;
3440fdc7489SMatthew Knepley       in.edgelist[idx*2 + 1] = cone[1] - vStart;
3450fdc7489SMatthew Knepley 
3469566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, e, &val));
347af226901SMatthew G. Knepley       if (val != defVal) in.edgemarkerlist[idx] = (int) val;
3480fdc7489SMatthew Knepley     }
3490fdc7489SMatthew Knepley   }
3500fdc7489SMatthew Knepley 
3519566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
3520fdc7489SMatthew Knepley   in.numberoffacets = fEnd - fStart;
3530fdc7489SMatthew Knepley   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberoffacets > 0) {
3540fdc7489SMatthew Knepley     in.facetlist       = new tetgenio::facet[in.numberoffacets];
3550fdc7489SMatthew Knepley     in.facetmarkerlist = new int[in.numberoffacets];
3560fdc7489SMatthew Knepley     for (f = fStart; f < fEnd; ++f) {
3570fdc7489SMatthew Knepley       const PetscInt idx    = f - fStart;
3580fdc7489SMatthew Knepley       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, val;
3590fdc7489SMatthew Knepley 
3600fdc7489SMatthew Knepley       in.facetlist[idx].numberofpolygons = 1;
3610fdc7489SMatthew Knepley       in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
3620fdc7489SMatthew Knepley       in.facetlist[idx].numberofholes    = 0;
3630fdc7489SMatthew Knepley       in.facetlist[idx].holelist         = NULL;
3640fdc7489SMatthew Knepley 
3659566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points));
3660fdc7489SMatthew Knepley       for (p = 0; p < numPoints*2; p += 2) {
3670fdc7489SMatthew Knepley         const PetscInt point = points[p];
3680fdc7489SMatthew Knepley         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
3690fdc7489SMatthew Knepley       }
3700fdc7489SMatthew Knepley 
3710fdc7489SMatthew Knepley       tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
3720fdc7489SMatthew Knepley       poly->numberofvertices = numVertices;
3730fdc7489SMatthew Knepley       poly->vertexlist       = new int[poly->numberofvertices];
3740fdc7489SMatthew Knepley       for (v = 0; v < numVertices; ++v) {
3750fdc7489SMatthew Knepley         const PetscInt vIdx = points[v] - vStart;
3760fdc7489SMatthew Knepley         poly->vertexlist[v] = vIdx;
3770fdc7489SMatthew Knepley       }
3780fdc7489SMatthew Knepley 
3799566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(universal->label, f, &val));
380af226901SMatthew G. Knepley       if (val != defVal) in.facetmarkerlist[idx] = (int) val;
3810fdc7489SMatthew Knepley 
3829566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points));
3830fdc7489SMatthew Knepley     }
3840fdc7489SMatthew Knepley   }
3850fdc7489SMatthew Knepley 
3869566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
3873a074057SBarry Smith   in.numberofcorners       = 4;
3883a074057SBarry Smith   in.numberoftetrahedra    = cEnd - cStart;
3893a074057SBarry Smith   in.tetrahedronvolumelist = (double *) maxVolumes;
3903a074057SBarry Smith   if (in.numberoftetrahedra > 0) {
3913a074057SBarry Smith     in.tetrahedronlist = new int[in.numberoftetrahedra*in.numberofcorners];
3923a074057SBarry Smith     for (c = cStart; c < cEnd; ++c) {
3933a074057SBarry Smith       const PetscInt idx     = c - cStart;
3943a074057SBarry Smith       PetscInt      *closure = NULL;
3953a074057SBarry Smith       PetscInt       closureSize;
3963a074057SBarry Smith 
3979566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
39863a3b9bcSJacob 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);
3990fdc7489SMatthew Knepley       for (v = 0; v < 4; ++v) in.tetrahedronlist[idx*in.numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
4009566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
4013a074057SBarry Smith     }
4023a074057SBarry Smith   }
4030fdc7489SMatthew Knepley 
404dd400576SPatrick Sanan   if (rank == 0) {
4053a074057SBarry Smith     char args[32];
4063a074057SBarry Smith 
4073a074057SBarry Smith     /* Take away 'Q' for verbose output */
4089566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(args, "qezQra"));
4093a074057SBarry Smith     ::tetrahedralize(args, &in, &out);
4103a074057SBarry Smith   }
4113a074057SBarry Smith 
4120fdc7489SMatthew Knepley   in.tetrahedronvolumelist = NULL;
4133a074057SBarry Smith   {
4143a074057SBarry Smith     const PetscInt   numCorners  = 4;
4153a074057SBarry Smith     const PetscInt   numCells    = out.numberoftetrahedra;
4163a074057SBarry Smith     const PetscInt   numVertices = out.numberofpoints;
417a4a685f2SJacob Faibussowitsch     PetscReal        *meshCoords = NULL;
418a4a685f2SJacob Faibussowitsch     PetscInt         *cells      = NULL;
4190fdc7489SMatthew Knepley     PetscBool        interpolate = isInterpolated == DMPLEX_INTERPOLATED_FULL ? PETSC_TRUE : PETSC_FALSE;
4203a074057SBarry Smith 
421a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) == sizeof (out.pointlist[0])) {
422a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *) out.pointlist;
423a4a685f2SJacob Faibussowitsch     } else {
424a4a685f2SJacob Faibussowitsch       PetscInt i;
425a4a685f2SJacob Faibussowitsch 
426a4a685f2SJacob Faibussowitsch       meshCoords = new PetscReal[dim * numVertices];
4270fdc7489SMatthew Knepley       for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out.pointlist[i];
428a4a685f2SJacob Faibussowitsch     }
429a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) == sizeof (out.tetrahedronlist[0])) {
430a4a685f2SJacob Faibussowitsch       cells = (PetscInt *) out.tetrahedronlist;
431a4a685f2SJacob Faibussowitsch     } else {
432a4a685f2SJacob Faibussowitsch       PetscInt i;
433a4a685f2SJacob Faibussowitsch 
434a4a685f2SJacob Faibussowitsch       cells = new PetscInt[numCells * numCorners];
4350fdc7489SMatthew Knepley       for (i = 0; i < numCells * numCorners; ++i)cells[i] = (PetscInt) out.tetrahedronlist[i];
436a4a685f2SJacob Faibussowitsch     }
437a4a685f2SJacob Faibussowitsch 
4389566063dSJacob Faibussowitsch     PetscCall(DMPlexInvertCells_Tetgen(numCells, numCorners, cells));
4399566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined));
4400fdc7489SMatthew Knepley     if (sizeof (PetscReal) != sizeof (out.pointlist[0])) {delete [] meshCoords;}
4410fdc7489SMatthew Knepley     if (sizeof (PetscInt) != sizeof (out.tetrahedronlist[0])) {delete [] cells;}
4420fdc7489SMatthew Knepley 
4433a074057SBarry Smith     /* Set labels */
4449566063dSJacob Faibussowitsch     PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dmRefined));
4453a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
4463a074057SBarry Smith       if (out.pointmarkerlist[v]) {
4479566063dSJacob Faibussowitsch         PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, v+numCells, out.pointmarkerlist[v]));
4483a074057SBarry Smith       }
4493a074057SBarry Smith     }
4503a074057SBarry Smith     if (interpolate) {
4510fdc7489SMatthew Knepley       PetscInt e, f;
4523a074057SBarry Smith 
4530fdc7489SMatthew Knepley       for (e = 0; e < out.numberofedges; ++e) {
4543a074057SBarry Smith         if (out.edgemarkerlist[e]) {
4553a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
4563a074057SBarry Smith           const PetscInt *edges;
4573a074057SBarry Smith           PetscInt        numEdges;
4583a074057SBarry Smith 
4599566063dSJacob Faibussowitsch           PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges));
46063a3b9bcSJacob Faibussowitsch           PetscCheck(numEdges == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
4619566063dSJacob Faibussowitsch           PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, edges[0], out.edgemarkerlist[e]));
4629566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges));
4633a074057SBarry Smith         }
4643a074057SBarry Smith       }
4650fdc7489SMatthew Knepley       for (f = 0; f < out.numberoftrifaces; ++f) {
4663a074057SBarry Smith         if (out.trifacemarkerlist[f]) {
4673a074057SBarry Smith           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
4683a074057SBarry Smith           const PetscInt *faces;
4693a074057SBarry Smith           PetscInt        numFaces;
4703a074057SBarry Smith 
4719566063dSJacob Faibussowitsch           PetscCall(DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces));
47263a3b9bcSJacob Faibussowitsch           PetscCheck(numFaces == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces);
4739566063dSJacob Faibussowitsch           PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, faces[0], out.trifacemarkerlist[f]));
4749566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces));
4753a074057SBarry Smith         }
4763a074057SBarry Smith       }
4773a074057SBarry Smith     }
4780fdc7489SMatthew Knepley 
4799566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj));
4809318fe57SMatthew G. Knepley     if (modelObj) {
4810fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
4820fdc7489SMatthew Knepley       DMLabel        bodyLabel;
4830fdc7489SMatthew Knepley       PetscInt       cStart, cEnd, c, eStart, eEnd, fStart, fEnd;
484c1cad2e7SMatthew G. Knepley       PetscBool      islite = PETSC_FALSE;
4850fdc7489SMatthew Knepley       ego           *bodies;
4860fdc7489SMatthew Knepley       ego            model, geom;
4870fdc7489SMatthew Knepley       int            Nb, oclass, mtype, *senses;
4880fdc7489SMatthew Knepley 
4890fdc7489SMatthew Knepley       /* Get Attached EGADS Model from Original DMPlex */
4909566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj));
491c1cad2e7SMatthew G. Knepley       if (modelObj) {
4929566063dSJacob Faibussowitsch         PetscCall(PetscContainerGetPointer(modelObj, (void **) &model));
4939566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
4940fdc7489SMatthew Knepley         /* Transfer EGADS Model to Volumetric Mesh */
4959566063dSJacob Faibussowitsch         PetscCall(PetscObjectCompose((PetscObject) *dmRefined, "EGADS Model", (PetscObject) modelObj));
496c1cad2e7SMatthew G. Knepley       } else {
4979566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject) dm, "EGADSLite Model", (PetscObject *) &modelObj));
498c1cad2e7SMatthew G. Knepley         if (modelObj) {
4999566063dSJacob Faibussowitsch           PetscCall(PetscContainerGetPointer(modelObj, (void **) &model));
5009566063dSJacob Faibussowitsch           PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
501c1cad2e7SMatthew G. Knepley           /* Transfer EGADS Model to Volumetric Mesh */
5029566063dSJacob Faibussowitsch           PetscCall(PetscObjectCompose((PetscObject) *dmRefined, "EGADSLite Model", (PetscObject) modelObj));
503c1cad2e7SMatthew G. Knepley           islite = PETSC_TRUE;
504c1cad2e7SMatthew G. Knepley         }
505c1cad2e7SMatthew G. Knepley       }
506c1cad2e7SMatthew G. Knepley       if (!modelObj) goto skip_egads;
5070fdc7489SMatthew Knepley 
5080fdc7489SMatthew Knepley       /* Set Cell Labels */
5099566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dmRefined, "EGADS Body ID", &bodyLabel));
5109566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(*dmRefined, 0, &cStart, &cEnd));
5119566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(*dmRefined, 1, &fStart, &fEnd));
5129566063dSJacob Faibussowitsch       PetscCall(DMPlexGetDepthStratum(*dmRefined, 1, &eStart, &eEnd));
5130fdc7489SMatthew Knepley 
5140fdc7489SMatthew Knepley       for (c = cStart; c < cEnd; ++c) {
5150fdc7489SMatthew Knepley         PetscReal centroid[3] = {0., 0., 0.};
5160fdc7489SMatthew Knepley         PetscInt  b;
5170fdc7489SMatthew Knepley 
5180fdc7489SMatthew Knepley         /* Deterimine what body the cell's centroid is located in */
5190fdc7489SMatthew Knepley         if (!interpolate) {
5200fdc7489SMatthew Knepley           PetscSection   coordSection;
5210fdc7489SMatthew Knepley           Vec            coordinates;
5220fdc7489SMatthew Knepley           PetscScalar   *coords = NULL;
5230fdc7489SMatthew Knepley           PetscInt       coordSize, s, d;
5240fdc7489SMatthew Knepley 
5259566063dSJacob Faibussowitsch           PetscCall(DMGetCoordinatesLocal(*dmRefined, &coordinates));
5269566063dSJacob Faibussowitsch           PetscCall(DMGetCoordinateSection(*dmRefined, &coordSection));
5279566063dSJacob Faibussowitsch           PetscCall(DMPlexVecGetClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords));
5280fdc7489SMatthew Knepley           for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d];
5299566063dSJacob Faibussowitsch           PetscCall(DMPlexVecRestoreClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords));
5301baa6e33SBarry Smith         } else PetscCall(DMPlexComputeCellGeometryFVM(*dmRefined, c, NULL, centroid, NULL));
5310fdc7489SMatthew Knepley         for (b = 0; b < Nb; ++b) {
532c1cad2e7SMatthew G. Knepley           if (islite) {if (EGlite_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
533c1cad2e7SMatthew G. Knepley           else        {if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
5340fdc7489SMatthew Knepley         }
5350fdc7489SMatthew Knepley         if (b < Nb) {
5360fdc7489SMatthew Knepley           PetscInt   cval = b, eVal, fVal;
5370fdc7489SMatthew Knepley           PetscInt *closure = NULL, Ncl, cl;
5380fdc7489SMatthew Knepley 
5399566063dSJacob Faibussowitsch           PetscCall(DMLabelSetValue(bodyLabel, c, cval));
5409566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure));
5410fdc7489SMatthew Knepley           for (cl = 0; cl < Ncl; cl += 2) {
5420fdc7489SMatthew Knepley             const PetscInt p = closure[cl];
5430fdc7489SMatthew Knepley 
5440fdc7489SMatthew Knepley             if (p >= eStart && p < eEnd) {
5459566063dSJacob Faibussowitsch               PetscCall(DMLabelGetValue(bodyLabel, p, &eVal));
5469566063dSJacob Faibussowitsch               if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
5470fdc7489SMatthew Knepley             }
5480fdc7489SMatthew Knepley             if (p >= fStart && p < fEnd) {
5499566063dSJacob Faibussowitsch               PetscCall(DMLabelGetValue(bodyLabel, p, &fVal));
5509566063dSJacob Faibussowitsch               if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval));
5510fdc7489SMatthew Knepley             }
5520fdc7489SMatthew Knepley           }
5539566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure));
5540fdc7489SMatthew Knepley         }
5550fdc7489SMatthew Knepley       }
556c1cad2e7SMatthew G. Knepley skip_egads: ;
5570fdc7489SMatthew Knepley #endif
5589318fe57SMatthew G. Knepley     }
5599566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE));
5603a074057SBarry Smith   }
5613a074057SBarry Smith   PetscFunctionReturn(0);
5623a074057SBarry Smith }
563