xref: /petsc/src/dm/impls/plex/generators/tetgen/tetgenerate.cxx (revision c1cad2e7c584d112105df06cee8bea5f275d3d5d)
13a074057SBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
23a074057SBarry Smith 
30fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
40fdc7489SMatthew Knepley #include <egads.h>
5*c1cad2e7SMatthew G. Knepley /* Need to make EGADSLite header compatible */
6*c1cad2e7SMatthew G. Knepley extern "C" int EGlite_getTopology(const ego, ego *, int *, int *, double *, int *, ego **, int **);
7*c1cad2e7SMatthew 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;
350fdc7489SMatthew Knepley   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f;
360fdc7489SMatthew Knepley   DMPlexInterpolatedFlag isInterpolated;
373a074057SBarry Smith   PetscMPIInt            rank;
383a074057SBarry Smith   PetscErrorCode         ierr;
393a074057SBarry Smith 
403a074057SBarry Smith   PetscFunctionBegin;
413a074057SBarry Smith   ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr);
42ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
430fdc7489SMatthew Knepley   ierr = DMPlexIsInterpolatedCollective(boundary, &isInterpolated);CHKERRQ(ierr);
440fdc7489SMatthew Knepley   ierr = DMUniversalLabelCreate(boundary, &universal);CHKERRQ(ierr);
453a074057SBarry Smith 
460fdc7489SMatthew Knepley   ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr);
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 
563a074057SBarry Smith     ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr);
573a074057SBarry Smith     ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr);
580fdc7489SMatthew Knepley     ierr = VecGetArrayRead(coordinates, &array);CHKERRQ(ierr);
593a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
603a074057SBarry Smith       const PetscInt idx = v - vStart;
610fdc7489SMatthew Knepley       PetscInt       off, d, val;
623a074057SBarry Smith 
633a074057SBarry Smith       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
643a074057SBarry Smith       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
650fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, v, &val);CHKERRQ(ierr);
663a074057SBarry Smith       in.pointmarkerlist[idx] = (int) val;
673a074057SBarry Smith     }
680fdc7489SMatthew Knepley     ierr = VecRestoreArrayRead(coordinates, &array);CHKERRQ(ierr);
693a074057SBarry Smith   }
703a074057SBarry Smith 
710fdc7489SMatthew Knepley   ierr = DMPlexGetHeightStratum(boundary, 1, &eStart, &eEnd);CHKERRQ(ierr);
720fdc7489SMatthew Knepley   in.numberofedges = eEnd - eStart;
730fdc7489SMatthew Knepley   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberofedges > 0) {
740fdc7489SMatthew Knepley     in.edgelist       = new int[in.numberofedges * 2];
750fdc7489SMatthew Knepley     in.edgemarkerlist = new int[in.numberofedges];
760fdc7489SMatthew Knepley     for (e = eStart; e < eEnd; ++e) {
770fdc7489SMatthew Knepley       const PetscInt  idx = e - eStart;
780fdc7489SMatthew Knepley       const PetscInt *cone;
790fdc7489SMatthew Knepley       PetscInt        coneSize, val;
800fdc7489SMatthew Knepley 
810fdc7489SMatthew Knepley       ierr = DMPlexGetConeSize(boundary, e, &coneSize);CHKERRQ(ierr);
820fdc7489SMatthew Knepley       ierr = DMPlexGetCone(boundary, e, &cone);CHKERRQ(ierr);
830fdc7489SMatthew Knepley       in.edgelist[idx*2]     = cone[0] - vStart;
840fdc7489SMatthew Knepley       in.edgelist[idx*2 + 1] = cone[1] - vStart;
850fdc7489SMatthew Knepley 
860fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, e, &val);CHKERRQ(ierr);
870fdc7489SMatthew Knepley       in.edgemarkerlist[idx] = (int) val;
880fdc7489SMatthew Knepley     }
890fdc7489SMatthew Knepley   }
900fdc7489SMatthew Knepley 
910fdc7489SMatthew Knepley   ierr = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr);
923a074057SBarry Smith   in.numberoffacets = fEnd - fStart;
933a074057SBarry Smith   if (in.numberoffacets > 0) {
943a074057SBarry Smith     in.facetlist       = new tetgenio::facet[in.numberoffacets];
953a074057SBarry Smith     in.facetmarkerlist = new int[in.numberoffacets];
963a074057SBarry Smith     for (f = fStart; f < fEnd; ++f) {
973a074057SBarry Smith       const PetscInt idx    = f - fStart;
980fdc7489SMatthew Knepley       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, val = -1;
993a074057SBarry Smith 
1003a074057SBarry Smith       in.facetlist[idx].numberofpolygons = 1;
1013a074057SBarry Smith       in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
1023a074057SBarry Smith       in.facetlist[idx].numberofholes    = 0;
1033a074057SBarry Smith       in.facetlist[idx].holelist         = NULL;
1043a074057SBarry Smith 
1053a074057SBarry Smith       ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
1063a074057SBarry Smith       for (p = 0; p < numPoints*2; p += 2) {
1073a074057SBarry Smith         const PetscInt point = points[p];
1083a074057SBarry Smith         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
1093a074057SBarry Smith       }
1103a074057SBarry Smith 
1113a074057SBarry Smith       tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
1123a074057SBarry Smith       poly->numberofvertices = numVertices;
1133a074057SBarry Smith       poly->vertexlist       = new int[poly->numberofvertices];
1143a074057SBarry Smith       for (v = 0; v < numVertices; ++v) {
1153a074057SBarry Smith         const PetscInt vIdx = points[v] - vStart;
1163a074057SBarry Smith         poly->vertexlist[v] = vIdx;
1173a074057SBarry Smith       }
1180fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, f, &val);CHKERRQ(ierr);
1193a074057SBarry Smith       in.facetmarkerlist[idx] = (int) val;
1203a074057SBarry Smith       ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
1213a074057SBarry Smith     }
1223a074057SBarry Smith   }
1233a074057SBarry Smith   if (!rank) {
1240fdc7489SMatthew Knepley     DM_Plex *mesh = (DM_Plex *) boundary->data;
1253a074057SBarry Smith     char     args[32];
1263a074057SBarry Smith 
1273a074057SBarry Smith     /* Take away 'Q' for verbose output */
1280fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
1290fdc7489SMatthew Knepley     ierr = PetscStrcpy(args, "pqezQY");CHKERRQ(ierr);
1300fdc7489SMatthew Knepley #else
1313a074057SBarry Smith     ierr = PetscStrcpy(args, "pqezQ");CHKERRQ(ierr);
1320fdc7489SMatthew Knepley #endif
1333a074057SBarry Smith     if (mesh->tetgenOpts) {::tetrahedralize(mesh->tetgenOpts, &in, &out);}
1343a074057SBarry Smith     else                  {::tetrahedralize(args, &in, &out);}
1353a074057SBarry Smith   }
1363a074057SBarry Smith   {
1373a074057SBarry Smith     const PetscInt   numCorners  = 4;
1383a074057SBarry Smith     const PetscInt   numCells    = out.numberoftetrahedra;
1393a074057SBarry Smith     const PetscInt   numVertices = out.numberofpoints;
140a4a685f2SJacob Faibussowitsch     PetscReal        *meshCoords = NULL;
141a4a685f2SJacob Faibussowitsch     PetscInt         *cells      = NULL;
142a4a685f2SJacob Faibussowitsch 
143a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) == sizeof (out.pointlist[0])) {
144a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *) out.pointlist;
145a4a685f2SJacob Faibussowitsch     } else {
146a4a685f2SJacob Faibussowitsch       PetscInt i;
147a4a685f2SJacob Faibussowitsch 
148a4a685f2SJacob Faibussowitsch       meshCoords = new PetscReal[dim * numVertices];
1490fdc7489SMatthew Knepley       for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out.pointlist[i];
150a4a685f2SJacob Faibussowitsch     }
151a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) == sizeof (out.tetrahedronlist[0])) {
152a4a685f2SJacob Faibussowitsch       cells = (PetscInt *) out.tetrahedronlist;
153a4a685f2SJacob Faibussowitsch     } else {
154a4a685f2SJacob Faibussowitsch       PetscInt i;
155a4a685f2SJacob Faibussowitsch 
156a4a685f2SJacob Faibussowitsch       cells = new PetscInt[numCells * numCorners];
1570fdc7489SMatthew Knepley       for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt) out.tetrahedronlist[i];
158a4a685f2SJacob Faibussowitsch     }
1593a074057SBarry Smith 
16096ca5757SLisandro Dalcin     ierr = DMPlexInvertCells_Tetgen(numCells, numCorners, cells);CHKERRQ(ierr);
161a4a685f2SJacob Faibussowitsch     ierr = DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr);
1620fdc7489SMatthew Knepley 
1633a074057SBarry Smith     /* Set labels */
1640fdc7489SMatthew Knepley     ierr = DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dm);CHKERRQ(ierr);
1653a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
1663a074057SBarry Smith       if (out.pointmarkerlist[v]) {
1670fdc7489SMatthew Knepley         ierr = DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);
1683a074057SBarry Smith       }
1693a074057SBarry Smith     }
1703a074057SBarry Smith     if (interpolate) {
1713a074057SBarry Smith       PetscInt e;
1723a074057SBarry Smith 
1733a074057SBarry Smith       for (e = 0; e < out.numberofedges; e++) {
1743a074057SBarry Smith         if (out.edgemarkerlist[e]) {
1753a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
1763a074057SBarry Smith           const PetscInt *edges;
1773a074057SBarry Smith           PetscInt        numEdges;
1783a074057SBarry Smith 
1793a074057SBarry Smith           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1803a074057SBarry Smith           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
1810fdc7489SMatthew Knepley           ierr = DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);
1823a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1833a074057SBarry Smith         }
1843a074057SBarry Smith       }
1853a074057SBarry Smith       for (f = 0; f < out.numberoftrifaces; f++) {
1863a074057SBarry Smith         if (out.trifacemarkerlist[f]) {
1873a074057SBarry Smith           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
1883a074057SBarry Smith           const PetscInt *faces;
1893a074057SBarry Smith           PetscInt        numFaces;
1903a074057SBarry Smith 
1913a074057SBarry Smith           ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1923a074057SBarry Smith           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
1930fdc7489SMatthew Knepley           ierr = DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);
1943a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1953a074057SBarry Smith         }
1963a074057SBarry Smith       }
1973a074057SBarry Smith     }
1980fdc7489SMatthew Knepley 
1999318fe57SMatthew G. Knepley     ierr = PetscObjectQuery((PetscObject) boundary, "EGADS Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
2009318fe57SMatthew G. Knepley     if (modelObj) {
2010fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
2020fdc7489SMatthew Knepley       DMLabel        bodyLabel;
2030fdc7489SMatthew Knepley       PetscInt       cStart, cEnd, c, eStart, eEnd, fStart, fEnd;
204*c1cad2e7SMatthew G. Knepley       PetscBool      islite = PETSC_FALSE;
2050fdc7489SMatthew Knepley       ego           *bodies;
2060fdc7489SMatthew Knepley       ego            model, geom;
2070fdc7489SMatthew Knepley       int            Nb, oclass, mtype, *senses;
2080fdc7489SMatthew Knepley 
2090fdc7489SMatthew Knepley       /* Get Attached EGADS Model from Original DMPlex */
210*c1cad2e7SMatthew G. Knepley       ierr = PetscObjectQuery((PetscObject) boundary, "EGADS Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
211*c1cad2e7SMatthew G. Knepley       if (modelObj) {
2120fdc7489SMatthew Knepley         ierr = PetscContainerGetPointer(modelObj, (void **) &model);CHKERRQ(ierr);
2130fdc7489SMatthew Knepley         ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
2140fdc7489SMatthew Knepley         /* Transfer EGADS Model to Volumetric Mesh */
2150fdc7489SMatthew Knepley         ierr = PetscObjectCompose((PetscObject) *dm, "EGADS Model", (PetscObject) modelObj);CHKERRQ(ierr);
216*c1cad2e7SMatthew G. Knepley       } else {
217*c1cad2e7SMatthew G. Knepley         ierr = PetscObjectQuery((PetscObject) boundary, "EGADSLite Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
218*c1cad2e7SMatthew G. Knepley         if (modelObj) {
219*c1cad2e7SMatthew G. Knepley           ierr = PetscContainerGetPointer(modelObj, (void **) &model);CHKERRQ(ierr);
220*c1cad2e7SMatthew G. Knepley           ierr = EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
221*c1cad2e7SMatthew G. Knepley           /* Transfer EGADS Model to Volumetric Mesh */
222*c1cad2e7SMatthew G. Knepley           ierr = PetscObjectCompose((PetscObject) *dm, "EGADSLite Model", (PetscObject) modelObj);CHKERRQ(ierr);
223*c1cad2e7SMatthew G. Knepley           islite = PETSC_TRUE;
224*c1cad2e7SMatthew G. Knepley         }
225*c1cad2e7SMatthew G. Knepley       }
226*c1cad2e7SMatthew G. Knepley       if (!modelObj) goto skip_egads;
2270fdc7489SMatthew Knepley 
2280fdc7489SMatthew Knepley       /* Set Cell Labels */
2290fdc7489SMatthew Knepley       ierr = DMGetLabel(*dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
2300fdc7489SMatthew Knepley       ierr = DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2310fdc7489SMatthew Knepley       ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
2320fdc7489SMatthew Knepley       ierr = DMPlexGetDepthStratum(*dm, 1, &eStart, &eEnd);CHKERRQ(ierr);
2330fdc7489SMatthew Knepley 
2340fdc7489SMatthew Knepley       for (c = cStart; c < cEnd; ++c) {
2350fdc7489SMatthew Knepley         PetscReal centroid[3] = {0., 0., 0.};
2360fdc7489SMatthew Knepley         PetscInt  b;
2370fdc7489SMatthew Knepley 
2380fdc7489SMatthew Knepley         /* Deterimine what body the cell's centroid is located in */
2390fdc7489SMatthew Knepley         if (!interpolate) {
2400fdc7489SMatthew Knepley           PetscSection   coordSection;
2410fdc7489SMatthew Knepley           Vec            coordinates;
2420fdc7489SMatthew Knepley           PetscScalar   *coords = NULL;
2430fdc7489SMatthew Knepley           PetscInt       coordSize, s, d;
2440fdc7489SMatthew Knepley 
2450fdc7489SMatthew Knepley           ierr = DMGetCoordinatesLocal(*dm, &coordinates);CHKERRQ(ierr);
2460fdc7489SMatthew Knepley           ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
2470fdc7489SMatthew Knepley           ierr = DMPlexVecGetClosure(*dm, coordSection, coordinates, c, &coordSize, &coords);CHKERRQ(ierr);
2480fdc7489SMatthew Knepley           for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d];
2490fdc7489SMatthew Knepley           ierr = DMPlexVecRestoreClosure(*dm, coordSection, coordinates, c, &coordSize, &coords);CHKERRQ(ierr);
2500fdc7489SMatthew Knepley         } else {
2510fdc7489SMatthew Knepley           ierr = DMPlexComputeCellGeometryFVM(*dm, c, NULL, centroid, NULL);CHKERRQ(ierr);
2520fdc7489SMatthew Knepley         }
2530fdc7489SMatthew Knepley         for (b = 0; b < Nb; ++b) {
254*c1cad2e7SMatthew G. Knepley           if (islite) {if (EGlite_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
255*c1cad2e7SMatthew G. Knepley           else        {if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
2560fdc7489SMatthew Knepley         }
2570fdc7489SMatthew Knepley         if (b < Nb) {
2580fdc7489SMatthew Knepley           PetscInt   cval = b, eVal, fVal;
2590fdc7489SMatthew Knepley           PetscInt *closure = NULL, Ncl, cl;
2600fdc7489SMatthew Knepley 
2610fdc7489SMatthew Knepley           ierr = DMLabelSetValue(bodyLabel, c, cval);CHKERRQ(ierr);
2620fdc7489SMatthew Knepley           ierr = DMPlexGetTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
2630fdc7489SMatthew Knepley           for (cl = 0; cl < Ncl; cl += 2) {
2640fdc7489SMatthew Knepley             const PetscInt p = closure[cl];
2650fdc7489SMatthew Knepley 
2660fdc7489SMatthew Knepley             if (p >= eStart && p < eEnd) {
2670fdc7489SMatthew Knepley               ierr = DMLabelGetValue(bodyLabel, p, &eVal);CHKERRQ(ierr);
2680fdc7489SMatthew Knepley               if (eVal < 0) {ierr = DMLabelSetValue(bodyLabel, p, cval);CHKERRQ(ierr);}
2690fdc7489SMatthew Knepley             }
2700fdc7489SMatthew Knepley             if (p >= fStart && p < fEnd) {
2710fdc7489SMatthew Knepley               ierr = DMLabelGetValue(bodyLabel, p, &fVal);CHKERRQ(ierr);
2720fdc7489SMatthew Knepley               if (fVal < 0) {ierr = DMLabelSetValue(bodyLabel, p, cval);CHKERRQ(ierr);}
2730fdc7489SMatthew Knepley             }
2740fdc7489SMatthew Knepley           }
2750fdc7489SMatthew Knepley           ierr = DMPlexRestoreTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
2760fdc7489SMatthew Knepley         }
2770fdc7489SMatthew Knepley       }
278*c1cad2e7SMatthew G. Knepley skip_egads: ;
2790fdc7489SMatthew Knepley #endif
2809318fe57SMatthew G. Knepley     }
2813a074057SBarry Smith     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
2823a074057SBarry Smith   }
283*c1cad2e7SMatthew G. Knepley   ierr = DMUniversalLabelDestroy(&universal);CHKERRQ(ierr);
2843a074057SBarry Smith   PetscFunctionReturn(0);
2853a074057SBarry Smith }
2863a074057SBarry Smith 
2873a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexRefine_Tetgen(DM dm, double *maxVolumes, DM *dmRefined)
2883a074057SBarry Smith {
2893a074057SBarry Smith   MPI_Comm               comm;
2903a074057SBarry Smith   const PetscInt         dim = 3;
2913a074057SBarry Smith   ::tetgenio             in;
2923a074057SBarry Smith   ::tetgenio             out;
2939318fe57SMatthew G. Knepley   PetscContainer         modelObj;
2940fdc7489SMatthew Knepley   DMUniversalLabel       universal;
2950fdc7489SMatthew Knepley   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, cStart, cEnd, c;
2960fdc7489SMatthew Knepley   DMPlexInterpolatedFlag isInterpolated;
2973a074057SBarry Smith   PetscMPIInt            rank;
2983a074057SBarry Smith   PetscErrorCode         ierr;
2993a074057SBarry Smith 
3003a074057SBarry Smith   PetscFunctionBegin;
3013a074057SBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
302ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
3030fdc7489SMatthew Knepley   ierr = DMPlexIsInterpolatedCollective(dm, &isInterpolated);CHKERRQ(ierr);
3040fdc7489SMatthew Knepley   ierr = DMUniversalLabelCreate(dm, &universal);CHKERRQ(ierr);
3053a074057SBarry Smith 
3060fdc7489SMatthew Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
3073a074057SBarry Smith   in.numberofpoints = vEnd - vStart;
3083a074057SBarry Smith   if (in.numberofpoints > 0) {
3093a074057SBarry Smith     PetscSection coordSection;
3103a074057SBarry Smith     Vec          coordinates;
3113a074057SBarry Smith     PetscScalar *array;
3123a074057SBarry Smith 
3133a074057SBarry Smith     in.pointlist       = new double[in.numberofpoints*dim];
3143a074057SBarry Smith     in.pointmarkerlist = new int[in.numberofpoints];
3153a074057SBarry Smith 
3163a074057SBarry Smith     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
3173a074057SBarry Smith     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
3183a074057SBarry Smith     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
3193a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
3203a074057SBarry Smith       const PetscInt idx = v - vStart;
3210fdc7489SMatthew Knepley       PetscInt       off, d, val;
3223a074057SBarry Smith 
3233a074057SBarry Smith       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
3243a074057SBarry Smith       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
3250fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, v, &val);CHKERRQ(ierr);
3263a074057SBarry Smith       in.pointmarkerlist[idx] = (int) val;
3273a074057SBarry Smith     }
3283a074057SBarry Smith     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
3293a074057SBarry Smith   }
3303a074057SBarry Smith 
3310fdc7489SMatthew Knepley   ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);
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 
3410fdc7489SMatthew Knepley       ierr = DMPlexGetConeSize(dm, e, &coneSize);CHKERRQ(ierr);
3420fdc7489SMatthew Knepley       ierr = DMPlexGetCone(dm, e, &cone);CHKERRQ(ierr);
3430fdc7489SMatthew Knepley       in.edgelist[idx*2]     = cone[0] - vStart;
3440fdc7489SMatthew Knepley       in.edgelist[idx*2 + 1] = cone[1] - vStart;
3450fdc7489SMatthew Knepley 
3460fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, e, &val);CHKERRQ(ierr);
3470fdc7489SMatthew Knepley       in.edgemarkerlist[idx] = (int) val;
3480fdc7489SMatthew Knepley     }
3490fdc7489SMatthew Knepley   }
3500fdc7489SMatthew Knepley 
3510fdc7489SMatthew Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
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 
3650fdc7489SMatthew Knepley       ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
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 
3790fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, f, &val);CHKERRQ(ierr);
3800fdc7489SMatthew Knepley       in.facetmarkerlist[idx] = (int) val;
3810fdc7489SMatthew Knepley 
3820fdc7489SMatthew Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
3830fdc7489SMatthew Knepley     }
3840fdc7489SMatthew Knepley   }
3850fdc7489SMatthew Knepley 
3860fdc7489SMatthew Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
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 
3973a074057SBarry Smith       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
3983a074057SBarry Smith       if ((closureSize != 5) && (closureSize != 15)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %D vertices in closure", closureSize);
3990fdc7489SMatthew Knepley       for (v = 0; v < 4; ++v) in.tetrahedronlist[idx*in.numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
4003a074057SBarry Smith       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
4013a074057SBarry Smith     }
4023a074057SBarry Smith   }
4030fdc7489SMatthew Knepley 
4043a074057SBarry Smith   if (!rank) {
4053a074057SBarry Smith     char args[32];
4063a074057SBarry Smith 
4073a074057SBarry Smith     /* Take away 'Q' for verbose output */
4083a074057SBarry Smith     ierr = PetscStrcpy(args, "qezQra");CHKERRQ(ierr);
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 
43896ca5757SLisandro Dalcin     ierr = DMPlexInvertCells_Tetgen(numCells, numCorners, cells);CHKERRQ(ierr);
439a4a685f2SJacob Faibussowitsch     ierr = DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
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 */
4440fdc7489SMatthew Knepley     ierr = DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dmRefined);CHKERRQ(ierr);
4453a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
4463a074057SBarry Smith       if (out.pointmarkerlist[v]) {
4470fdc7489SMatthew Knepley         ierr = DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);
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 
4593a074057SBarry Smith           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
4603a074057SBarry Smith           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
4610fdc7489SMatthew Knepley           ierr = DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);
4623a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
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 
4713a074057SBarry Smith           ierr = DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
4723a074057SBarry Smith           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
4730fdc7489SMatthew Knepley           ierr = DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);
4743a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
4753a074057SBarry Smith         }
4763a074057SBarry Smith       }
4773a074057SBarry Smith     }
4780fdc7489SMatthew Knepley 
4799318fe57SMatthew G. Knepley     ierr = PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
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;
484*c1cad2e7SMatthew 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 */
490*c1cad2e7SMatthew G. Knepley       ierr = PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
491*c1cad2e7SMatthew G. Knepley       if (modelObj) {
4920fdc7489SMatthew Knepley         ierr = PetscContainerGetPointer(modelObj, (void **) &model);CHKERRQ(ierr);
4930fdc7489SMatthew Knepley         ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
4940fdc7489SMatthew Knepley         /* Transfer EGADS Model to Volumetric Mesh */
4950fdc7489SMatthew Knepley         ierr = PetscObjectCompose((PetscObject) *dmRefined, "EGADS Model", (PetscObject) modelObj);CHKERRQ(ierr);
496*c1cad2e7SMatthew G. Knepley       } else {
497*c1cad2e7SMatthew G. Knepley         ierr = PetscObjectQuery((PetscObject) dm, "EGADSLite Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
498*c1cad2e7SMatthew G. Knepley         if (modelObj) {
499*c1cad2e7SMatthew G. Knepley           ierr = PetscContainerGetPointer(modelObj, (void **) &model);CHKERRQ(ierr);
500*c1cad2e7SMatthew G. Knepley           ierr = EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
501*c1cad2e7SMatthew G. Knepley           /* Transfer EGADS Model to Volumetric Mesh */
502*c1cad2e7SMatthew G. Knepley           ierr = PetscObjectCompose((PetscObject) *dmRefined, "EGADSLite Model", (PetscObject) modelObj);CHKERRQ(ierr);
503*c1cad2e7SMatthew G. Knepley           islite = PETSC_TRUE;
504*c1cad2e7SMatthew G. Knepley         }
505*c1cad2e7SMatthew G. Knepley       }
506*c1cad2e7SMatthew G. Knepley       if (!modelObj) goto skip_egads;
5070fdc7489SMatthew Knepley 
5080fdc7489SMatthew Knepley       /* Set Cell Labels */
5090fdc7489SMatthew Knepley       ierr = DMGetLabel(*dmRefined, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
5100fdc7489SMatthew Knepley       ierr = DMPlexGetHeightStratum(*dmRefined, 0, &cStart, &cEnd);CHKERRQ(ierr);
5110fdc7489SMatthew Knepley       ierr = DMPlexGetHeightStratum(*dmRefined, 1, &fStart, &fEnd);CHKERRQ(ierr);
5120fdc7489SMatthew Knepley       ierr = DMPlexGetDepthStratum(*dmRefined, 1, &eStart, &eEnd);CHKERRQ(ierr);
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 
5250fdc7489SMatthew Knepley           ierr = DMGetCoordinatesLocal(*dmRefined, &coordinates);CHKERRQ(ierr);
5260fdc7489SMatthew Knepley           ierr = DMGetCoordinateSection(*dmRefined, &coordSection);CHKERRQ(ierr);
5270fdc7489SMatthew Knepley           ierr = DMPlexVecGetClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords);CHKERRQ(ierr);
5280fdc7489SMatthew Knepley           for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d];
5290fdc7489SMatthew Knepley           ierr = DMPlexVecRestoreClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords);CHKERRQ(ierr);
5300fdc7489SMatthew Knepley         } else {
5310fdc7489SMatthew Knepley           ierr = DMPlexComputeCellGeometryFVM(*dmRefined, c, NULL, centroid, NULL);CHKERRQ(ierr);
5320fdc7489SMatthew Knepley         }
5330fdc7489SMatthew Knepley         for (b = 0; b < Nb; ++b) {
534*c1cad2e7SMatthew G. Knepley           if (islite) {if (EGlite_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
535*c1cad2e7SMatthew G. Knepley           else        {if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;}
5360fdc7489SMatthew Knepley         }
5370fdc7489SMatthew Knepley         if (b < Nb) {
5380fdc7489SMatthew Knepley           PetscInt   cval = b, eVal, fVal;
5390fdc7489SMatthew Knepley           PetscInt *closure = NULL, Ncl, cl;
5400fdc7489SMatthew Knepley 
5410fdc7489SMatthew Knepley           ierr = DMLabelSetValue(bodyLabel, c, cval);CHKERRQ(ierr);
5420fdc7489SMatthew Knepley           ierr = DMPlexGetTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
5430fdc7489SMatthew Knepley           for (cl = 0; cl < Ncl; cl += 2) {
5440fdc7489SMatthew Knepley             const PetscInt p = closure[cl];
5450fdc7489SMatthew Knepley 
5460fdc7489SMatthew Knepley             if (p >= eStart && p < eEnd) {
5470fdc7489SMatthew Knepley               ierr = DMLabelGetValue(bodyLabel, p, &eVal);CHKERRQ(ierr);
5480fdc7489SMatthew Knepley               if (eVal < 0) {ierr = DMLabelSetValue(bodyLabel, p, cval);CHKERRQ(ierr);}
5490fdc7489SMatthew Knepley             }
5500fdc7489SMatthew Knepley             if (p >= fStart && p < fEnd) {
5510fdc7489SMatthew Knepley               ierr = DMLabelGetValue(bodyLabel, p, &fVal);CHKERRQ(ierr);
5520fdc7489SMatthew Knepley               if (fVal < 0) {ierr = DMLabelSetValue(bodyLabel, p, cval);CHKERRQ(ierr);}
5530fdc7489SMatthew Knepley             }
5540fdc7489SMatthew Knepley           }
5550fdc7489SMatthew Knepley           ierr = DMPlexRestoreTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
5560fdc7489SMatthew Knepley         }
5570fdc7489SMatthew Knepley       }
558*c1cad2e7SMatthew G. Knepley skip_egads: ;
5590fdc7489SMatthew Knepley #endif
5609318fe57SMatthew G. Knepley     }
5613a074057SBarry Smith     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
5623a074057SBarry Smith   }
5633a074057SBarry Smith   PetscFunctionReturn(0);
5643a074057SBarry Smith }
565