xref: /petsc/src/dm/impls/plex/generators/ctetgen/ctetgenerate.c (revision 0fdc74893d08b7f731a3f927b6584bf9ea7d2444)
13a074057SBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
23a074057SBarry Smith 
3*0fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
4*0fdc7489SMatthew Knepley #include <egads.h>
5*0fdc7489SMatthew Knepley #endif
6*0fdc7489SMatthew Knepley 
73a074057SBarry Smith #include <ctetgen.h>
83a074057SBarry Smith 
93a074057SBarry Smith /* This is to fix the tetrahedron orientation from TetGen */
10a4a685f2SJacob Faibussowitsch static PetscErrorCode DMPlexInvertCells_CTetgen(PetscInt numCells, PetscInt numCorners, PetscInt cells[])
113a074057SBarry Smith {
123a074057SBarry Smith   PetscInt bound = numCells*numCorners, coff;
133a074057SBarry Smith 
143a074057SBarry Smith   PetscFunctionBegin;
15a4a685f2SJacob Faibussowitsch #define SWAP(a,b) do { PetscInt tmp = (a); (a) = (b); (b) = tmp; } while (0)
1696ca5757SLisandro Dalcin   for (coff = 0; coff < bound; coff += numCorners) SWAP(cells[coff],cells[coff+1]);
1796ca5757SLisandro Dalcin #undef SWAP
183a074057SBarry Smith   PetscFunctionReturn(0);
193a074057SBarry Smith }
203a074057SBarry Smith 
213a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerate_CTetgen(DM boundary, PetscBool interpolate, DM *dm)
223a074057SBarry Smith {
233a074057SBarry Smith   MPI_Comm               comm;
243a074057SBarry Smith   const PetscInt         dim = 3;
253a074057SBarry Smith   PLC                   *in, *out;
26*0fdc7489SMatthew Knepley   DMUniversalLabel       universal;
27*0fdc7489SMatthew Knepley   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, verbose = 0;
28*0fdc7489SMatthew Knepley   DMPlexInterpolatedFlag isInterpolated;
293a074057SBarry Smith   PetscMPIInt            rank;
303a074057SBarry Smith   PetscErrorCode         ierr;
313a074057SBarry Smith 
323a074057SBarry Smith   PetscFunctionBegin;
333a074057SBarry Smith   ierr = PetscOptionsGetInt(NULL,((PetscObject) boundary)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr);
34*0fdc7489SMatthew Knepley   ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr);
353a074057SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
36*0fdc7489SMatthew Knepley   ierr = DMPlexIsInterpolatedCollective(boundary, &isInterpolated);CHKERRQ(ierr);
37*0fdc7489SMatthew Knepley   ierr = DMUniversalLabelCreate(boundary, &universal);CHKERRQ(ierr);
38*0fdc7489SMatthew Knepley 
393a074057SBarry Smith   ierr = PLCCreate(&in);CHKERRQ(ierr);
403a074057SBarry Smith   ierr = PLCCreate(&out);CHKERRQ(ierr);
413a074057SBarry Smith 
42*0fdc7489SMatthew Knepley   ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr);
433a074057SBarry Smith   in->numberofpoints = vEnd - vStart;
443a074057SBarry Smith   if (in->numberofpoints > 0) {
453a074057SBarry Smith     PetscSection       coordSection;
463a074057SBarry Smith     Vec                coordinates;
47*0fdc7489SMatthew Knepley     const PetscScalar *array;
483a074057SBarry Smith 
493a074057SBarry Smith     ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr);
503a074057SBarry Smith     ierr = PetscMalloc1(in->numberofpoints,     &in->pointmarkerlist);CHKERRQ(ierr);
513a074057SBarry Smith     ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr);
523a074057SBarry Smith     ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr);
53*0fdc7489SMatthew Knepley     ierr = VecGetArrayRead(coordinates, &array);CHKERRQ(ierr);
543a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
553a074057SBarry Smith       const PetscInt idx = v - vStart;
56*0fdc7489SMatthew Knepley       PetscInt       off, d, m;
573a074057SBarry Smith 
583a074057SBarry Smith       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
59*0fdc7489SMatthew Knepley       for (d = 0; d < dim; ++d) in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
60*0fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, v, &m);CHKERRQ(ierr);
613a074057SBarry Smith       in->pointmarkerlist[idx] = (int) m;
623a074057SBarry Smith     }
63*0fdc7489SMatthew Knepley     ierr = VecRestoreArrayRead(coordinates, &array);CHKERRQ(ierr);
643a074057SBarry Smith   }
653a074057SBarry Smith 
66*0fdc7489SMatthew Knepley   ierr = DMPlexGetHeightStratum(boundary, 1, &eStart, &eEnd);CHKERRQ(ierr);
67*0fdc7489SMatthew Knepley   in->numberofedges = eEnd - eStart;
68*0fdc7489SMatthew Knepley   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberofedges > 0) {
69*0fdc7489SMatthew Knepley     ierr = PetscMalloc1(in->numberofedges*2, &in->edgelist);CHKERRQ(ierr);
70*0fdc7489SMatthew Knepley     ierr = PetscMalloc1(in->numberofedges,   &in->edgemarkerlist);CHKERRQ(ierr);
71*0fdc7489SMatthew Knepley     for (e = eStart; e < eEnd; ++e) {
72*0fdc7489SMatthew Knepley       const PetscInt  idx = e - eStart;
73*0fdc7489SMatthew Knepley       const PetscInt *cone;
74*0fdc7489SMatthew Knepley       PetscInt        coneSize, val;
75*0fdc7489SMatthew Knepley 
76*0fdc7489SMatthew Knepley       ierr = DMPlexGetConeSize(boundary, e, &coneSize);CHKERRQ(ierr);
77*0fdc7489SMatthew Knepley       ierr = DMPlexGetCone(boundary, e, &cone);CHKERRQ(ierr);
78*0fdc7489SMatthew Knepley       in->edgelist[idx*2]     = cone[0] - vStart;
79*0fdc7489SMatthew Knepley       in->edgelist[idx*2 + 1] = cone[1] - vStart;
80*0fdc7489SMatthew Knepley 
81*0fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, e, &val);CHKERRQ(ierr);
82*0fdc7489SMatthew Knepley       in->edgemarkerlist[idx] = (int) val;
83*0fdc7489SMatthew Knepley     }
84*0fdc7489SMatthew Knepley   }
85*0fdc7489SMatthew Knepley 
86*0fdc7489SMatthew Knepley   ierr = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr);
873a074057SBarry Smith   in->numberoffacets = fEnd - fStart;
883a074057SBarry Smith   if (in->numberoffacets > 0) {
893a074057SBarry Smith     ierr = PetscMalloc1(in->numberoffacets, &in->facetlist);CHKERRQ(ierr);
903a074057SBarry Smith     ierr = PetscMalloc1(in->numberoffacets, &in->facetmarkerlist);CHKERRQ(ierr);
913a074057SBarry Smith     for (f = fStart; f < fEnd; ++f) {
923a074057SBarry Smith       const PetscInt idx    = f - fStart;
933a074057SBarry Smith       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, m = -1;
943a074057SBarry Smith       polygon       *poly;
953a074057SBarry Smith 
963a074057SBarry Smith       in->facetlist[idx].numberofpolygons = 1;
973a074057SBarry Smith       ierr = PetscMalloc1(in->facetlist[idx].numberofpolygons, &in->facetlist[idx].polygonlist);CHKERRQ(ierr);
983a074057SBarry Smith       in->facetlist[idx].numberofholes    = 0;
993a074057SBarry Smith       in->facetlist[idx].holelist         = NULL;
1003a074057SBarry Smith 
1013a074057SBarry Smith       ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
1023a074057SBarry Smith       for (p = 0; p < numPoints*2; p += 2) {
1033a074057SBarry Smith         const PetscInt point = points[p];
1043a074057SBarry Smith         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
1053a074057SBarry Smith       }
1063a074057SBarry Smith 
1073a074057SBarry Smith       poly                   = in->facetlist[idx].polygonlist;
1083a074057SBarry Smith       poly->numberofvertices = numVertices;
1093a074057SBarry Smith       ierr = PetscMalloc1(poly->numberofvertices, &poly->vertexlist);CHKERRQ(ierr);
1103a074057SBarry Smith       for (v = 0; v < numVertices; ++v) {
1113a074057SBarry Smith         const PetscInt vIdx = points[v] - vStart;
1123a074057SBarry Smith         poly->vertexlist[v] = vIdx;
1133a074057SBarry Smith       }
114*0fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, f, &m);CHKERRQ(ierr);
1153a074057SBarry Smith       in->facetmarkerlist[idx] = (int) m;
1163a074057SBarry Smith       ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
1173a074057SBarry Smith     }
1183a074057SBarry Smith   }
1193a074057SBarry Smith   if (!rank) {
1203a074057SBarry Smith     TetGenOpts t;
1213a074057SBarry Smith 
1223a074057SBarry Smith     ierr        = TetGenOptsInitialize(&t);CHKERRQ(ierr);
1233a074057SBarry Smith     t.in        = boundary; /* Should go away */
1243a074057SBarry Smith     t.plc       = 1;
1253a074057SBarry Smith     t.quality   = 1;
1263a074057SBarry Smith     t.edgesout  = 1;
1273a074057SBarry Smith     t.zeroindex = 1;
1283a074057SBarry Smith     t.quiet     = 1;
1293a074057SBarry Smith     t.verbose   = verbose;
130*0fdc7489SMatthew Knepley #if 0
131*0fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
132*0fdc7489SMatthew Knepley     /* Need to add in more TetGen code */
133*0fdc7489SMatthew Knepley     t.nobisect  = 1; /* Add Y to preserve Surface Mesh for EGADS */
134*0fdc7489SMatthew Knepley #endif
135*0fdc7489SMatthew Knepley #endif
136*0fdc7489SMatthew Knepley 
1373a074057SBarry Smith     ierr        = TetGenCheckOpts(&t);CHKERRQ(ierr);
1383a074057SBarry Smith     ierr        = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr);
1393a074057SBarry Smith   }
1403a074057SBarry Smith   {
1413a074057SBarry Smith     const PetscInt numCorners  = 4;
1423a074057SBarry Smith     const PetscInt numCells    = out->numberoftetrahedra;
1433a074057SBarry Smith     const PetscInt numVertices = out->numberofpoints;
144*0fdc7489SMatthew Knepley     PetscReal      *meshCoords = NULL;
145*0fdc7489SMatthew Knepley     PetscInt       *cells      = NULL;
1463a074057SBarry Smith 
147a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) == sizeof (out->pointlist[0])) {
148a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *) out->pointlist;
149a4a685f2SJacob Faibussowitsch     } else {
1503a074057SBarry Smith       PetscInt i;
1513a074057SBarry Smith 
152a4a685f2SJacob Faibussowitsch       ierr = PetscMalloc1(dim * numVertices, &meshCoords);CHKERRQ(ierr);
153*0fdc7489SMatthew Knepley       for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out->pointlist[i];
1543a074057SBarry Smith     }
155a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) == sizeof (out->tetrahedronlist[0])) {
156a4a685f2SJacob Faibussowitsch       cells = (PetscInt *) out->tetrahedronlist;
157a4a685f2SJacob Faibussowitsch     } else {
158a4a685f2SJacob Faibussowitsch       PetscInt i;
159a4a685f2SJacob Faibussowitsch 
160a4a685f2SJacob Faibussowitsch       ierr = PetscMalloc1(numCells * numCorners, &cells);CHKERRQ(ierr);
161*0fdc7489SMatthew Knepley       for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt) out->tetrahedronlist[i];
162a4a685f2SJacob Faibussowitsch     }
1633a074057SBarry Smith 
16496ca5757SLisandro Dalcin     ierr = DMPlexInvertCells_CTetgen(numCells, numCorners, cells);CHKERRQ(ierr);
165a4a685f2SJacob Faibussowitsch     ierr = DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr);
166a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) != sizeof (out->pointlist[0])) {
1673a074057SBarry Smith       ierr = PetscFree(meshCoords);CHKERRQ(ierr);
1683a074057SBarry Smith     }
169a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) != sizeof (out->tetrahedronlist[0])) {
170a4a685f2SJacob Faibussowitsch       ierr = PetscFree(cells);CHKERRQ(ierr);
171a4a685f2SJacob Faibussowitsch     }
172*0fdc7489SMatthew Knepley 
1733a074057SBarry Smith     /* Set labels */
174*0fdc7489SMatthew Knepley     ierr = DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dm);CHKERRQ(ierr);
1753a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
1763a074057SBarry Smith       if (out->pointmarkerlist[v]) {
177*0fdc7489SMatthew Knepley         ierr = DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);
1783a074057SBarry Smith       }
1793a074057SBarry Smith     }
1803a074057SBarry Smith     if (interpolate) {
1813a074057SBarry Smith       PetscInt e;
1823a074057SBarry Smith 
1833a074057SBarry Smith       for (e = 0; e < out->numberofedges; e++) {
1843a074057SBarry Smith         if (out->edgemarkerlist[e]) {
1853a074057SBarry Smith           const PetscInt  vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells};
1863a074057SBarry Smith           const PetscInt *edges;
1873a074057SBarry Smith           PetscInt        numEdges;
1883a074057SBarry Smith 
1893a074057SBarry Smith           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1903a074057SBarry Smith           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
191*0fdc7489SMatthew Knepley           ierr = DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);
1923a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1933a074057SBarry Smith         }
1943a074057SBarry Smith       }
1953a074057SBarry Smith       for (f = 0; f < out->numberoftrifaces; f++) {
1963a074057SBarry Smith         if (out->trifacemarkerlist[f]) {
1973a074057SBarry Smith           const PetscInt  vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells};
1983a074057SBarry Smith           const PetscInt *faces;
1993a074057SBarry Smith           PetscInt        numFaces;
2003a074057SBarry Smith 
2013a074057SBarry Smith           ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
2023a074057SBarry Smith           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
203*0fdc7489SMatthew Knepley           ierr = DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);
2043a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
2053a074057SBarry Smith         }
2063a074057SBarry Smith       }
2073a074057SBarry Smith     }
208*0fdc7489SMatthew Knepley 
209*0fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
210*0fdc7489SMatthew Knepley     {
211*0fdc7489SMatthew Knepley       DMLabel        bodyLabel;
212*0fdc7489SMatthew Knepley       PetscContainer modelObj;
213*0fdc7489SMatthew Knepley       PetscInt       cStart, cEnd, c, eStart, eEnd, fStart, fEnd;
214*0fdc7489SMatthew Knepley       ego           *bodies;
215*0fdc7489SMatthew Knepley       ego            model, geom;
216*0fdc7489SMatthew Knepley       int            Nb, oclass, mtype, *senses;
217*0fdc7489SMatthew Knepley 
218*0fdc7489SMatthew Knepley       /* Get Attached EGADS Model from Original DMPlex */
219*0fdc7489SMatthew Knepley       ierr = PetscObjectQuery((PetscObject) boundary, "EGADS Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
220*0fdc7489SMatthew Knepley       if (modelObj) {
221*0fdc7489SMatthew Knepley         ierr = PetscContainerGetPointer(modelObj, (void **) &model);CHKERRQ(ierr);
222*0fdc7489SMatthew Knepley         ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
223*0fdc7489SMatthew Knepley         /* Transfer EGADS Model to Volumetric Mesh */
224*0fdc7489SMatthew Knepley         ierr = PetscObjectCompose((PetscObject) *dm, "EGADS Model", (PetscObject) modelObj);CHKERRQ(ierr);
225*0fdc7489SMatthew Knepley 
226*0fdc7489SMatthew Knepley         /* Set Cell Labels */
227*0fdc7489SMatthew Knepley         ierr = DMGetLabel(*dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
228*0fdc7489SMatthew Knepley         ierr = DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
229*0fdc7489SMatthew Knepley         ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
230*0fdc7489SMatthew Knepley         ierr = DMPlexGetDepthStratum(*dm, 1, &eStart, &eEnd);CHKERRQ(ierr);
231*0fdc7489SMatthew Knepley 
232*0fdc7489SMatthew Knepley         for (c = cStart; c < cEnd; ++c) {
233*0fdc7489SMatthew Knepley           PetscReal centroid[3] = {0., 0., 0.};
234*0fdc7489SMatthew Knepley           PetscInt  b;
235*0fdc7489SMatthew Knepley 
236*0fdc7489SMatthew Knepley           /* Deterimine what body the cell's centroid is located in */
237*0fdc7489SMatthew Knepley           if (!interpolate) {
238*0fdc7489SMatthew Knepley             PetscSection   coordSection;
239*0fdc7489SMatthew Knepley             Vec            coordinates;
240*0fdc7489SMatthew Knepley             PetscScalar   *coords = NULL;
241*0fdc7489SMatthew Knepley             PetscInt       coordSize, s, d;
242*0fdc7489SMatthew Knepley 
243*0fdc7489SMatthew Knepley             ierr = DMGetCoordinatesLocal(*dm, &coordinates);CHKERRQ(ierr);
244*0fdc7489SMatthew Knepley             ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
245*0fdc7489SMatthew Knepley             ierr = DMPlexVecGetClosure(*dm, coordSection, coordinates, c, &coordSize, &coords);CHKERRQ(ierr);
246*0fdc7489SMatthew Knepley             for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d];
247*0fdc7489SMatthew Knepley             ierr = DMPlexVecRestoreClosure(*dm, coordSection, coordinates, c, &coordSize, &coords);CHKERRQ(ierr);
248*0fdc7489SMatthew Knepley           } else {
249*0fdc7489SMatthew Knepley             ierr = DMPlexComputeCellGeometryFVM(*dm, c, NULL, centroid, NULL);CHKERRQ(ierr);
250*0fdc7489SMatthew Knepley           }
251*0fdc7489SMatthew Knepley           for (b = 0; b < Nb; ++b) {
252*0fdc7489SMatthew Knepley             if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;
253*0fdc7489SMatthew Knepley           }
254*0fdc7489SMatthew Knepley           if (b < Nb) {
255*0fdc7489SMatthew Knepley             PetscInt   cval = b, eVal, fVal;
256*0fdc7489SMatthew Knepley             PetscInt *closure = NULL, Ncl, cl;
257*0fdc7489SMatthew Knepley 
258*0fdc7489SMatthew Knepley             ierr = DMLabelSetValue(bodyLabel, c, cval);CHKERRQ(ierr);
259*0fdc7489SMatthew Knepley             ierr = DMPlexGetTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
260*0fdc7489SMatthew Knepley             for (cl = 0; cl < Ncl; cl += 2) {
261*0fdc7489SMatthew Knepley               const PetscInt p = closure[cl];
262*0fdc7489SMatthew Knepley 
263*0fdc7489SMatthew Knepley               if (p >= eStart && p < eEnd) {
264*0fdc7489SMatthew Knepley                 ierr = DMLabelGetValue(bodyLabel, p, &eVal);CHKERRQ(ierr);
265*0fdc7489SMatthew Knepley                 if (eVal < 0) {ierr = DMLabelSetValue(bodyLabel, p, cval);CHKERRQ(ierr);}
266*0fdc7489SMatthew Knepley               }
267*0fdc7489SMatthew Knepley               if (p >= fStart && p < fEnd) {
268*0fdc7489SMatthew Knepley                 ierr = DMLabelGetValue(bodyLabel, p, &fVal);CHKERRQ(ierr);
269*0fdc7489SMatthew Knepley                 if (fVal < 0) {ierr = DMLabelSetValue(bodyLabel, p, cval);CHKERRQ(ierr);}
270*0fdc7489SMatthew Knepley               }
271*0fdc7489SMatthew Knepley             }
272*0fdc7489SMatthew Knepley             ierr = DMPlexRestoreTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
273*0fdc7489SMatthew Knepley           }
274*0fdc7489SMatthew Knepley         }
275*0fdc7489SMatthew Knepley       }
276*0fdc7489SMatthew Knepley     }
277*0fdc7489SMatthew Knepley #endif
2783a074057SBarry Smith     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
2793a074057SBarry Smith   }
2803a074057SBarry Smith 
281*0fdc7489SMatthew Knepley   ierr = DMUniversalLabelDestroy(&universal);CHKERRQ(ierr);
2823a074057SBarry Smith   ierr = PLCDestroy(&in);CHKERRQ(ierr);
2833a074057SBarry Smith   ierr = PLCDestroy(&out);CHKERRQ(ierr);
2843a074057SBarry Smith   PetscFunctionReturn(0);
2853a074057SBarry Smith }
2863a074057SBarry Smith 
2873a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined)
2883a074057SBarry Smith {
2893a074057SBarry Smith   MPI_Comm               comm;
2903a074057SBarry Smith   const PetscInt         dim = 3;
2913a074057SBarry Smith   PLC                   *in, *out;
292*0fdc7489SMatthew Knepley   DMUniversalLabel       universal;
293*0fdc7489SMatthew Knepley   PetscInt               vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, cStart, cEnd, c, verbose = 0;
294*0fdc7489SMatthew Knepley   DMPlexInterpolatedFlag isInterpolated;
2953a074057SBarry Smith   PetscMPIInt            rank;
2963a074057SBarry Smith   PetscErrorCode         ierr;
2973a074057SBarry Smith 
2983a074057SBarry Smith   PetscFunctionBegin;
2993a074057SBarry Smith   ierr = PetscOptionsGetInt(NULL,((PetscObject) dm)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr);
300*0fdc7489SMatthew Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
3013a074057SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
302*0fdc7489SMatthew Knepley   ierr = DMPlexIsInterpolatedCollective(dm, &isInterpolated);CHKERRQ(ierr);
303*0fdc7489SMatthew Knepley   ierr = DMUniversalLabelCreate(dm, &universal);CHKERRQ(ierr);
304*0fdc7489SMatthew Knepley 
3053a074057SBarry Smith   ierr = PLCCreate(&in);CHKERRQ(ierr);
3063a074057SBarry Smith   ierr = PLCCreate(&out);CHKERRQ(ierr);
3073a074057SBarry Smith 
308*0fdc7489SMatthew Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
3093a074057SBarry Smith   in->numberofpoints = vEnd - vStart;
3103a074057SBarry Smith   if (in->numberofpoints > 0) {
3113a074057SBarry Smith     PetscSection coordSection;
3123a074057SBarry Smith     Vec          coordinates;
3133a074057SBarry Smith     PetscScalar *array;
3143a074057SBarry Smith 
3153a074057SBarry Smith     ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr);
3163a074057SBarry Smith     ierr = PetscMalloc1(in->numberofpoints, &in->pointmarkerlist);CHKERRQ(ierr);
3173a074057SBarry Smith     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
3183a074057SBarry Smith     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
3193a074057SBarry Smith     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
3203a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
3213a074057SBarry Smith       const PetscInt idx = v - vStart;
322*0fdc7489SMatthew Knepley       PetscInt       off, d, m;
3233a074057SBarry Smith 
3243a074057SBarry Smith       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
325*0fdc7489SMatthew Knepley       for (d = 0; d < dim; ++d) in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
326*0fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, v, &m);CHKERRQ(ierr);
3273a074057SBarry Smith       in->pointmarkerlist[idx] = (int) m;
3283a074057SBarry Smith     }
3293a074057SBarry Smith     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
3303a074057SBarry Smith   }
3313a074057SBarry Smith 
332*0fdc7489SMatthew Knepley   ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);
333*0fdc7489SMatthew Knepley   in->numberofedges = eEnd - eStart;
334*0fdc7489SMatthew Knepley   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberofedges > 0) {
335*0fdc7489SMatthew Knepley     ierr = PetscMalloc1(in->numberofedges * 2, &in->edgelist);CHKERRQ(ierr);
336*0fdc7489SMatthew Knepley     ierr = PetscMalloc1(in->numberofedges,     &in->edgemarkerlist);CHKERRQ(ierr);
337*0fdc7489SMatthew Knepley     for (e = eStart; e < eEnd; ++e) {
338*0fdc7489SMatthew Knepley       const PetscInt  idx = e - eStart;
339*0fdc7489SMatthew Knepley       const PetscInt *cone;
340*0fdc7489SMatthew Knepley       PetscInt        coneSize, val;
341*0fdc7489SMatthew Knepley 
342*0fdc7489SMatthew Knepley       ierr = DMPlexGetConeSize(dm, e, &coneSize);CHKERRQ(ierr);
343*0fdc7489SMatthew Knepley       ierr = DMPlexGetCone(dm, e, &cone);CHKERRQ(ierr);
344*0fdc7489SMatthew Knepley       in->edgelist[idx*2]     = cone[0] - vStart;
345*0fdc7489SMatthew Knepley       in->edgelist[idx*2 + 1] = cone[1] - vStart;
346*0fdc7489SMatthew Knepley 
347*0fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, e, &val);CHKERRQ(ierr);
348*0fdc7489SMatthew Knepley       in->edgemarkerlist[idx] = (int) val;
349*0fdc7489SMatthew Knepley     }
350*0fdc7489SMatthew Knepley   }
351*0fdc7489SMatthew Knepley 
352*0fdc7489SMatthew Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
353*0fdc7489SMatthew Knepley   in->numberoftrifaces = 0;
354*0fdc7489SMatthew Knepley   for (f = fStart; f < fEnd; ++f) {
355*0fdc7489SMatthew Knepley     PetscInt supportSize;
356*0fdc7489SMatthew Knepley 
357*0fdc7489SMatthew Knepley     ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
358*0fdc7489SMatthew Knepley     if (supportSize == 1) ++in->numberoftrifaces;
359*0fdc7489SMatthew Knepley   }
360*0fdc7489SMatthew Knepley   if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberoftrifaces > 0) {
361*0fdc7489SMatthew Knepley     PetscInt tf = 0;
362*0fdc7489SMatthew Knepley 
363*0fdc7489SMatthew Knepley     ierr = PetscMalloc1(in->numberoftrifaces*3, &in->trifacelist);CHKERRQ(ierr);
364*0fdc7489SMatthew Knepley     ierr = PetscMalloc1(in->numberoftrifaces, &in->trifacemarkerlist);CHKERRQ(ierr);
365*0fdc7489SMatthew Knepley     for (f = fStart; f < fEnd; ++f) {
366*0fdc7489SMatthew Knepley       PetscInt *points = NULL;
367*0fdc7489SMatthew Knepley       PetscInt supportSize, numPoints, p, Nv = 0, val;
368*0fdc7489SMatthew Knepley 
369*0fdc7489SMatthew Knepley       ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
370*0fdc7489SMatthew Knepley       if (supportSize != 1) continue;
371*0fdc7489SMatthew Knepley       ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
372*0fdc7489SMatthew Knepley       for (p = 0; p < numPoints*2; p += 2) {
373*0fdc7489SMatthew Knepley         const PetscInt point = points[p];
374*0fdc7489SMatthew Knepley         if ((point >= vStart) && (point < vEnd)) in->trifacelist[tf*3 + Nv++] = point - vStart;
375*0fdc7489SMatthew Knepley       }
376*0fdc7489SMatthew Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
377*0fdc7489SMatthew Knepley       if (Nv != 3) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %D has %D vertices, not 3", f, Nv);
378*0fdc7489SMatthew Knepley       ierr = DMLabelGetValue(universal->label, f, &val);CHKERRQ(ierr);
379*0fdc7489SMatthew Knepley       in->trifacemarkerlist[tf] = (int) val;
380*0fdc7489SMatthew Knepley       ++tf;
381*0fdc7489SMatthew Knepley     }
382*0fdc7489SMatthew Knepley   }
383*0fdc7489SMatthew Knepley 
384*0fdc7489SMatthew Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
3853a074057SBarry Smith   in->numberofcorners       = 4;
3863a074057SBarry Smith   in->numberoftetrahedra    = cEnd - cStart;
3873a074057SBarry Smith   in->tetrahedronvolumelist = maxVolumes;
3883a074057SBarry Smith   if (in->numberoftetrahedra > 0) {
3893a074057SBarry Smith     ierr = PetscMalloc1(in->numberoftetrahedra*in->numberofcorners, &in->tetrahedronlist);CHKERRQ(ierr);
3903a074057SBarry Smith     for (c = cStart; c < cEnd; ++c) {
3913a074057SBarry Smith       const PetscInt idx     = c - cStart;
3923a074057SBarry Smith       PetscInt      *closure = NULL;
3933a074057SBarry Smith       PetscInt       closureSize;
3943a074057SBarry Smith 
3953a074057SBarry Smith       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
3963a074057SBarry 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);
397*0fdc7489SMatthew Knepley       for (v = 0; v < 4; ++v) in->tetrahedronlist[idx*in->numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
3983a074057SBarry Smith       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
3993a074057SBarry Smith     }
4003a074057SBarry Smith   }
401*0fdc7489SMatthew Knepley 
4023a074057SBarry Smith   if (!rank) {
4033a074057SBarry Smith     TetGenOpts t;
4043a074057SBarry Smith 
4053a074057SBarry Smith     ierr = TetGenOptsInitialize(&t);CHKERRQ(ierr);
4063a074057SBarry Smith 
4073a074057SBarry Smith     t.in        = dm; /* Should go away */
4083a074057SBarry Smith     t.refine    = 1;
4093a074057SBarry Smith     t.varvolume = 1;
4103a074057SBarry Smith     t.quality   = 1;
4113a074057SBarry Smith     t.edgesout  = 1;
4123a074057SBarry Smith     t.zeroindex = 1;
4133a074057SBarry Smith     t.quiet     = 1;
4143a074057SBarry Smith     t.verbose   = verbose; /* Change this */
4153a074057SBarry Smith 
4163a074057SBarry Smith     ierr = TetGenCheckOpts(&t);CHKERRQ(ierr);
4173a074057SBarry Smith     ierr = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr);
4183a074057SBarry Smith   }
419*0fdc7489SMatthew Knepley 
4203a074057SBarry Smith   in->tetrahedronvolumelist = NULL;
4213a074057SBarry Smith   {
4223a074057SBarry Smith     const PetscInt numCorners  = 4;
4233a074057SBarry Smith     const PetscInt numCells    = out->numberoftetrahedra;
4243a074057SBarry Smith     const PetscInt numVertices = out->numberofpoints;
425*0fdc7489SMatthew Knepley     PetscReal      *meshCoords = NULL;
426*0fdc7489SMatthew Knepley     PetscInt       *cells      = NULL;
427*0fdc7489SMatthew Knepley     PetscBool      interpolate = isInterpolated == DMPLEX_INTERPOLATED_FULL ? PETSC_TRUE : PETSC_FALSE;
4283a074057SBarry Smith 
429a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) == sizeof (out->pointlist[0])) {
430a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *) out->pointlist;
431a4a685f2SJacob Faibussowitsch     } else {
4323a074057SBarry Smith       PetscInt i;
4333a074057SBarry Smith 
434a4a685f2SJacob Faibussowitsch       ierr = PetscMalloc1(dim * numVertices, &meshCoords);CHKERRQ(ierr);
435*0fdc7489SMatthew Knepley       for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out->pointlist[i];
4363a074057SBarry Smith     }
437a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) == sizeof (out->tetrahedronlist[0])) {
438a4a685f2SJacob Faibussowitsch       cells = (PetscInt *) out->tetrahedronlist;
439a4a685f2SJacob Faibussowitsch     } else {
440a4a685f2SJacob Faibussowitsch       PetscInt i;
441a4a685f2SJacob Faibussowitsch 
442a4a685f2SJacob Faibussowitsch       ierr = PetscMalloc1(numCells * numCorners, &cells);CHKERRQ(ierr);
443*0fdc7489SMatthew Knepley       for (i = 0; i < numCells * numCorners; ++i) cells[i] = (PetscInt) out->tetrahedronlist[i];
444a4a685f2SJacob Faibussowitsch     }
4453a074057SBarry Smith 
44696ca5757SLisandro Dalcin     ierr = DMPlexInvertCells_CTetgen(numCells, numCorners, cells);CHKERRQ(ierr);
447a4a685f2SJacob Faibussowitsch     ierr = DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
448*0fdc7489SMatthew Knepley     if (sizeof (PetscReal) != sizeof (out->pointlist[0])) {ierr = PetscFree(meshCoords);CHKERRQ(ierr);}
449*0fdc7489SMatthew Knepley     if (sizeof (PetscInt) != sizeof (out->tetrahedronlist[0])) {ierr = PetscFree(cells);CHKERRQ(ierr);}
450*0fdc7489SMatthew Knepley 
4513a074057SBarry Smith     /* Set labels */
452*0fdc7489SMatthew Knepley     ierr = DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dmRefined);CHKERRQ(ierr);
4533a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
4543a074057SBarry Smith       if (out->pointmarkerlist[v]) {
455*0fdc7489SMatthew Knepley         ierr = DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);
4563a074057SBarry Smith       }
4573a074057SBarry Smith     }
4583a074057SBarry Smith     if (interpolate) {
4593a074057SBarry Smith       PetscInt e, f;
4603a074057SBarry Smith 
4613a074057SBarry Smith       for (e = 0; e < out->numberofedges; e++) {
4623a074057SBarry Smith         if (out->edgemarkerlist[e]) {
4633a074057SBarry Smith           const PetscInt  vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells};
4643a074057SBarry Smith           const PetscInt *edges;
4653a074057SBarry Smith           PetscInt        numEdges;
4663a074057SBarry Smith 
4673a074057SBarry Smith           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
4683a074057SBarry Smith           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
469*0fdc7489SMatthew Knepley           ierr = DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);
4703a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
4713a074057SBarry Smith         }
4723a074057SBarry Smith       }
4733a074057SBarry Smith       for (f = 0; f < out->numberoftrifaces; f++) {
4743a074057SBarry Smith         if (out->trifacemarkerlist[f]) {
4753a074057SBarry Smith           const PetscInt  vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells};
4763a074057SBarry Smith           const PetscInt *faces;
4773a074057SBarry Smith           PetscInt        numFaces;
4783a074057SBarry Smith 
4793a074057SBarry Smith           ierr = DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
4803a074057SBarry Smith           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
481*0fdc7489SMatthew Knepley           ierr = DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);
4823a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
4833a074057SBarry Smith         }
4843a074057SBarry Smith       }
4853a074057SBarry Smith     }
486*0fdc7489SMatthew Knepley 
487*0fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS
488*0fdc7489SMatthew Knepley     {
489*0fdc7489SMatthew Knepley       DMLabel        bodyLabel;
490*0fdc7489SMatthew Knepley       PetscContainer modelObj;
491*0fdc7489SMatthew Knepley       PetscInt       cStart, cEnd, c, eStart, eEnd, fStart, fEnd;
492*0fdc7489SMatthew Knepley       ego           *bodies;
493*0fdc7489SMatthew Knepley       ego            model, geom;
494*0fdc7489SMatthew Knepley       int            Nb, oclass, mtype, *senses;
495*0fdc7489SMatthew Knepley 
496*0fdc7489SMatthew Knepley       /* Get Attached EGADS Model from Original DMPlex */
497*0fdc7489SMatthew Knepley       ierr = PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
498*0fdc7489SMatthew Knepley       if (modelObj) {
499*0fdc7489SMatthew Knepley         ierr = PetscContainerGetPointer(modelObj, (void **) &model);CHKERRQ(ierr);
500*0fdc7489SMatthew Knepley         ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
501*0fdc7489SMatthew Knepley         /* Transfer EGADS Model to Volumetric Mesh */
502*0fdc7489SMatthew Knepley         ierr = PetscObjectCompose((PetscObject) *dmRefined, "EGADS Model", (PetscObject) modelObj);CHKERRQ(ierr);
503*0fdc7489SMatthew Knepley 
504*0fdc7489SMatthew Knepley         /* Set Cell Labels */
505*0fdc7489SMatthew Knepley         ierr = DMGetLabel(*dmRefined, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
506*0fdc7489SMatthew Knepley         ierr = DMPlexGetHeightStratum(*dmRefined, 0, &cStart, &cEnd);CHKERRQ(ierr);
507*0fdc7489SMatthew Knepley         ierr = DMPlexGetHeightStratum(*dmRefined, 1, &fStart, &fEnd);CHKERRQ(ierr);
508*0fdc7489SMatthew Knepley         ierr = DMPlexGetDepthStratum(*dmRefined, 1, &eStart, &eEnd);CHKERRQ(ierr);
509*0fdc7489SMatthew Knepley 
510*0fdc7489SMatthew Knepley         for (c = cStart; c < cEnd; ++c) {
511*0fdc7489SMatthew Knepley           PetscReal centroid[3] = {0., 0., 0.};
512*0fdc7489SMatthew Knepley           PetscInt  b;
513*0fdc7489SMatthew Knepley 
514*0fdc7489SMatthew Knepley           /* Deterimine what body the cell's centroid is located in */
515*0fdc7489SMatthew Knepley           if (!interpolate) {
516*0fdc7489SMatthew Knepley             PetscSection   coordSection;
517*0fdc7489SMatthew Knepley             Vec            coordinates;
518*0fdc7489SMatthew Knepley             PetscScalar   *coords = NULL;
519*0fdc7489SMatthew Knepley             PetscInt       coordSize, s, d;
520*0fdc7489SMatthew Knepley 
521*0fdc7489SMatthew Knepley             ierr = DMGetCoordinatesLocal(*dmRefined, &coordinates);CHKERRQ(ierr);
522*0fdc7489SMatthew Knepley             ierr = DMGetCoordinateSection(*dmRefined, &coordSection);CHKERRQ(ierr);
523*0fdc7489SMatthew Knepley             ierr = DMPlexVecGetClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords);CHKERRQ(ierr);
524*0fdc7489SMatthew Knepley             for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d];
525*0fdc7489SMatthew Knepley             ierr = DMPlexVecRestoreClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords);CHKERRQ(ierr);
526*0fdc7489SMatthew Knepley           } else {
527*0fdc7489SMatthew Knepley             ierr = DMPlexComputeCellGeometryFVM(*dmRefined, c, NULL, centroid, NULL);CHKERRQ(ierr);
528*0fdc7489SMatthew Knepley           }
529*0fdc7489SMatthew Knepley           for (b = 0; b < Nb; ++b) {
530*0fdc7489SMatthew Knepley             if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break;
531*0fdc7489SMatthew Knepley           }
532*0fdc7489SMatthew Knepley           if (b < Nb) {
533*0fdc7489SMatthew Knepley             PetscInt   cval = b, eVal, fVal;
534*0fdc7489SMatthew Knepley             PetscInt *closure = NULL, Ncl, cl;
535*0fdc7489SMatthew Knepley 
536*0fdc7489SMatthew Knepley             ierr = DMLabelSetValue(bodyLabel, c, cval);CHKERRQ(ierr);
537*0fdc7489SMatthew Knepley             ierr = DMPlexGetTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
538*0fdc7489SMatthew Knepley             for (cl = 0; cl < Ncl; cl += 2) {
539*0fdc7489SMatthew Knepley               const PetscInt p = closure[cl];
540*0fdc7489SMatthew Knepley 
541*0fdc7489SMatthew Knepley               if (p >= eStart && p < eEnd) {
542*0fdc7489SMatthew Knepley                 ierr = DMLabelGetValue(bodyLabel, p, &eVal);CHKERRQ(ierr);
543*0fdc7489SMatthew Knepley                 if (eVal < 0) {ierr = DMLabelSetValue(bodyLabel, p, cval);CHKERRQ(ierr);}
544*0fdc7489SMatthew Knepley               }
545*0fdc7489SMatthew Knepley               if (p >= fStart && p < fEnd) {
546*0fdc7489SMatthew Knepley                 ierr = DMLabelGetValue(bodyLabel, p, &fVal);CHKERRQ(ierr);
547*0fdc7489SMatthew Knepley                 if (fVal < 0) {ierr = DMLabelSetValue(bodyLabel, p, cval);CHKERRQ(ierr);}
548*0fdc7489SMatthew Knepley               }
549*0fdc7489SMatthew Knepley             }
550*0fdc7489SMatthew Knepley             ierr = DMPlexRestoreTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
551*0fdc7489SMatthew Knepley           }
552*0fdc7489SMatthew Knepley         }
553*0fdc7489SMatthew Knepley       }
554*0fdc7489SMatthew Knepley     }
555*0fdc7489SMatthew Knepley #endif
5563a074057SBarry Smith     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
5573a074057SBarry Smith   }
558*0fdc7489SMatthew Knepley   ierr = DMUniversalLabelDestroy(&universal);CHKERRQ(ierr);
5593a074057SBarry Smith   ierr = PLCDestroy(&in);CHKERRQ(ierr);
5603a074057SBarry Smith   ierr = PLCDestroy(&out);CHKERRQ(ierr);
5613a074057SBarry Smith   PetscFunctionReturn(0);
5623a074057SBarry Smith }
563