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 13*f22e26b7SPierre Jolivet #if defined(__clang__) 14*f22e26b7SPierre Jolivet #pragma clang diagnostic push 15*f22e26b7SPierre Jolivet #pragma clang diagnostic ignored "-Wunused-parameter" 16*f22e26b7SPierre Jolivet #pragma clang diagnostic ignored "-Wzero-as-null-pointer-constant" 17*f22e26b7SPierre Jolivet #elif defined(__GNUC__) || defined(__GNUG__) 18*f22e26b7SPierre Jolivet #pragma GCC diagnostic push 19*f22e26b7SPierre Jolivet #pragma GCC diagnostic ignored "-Wunused-parameter" 20*f22e26b7SPierre Jolivet #endif 213a074057SBarry Smith #include <tetgen.h> 22*f22e26b7SPierre Jolivet #if defined(__clang__) 23*f22e26b7SPierre Jolivet #pragma clang diagnostic pop 24*f22e26b7SPierre Jolivet #elif defined(__GNUC__) || defined(__GNUG__) 25*f22e26b7SPierre Jolivet #pragma GCC diagnostic pop 26*f22e26b7SPierre Jolivet #endif 273a074057SBarry Smith 283a074057SBarry Smith /* This is to fix the tetrahedron orientation from TetGen */ 29d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexInvertCells_Tetgen(PetscInt numCells, PetscInt numCorners, PetscInt cells[]) 30d71ae5a4SJacob Faibussowitsch { 313a074057SBarry Smith PetscInt bound = numCells * numCorners, coff; 323a074057SBarry Smith 333a074057SBarry Smith PetscFunctionBegin; 349371c9d4SSatish Balay #define SWAP(a, b) \ 359371c9d4SSatish Balay do { \ 369371c9d4SSatish Balay PetscInt tmp = (a); \ 379371c9d4SSatish Balay (a) = (b); \ 389371c9d4SSatish Balay (b) = tmp; \ 399371c9d4SSatish Balay } while (0) 4096ca5757SLisandro Dalcin for (coff = 0; coff < bound; coff += numCorners) SWAP(cells[coff], cells[coff + 1]); 4196ca5757SLisandro Dalcin #undef SWAP 423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 433a074057SBarry Smith } 443a074057SBarry Smith 45d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexGenerate_Tetgen(DM boundary, PetscBool interpolate, DM *dm) 46d71ae5a4SJacob Faibussowitsch { 473a074057SBarry Smith MPI_Comm comm; 483a074057SBarry Smith const PetscInt dim = 3; 493a074057SBarry Smith ::tetgenio in; 503a074057SBarry Smith ::tetgenio out; 519318fe57SMatthew G. Knepley PetscContainer modelObj; 520fdc7489SMatthew Knepley DMUniversalLabel universal; 53af226901SMatthew G. Knepley PetscInt vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, defVal; 540fdc7489SMatthew Knepley DMPlexInterpolatedFlag isInterpolated; 553a074057SBarry Smith PetscMPIInt rank; 563a074057SBarry Smith 573a074057SBarry Smith PetscFunctionBegin; 589566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)boundary, &comm)); 599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 609566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(boundary, &isInterpolated)); 619566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelCreate(boundary, &universal)); 62af226901SMatthew G. Knepley PetscCall(DMLabelGetDefaultValue(universal->label, &defVal)); 633a074057SBarry Smith 649566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd)); 653a074057SBarry Smith in.numberofpoints = vEnd - vStart; 663a074057SBarry Smith if (in.numberofpoints > 0) { 673a074057SBarry Smith PetscSection coordSection; 683a074057SBarry Smith Vec coordinates; 690fdc7489SMatthew Knepley const PetscScalar *array; 703a074057SBarry Smith 713a074057SBarry Smith in.pointlist = new double[in.numberofpoints * dim]; 723a074057SBarry Smith in.pointmarkerlist = new int[in.numberofpoints]; 733a074057SBarry Smith 7441e9d8b5SMatthew G. Knepley PetscCall(PetscArrayzero(in.pointmarkerlist, (size_t)in.numberofpoints)); 759566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(boundary, &coordinates)); 769566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(boundary, &coordSection)); 779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &array)); 783a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 793a074057SBarry Smith const PetscInt idx = v - vStart; 800fdc7489SMatthew Knepley PetscInt off, d, val; 813a074057SBarry Smith 829566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 833a074057SBarry Smith for (d = 0; d < dim; ++d) in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]); 849566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, v, &val)); 85af226901SMatthew G. Knepley if (val != defVal) in.pointmarkerlist[idx] = (int)val; 863a074057SBarry Smith } 879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &array)); 883a074057SBarry Smith } 893a074057SBarry Smith 909566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(boundary, 1, &eStart, &eEnd)); 910fdc7489SMatthew Knepley in.numberofedges = eEnd - eStart; 920fdc7489SMatthew Knepley if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberofedges > 0) { 930fdc7489SMatthew Knepley in.edgelist = new int[in.numberofedges * 2]; 940fdc7489SMatthew Knepley in.edgemarkerlist = new int[in.numberofedges]; 950fdc7489SMatthew Knepley for (e = eStart; e < eEnd; ++e) { 960fdc7489SMatthew Knepley const PetscInt idx = e - eStart; 970fdc7489SMatthew Knepley const PetscInt *cone; 980fdc7489SMatthew Knepley PetscInt coneSize, val; 990fdc7489SMatthew Knepley 1009566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(boundary, e, &coneSize)); 1019566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(boundary, e, &cone)); 1020fdc7489SMatthew Knepley in.edgelist[idx * 2] = cone[0] - vStart; 1030fdc7489SMatthew Knepley in.edgelist[idx * 2 + 1] = cone[1] - vStart; 1040fdc7489SMatthew Knepley 1059566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, e, &val)); 106af226901SMatthew G. Knepley if (val != defVal) in.edgemarkerlist[idx] = (int)val; 1070fdc7489SMatthew Knepley } 1080fdc7489SMatthew Knepley } 1090fdc7489SMatthew Knepley 1109566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd)); 1113a074057SBarry Smith in.numberoffacets = fEnd - fStart; 1123a074057SBarry Smith if (in.numberoffacets > 0) { 1133a074057SBarry Smith in.facetlist = new tetgenio::facet[in.numberoffacets]; 1143a074057SBarry Smith in.facetmarkerlist = new int[in.numberoffacets]; 1153a074057SBarry Smith for (f = fStart; f < fEnd; ++f) { 1163a074057SBarry Smith const PetscInt idx = f - fStart; 117*f22e26b7SPierre Jolivet PetscInt *points = nullptr, numPoints, p, numVertices = 0, v, val = -1; 1183a074057SBarry Smith 1193a074057SBarry Smith in.facetlist[idx].numberofpolygons = 1; 1203a074057SBarry Smith in.facetlist[idx].polygonlist = new tetgenio::polygon[in.facetlist[idx].numberofpolygons]; 1213a074057SBarry Smith in.facetlist[idx].numberofholes = 0; 122*f22e26b7SPierre Jolivet in.facetlist[idx].holelist = nullptr; 1233a074057SBarry Smith 1249566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points)); 1253a074057SBarry Smith for (p = 0; p < numPoints * 2; p += 2) { 1263a074057SBarry Smith const PetscInt point = points[p]; 1273a074057SBarry Smith if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point; 1283a074057SBarry Smith } 1293a074057SBarry Smith 1303a074057SBarry Smith tetgenio::polygon *poly = in.facetlist[idx].polygonlist; 1313a074057SBarry Smith poly->numberofvertices = numVertices; 1323a074057SBarry Smith poly->vertexlist = new int[poly->numberofvertices]; 1333a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 1343a074057SBarry Smith const PetscInt vIdx = points[v] - vStart; 1353a074057SBarry Smith poly->vertexlist[v] = vIdx; 1363a074057SBarry Smith } 1379566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, f, &val)); 138af226901SMatthew G. Knepley if (val != defVal) in.facetmarkerlist[idx] = (int)val; 1399566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points)); 1403a074057SBarry Smith } 1413a074057SBarry Smith } 142dd400576SPatrick Sanan if (rank == 0) { 1430fdc7489SMatthew Knepley DM_Plex *mesh = (DM_Plex *)boundary->data; 1443a074057SBarry Smith char args[32]; 1453a074057SBarry Smith 1463a074057SBarry Smith /* Take away 'Q' for verbose output */ 1470fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS 148c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(args, "pqezQY", sizeof(args))); 1490fdc7489SMatthew Knepley #else 150c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(args, "pqezQ", sizeof(args))); 1510fdc7489SMatthew Knepley #endif 1529371c9d4SSatish Balay if (mesh->tetgenOpts) { 1539371c9d4SSatish Balay ::tetrahedralize(mesh->tetgenOpts, &in, &out); 1549371c9d4SSatish Balay } else { 1559371c9d4SSatish Balay ::tetrahedralize(args, &in, &out); 1569371c9d4SSatish Balay } 1573a074057SBarry Smith } 1583a074057SBarry Smith { 1593a074057SBarry Smith const PetscInt numCorners = 4; 1603a074057SBarry Smith const PetscInt numCells = out.numberoftetrahedra; 1613a074057SBarry Smith const PetscInt numVertices = out.numberofpoints; 162*f22e26b7SPierre Jolivet PetscReal *meshCoords = nullptr; 163*f22e26b7SPierre Jolivet PetscInt *cells = nullptr; 164a4a685f2SJacob Faibussowitsch 165a4a685f2SJacob Faibussowitsch if (sizeof(PetscReal) == sizeof(out.pointlist[0])) { 166a4a685f2SJacob Faibussowitsch meshCoords = (PetscReal *)out.pointlist; 167a4a685f2SJacob Faibussowitsch } else { 168a4a685f2SJacob Faibussowitsch PetscInt i; 169a4a685f2SJacob Faibussowitsch 170a4a685f2SJacob Faibussowitsch meshCoords = new PetscReal[dim * numVertices]; 1710fdc7489SMatthew Knepley for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal)out.pointlist[i]; 172a4a685f2SJacob Faibussowitsch } 173a4a685f2SJacob Faibussowitsch if (sizeof(PetscInt) == sizeof(out.tetrahedronlist[0])) { 174a4a685f2SJacob Faibussowitsch cells = (PetscInt *)out.tetrahedronlist; 175a4a685f2SJacob Faibussowitsch } else { 176a4a685f2SJacob Faibussowitsch PetscInt i; 177a4a685f2SJacob Faibussowitsch 178a4a685f2SJacob Faibussowitsch cells = new PetscInt[numCells * numCorners]; 1790fdc7489SMatthew Knepley for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out.tetrahedronlist[i]; 180a4a685f2SJacob Faibussowitsch } 1813a074057SBarry Smith 1829566063dSJacob Faibussowitsch PetscCall(DMPlexInvertCells_Tetgen(numCells, numCorners, cells)); 1839566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm)); 1840fdc7489SMatthew Knepley 1853a074057SBarry Smith /* Set labels */ 1869566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dm)); 1873a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 18848a46eb9SPierre Jolivet if (out.pointmarkerlist[v]) PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, v + numCells, out.pointmarkerlist[v])); 1893a074057SBarry Smith } 1903a074057SBarry Smith if (interpolate) { 1913a074057SBarry Smith PetscInt e; 1923a074057SBarry Smith 1933a074057SBarry Smith for (e = 0; e < out.numberofedges; e++) { 1943a074057SBarry Smith if (out.edgemarkerlist[e]) { 1953a074057SBarry Smith const PetscInt vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells}; 1963a074057SBarry Smith const PetscInt *edges; 1973a074057SBarry Smith PetscInt numEdges; 1983a074057SBarry Smith 1999566063dSJacob Faibussowitsch PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges)); 20063a3b9bcSJacob Faibussowitsch PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges); 2019566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, edges[0], out.edgemarkerlist[e])); 2029566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges)); 2033a074057SBarry Smith } 2043a074057SBarry Smith } 2053a074057SBarry Smith for (f = 0; f < out.numberoftrifaces; f++) { 2063a074057SBarry Smith if (out.trifacemarkerlist[f]) { 2073a074057SBarry Smith const PetscInt vertices[3] = {out.trifacelist[f * 3 + 0] + numCells, out.trifacelist[f * 3 + 1] + numCells, out.trifacelist[f * 3 + 2] + numCells}; 2083a074057SBarry Smith const PetscInt *faces; 2093a074057SBarry Smith PetscInt numFaces; 2103a074057SBarry Smith 2119566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces)); 21263a3b9bcSJacob Faibussowitsch PetscCheck(numFaces == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces); 2139566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, faces[0], out.trifacemarkerlist[f])); 2149566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces)); 2153a074057SBarry Smith } 2163a074057SBarry Smith } 2173a074057SBarry Smith } 2180fdc7489SMatthew Knepley 2199566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)boundary, "EGADS Model", (PetscObject *)&modelObj)); 2209318fe57SMatthew G. Knepley if (modelObj) { 2210fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS 2220fdc7489SMatthew Knepley DMLabel bodyLabel; 2230fdc7489SMatthew Knepley PetscInt cStart, cEnd, c, eStart, eEnd, fStart, fEnd; 224c1cad2e7SMatthew G. Knepley PetscBool islite = PETSC_FALSE; 2250fdc7489SMatthew Knepley ego *bodies; 2260fdc7489SMatthew Knepley ego model, geom; 2270fdc7489SMatthew Knepley int Nb, oclass, mtype, *senses; 2280fdc7489SMatthew Knepley 2290fdc7489SMatthew Knepley /* Get Attached EGADS Model from Original DMPlex */ 2309566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)boundary, "EGADS Model", (PetscObject *)&modelObj)); 231c1cad2e7SMatthew G. Knepley if (modelObj) { 2329566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(modelObj, (void **)&model)); 233*f22e26b7SPierre Jolivet PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, nullptr, &Nb, &bodies, &senses)); 2340fdc7489SMatthew Knepley /* Transfer EGADS Model to Volumetric Mesh */ 2359566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)*dm, "EGADS Model", (PetscObject)modelObj)); 236c1cad2e7SMatthew G. Knepley } else { 2379566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)boundary, "EGADSLite Model", (PetscObject *)&modelObj)); 238c1cad2e7SMatthew G. Knepley if (modelObj) { 2399566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(modelObj, (void **)&model)); 240*f22e26b7SPierre Jolivet PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, nullptr, &Nb, &bodies, &senses)); 241c1cad2e7SMatthew G. Knepley /* Transfer EGADS Model to Volumetric Mesh */ 2429566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)*dm, "EGADSLite Model", (PetscObject)modelObj)); 243c1cad2e7SMatthew G. Knepley islite = PETSC_TRUE; 244c1cad2e7SMatthew G. Knepley } 245c1cad2e7SMatthew G. Knepley } 246c1cad2e7SMatthew G. Knepley if (!modelObj) goto skip_egads; 2470fdc7489SMatthew Knepley 2480fdc7489SMatthew Knepley /* Set Cell Labels */ 2499566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "EGADS Body ID", &bodyLabel)); 2509566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 2519566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd)); 2529566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(*dm, 1, &eStart, &eEnd)); 2530fdc7489SMatthew Knepley 2540fdc7489SMatthew Knepley for (c = cStart; c < cEnd; ++c) { 2550fdc7489SMatthew Knepley PetscReal centroid[3] = {0., 0., 0.}; 2560fdc7489SMatthew Knepley PetscInt b; 2570fdc7489SMatthew Knepley 25835cb6cd3SPierre Jolivet /* Determine what body the cell's centroid is located in */ 2590fdc7489SMatthew Knepley if (!interpolate) { 2600fdc7489SMatthew Knepley PetscSection coordSection; 2610fdc7489SMatthew Knepley Vec coordinates; 262*f22e26b7SPierre Jolivet PetscScalar *coords = nullptr; 2630fdc7489SMatthew Knepley PetscInt coordSize, s, d; 2640fdc7489SMatthew Knepley 2659566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(*dm, &coordinates)); 2669566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 2679566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(*dm, coordSection, coordinates, c, &coordSize, &coords)); 2689371c9d4SSatish Balay for (s = 0; s < coordSize; ++s) 2699371c9d4SSatish Balay for (d = 0; d < dim; ++d) centroid[d] += coords[s * dim + d]; 2709566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(*dm, coordSection, coordinates, c, &coordSize, &coords)); 271*f22e26b7SPierre Jolivet } else PetscCall(DMPlexComputeCellGeometryFVM(*dm, c, nullptr, centroid, nullptr)); 2720fdc7489SMatthew Knepley for (b = 0; b < Nb; ++b) { 2739371c9d4SSatish Balay if (islite) { 2749371c9d4SSatish Balay if (EGlite_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 2759371c9d4SSatish Balay } else { 2769371c9d4SSatish Balay if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 2779371c9d4SSatish Balay } 2780fdc7489SMatthew Knepley } 2790fdc7489SMatthew Knepley if (b < Nb) { 2800fdc7489SMatthew Knepley PetscInt cval = b, eVal, fVal; 281*f22e26b7SPierre Jolivet PetscInt *closure = nullptr, Ncl, cl; 2820fdc7489SMatthew Knepley 2839566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, c, cval)); 2849566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure)); 2850fdc7489SMatthew Knepley for (cl = 0; cl < Ncl; cl += 2) { 2860fdc7489SMatthew Knepley const PetscInt p = closure[cl]; 2870fdc7489SMatthew Knepley 2880fdc7489SMatthew Knepley if (p >= eStart && p < eEnd) { 2899566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, p, &eVal)); 2909566063dSJacob Faibussowitsch if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 2910fdc7489SMatthew Knepley } 2920fdc7489SMatthew Knepley if (p >= fStart && p < fEnd) { 2939566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, p, &fVal)); 2949566063dSJacob Faibussowitsch if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 2950fdc7489SMatthew Knepley } 2960fdc7489SMatthew Knepley } 2979566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure)); 2980fdc7489SMatthew Knepley } 2990fdc7489SMatthew Knepley } 300c1cad2e7SMatthew G. Knepley skip_egads:; 3010fdc7489SMatthew Knepley #endif 3029318fe57SMatthew G. Knepley } 3039566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE)); 3043a074057SBarry Smith } 3059566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelDestroy(&universal)); 3063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3073a074057SBarry Smith } 3083a074057SBarry Smith 309d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexRefine_Tetgen(DM dm, double *maxVolumes, DM *dmRefined) 310d71ae5a4SJacob Faibussowitsch { 3113a074057SBarry Smith MPI_Comm comm; 3123a074057SBarry Smith const PetscInt dim = 3; 3133a074057SBarry Smith ::tetgenio in; 3143a074057SBarry Smith ::tetgenio out; 3159318fe57SMatthew G. Knepley PetscContainer modelObj; 3160fdc7489SMatthew Knepley DMUniversalLabel universal; 317af226901SMatthew G. Knepley PetscInt vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, cStart, cEnd, c, defVal; 3180fdc7489SMatthew Knepley DMPlexInterpolatedFlag isInterpolated; 3193a074057SBarry Smith PetscMPIInt rank; 3203a074057SBarry Smith 3213a074057SBarry Smith PetscFunctionBegin; 3229566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 3249566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(dm, &isInterpolated)); 3259566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelCreate(dm, &universal)); 326af226901SMatthew G. Knepley PetscCall(DMLabelGetDefaultValue(universal->label, &defVal)); 3273a074057SBarry Smith 3289566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 3293a074057SBarry Smith in.numberofpoints = vEnd - vStart; 3303a074057SBarry Smith if (in.numberofpoints > 0) { 3313a074057SBarry Smith PetscSection coordSection; 3323a074057SBarry Smith Vec coordinates; 3333a074057SBarry Smith PetscScalar *array; 3343a074057SBarry Smith 3353a074057SBarry Smith in.pointlist = new double[in.numberofpoints * dim]; 3363a074057SBarry Smith in.pointmarkerlist = new int[in.numberofpoints]; 3373a074057SBarry Smith 33841e9d8b5SMatthew G. Knepley PetscCall(PetscArrayzero(in.pointmarkerlist, (size_t)in.numberofpoints)); 3399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 3409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 3419566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &array)); 3423a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 3433a074057SBarry Smith const PetscInt idx = v - vStart; 3440fdc7489SMatthew Knepley PetscInt off, d, val; 3453a074057SBarry Smith 3469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 3473a074057SBarry Smith for (d = 0; d < dim; ++d) in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]); 3489566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, v, &val)); 349af226901SMatthew G. Knepley if (val != defVal) in.pointmarkerlist[idx] = (int)val; 3503a074057SBarry Smith } 3519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &array)); 3523a074057SBarry Smith } 3533a074057SBarry Smith 3549566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd)); 3550fdc7489SMatthew Knepley in.numberofedges = eEnd - eStart; 3560fdc7489SMatthew Knepley if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberofedges > 0) { 3570fdc7489SMatthew Knepley in.edgelist = new int[in.numberofedges * 2]; 3580fdc7489SMatthew Knepley in.edgemarkerlist = new int[in.numberofedges]; 3590fdc7489SMatthew Knepley for (e = eStart; e < eEnd; ++e) { 3600fdc7489SMatthew Knepley const PetscInt idx = e - eStart; 3610fdc7489SMatthew Knepley const PetscInt *cone; 3620fdc7489SMatthew Knepley PetscInt coneSize, val; 3630fdc7489SMatthew Knepley 3649566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, e, &coneSize)); 3659566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, e, &cone)); 3660fdc7489SMatthew Knepley in.edgelist[idx * 2] = cone[0] - vStart; 3670fdc7489SMatthew Knepley in.edgelist[idx * 2 + 1] = cone[1] - vStart; 3680fdc7489SMatthew Knepley 3699566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, e, &val)); 370af226901SMatthew G. Knepley if (val != defVal) in.edgemarkerlist[idx] = (int)val; 3710fdc7489SMatthew Knepley } 3720fdc7489SMatthew Knepley } 3730fdc7489SMatthew Knepley 3749566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 3750fdc7489SMatthew Knepley in.numberoffacets = fEnd - fStart; 3760fdc7489SMatthew Knepley if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in.numberoffacets > 0) { 3770fdc7489SMatthew Knepley in.facetlist = new tetgenio::facet[in.numberoffacets]; 3780fdc7489SMatthew Knepley in.facetmarkerlist = new int[in.numberoffacets]; 3790fdc7489SMatthew Knepley for (f = fStart; f < fEnd; ++f) { 3800fdc7489SMatthew Knepley const PetscInt idx = f - fStart; 381*f22e26b7SPierre Jolivet PetscInt *points = nullptr, numPoints, p, numVertices = 0, v, val; 3820fdc7489SMatthew Knepley 3830fdc7489SMatthew Knepley in.facetlist[idx].numberofpolygons = 1; 3840fdc7489SMatthew Knepley in.facetlist[idx].polygonlist = new tetgenio::polygon[in.facetlist[idx].numberofpolygons]; 3850fdc7489SMatthew Knepley in.facetlist[idx].numberofholes = 0; 386*f22e26b7SPierre Jolivet in.facetlist[idx].holelist = nullptr; 3870fdc7489SMatthew Knepley 3889566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points)); 3890fdc7489SMatthew Knepley for (p = 0; p < numPoints * 2; p += 2) { 3900fdc7489SMatthew Knepley const PetscInt point = points[p]; 3910fdc7489SMatthew Knepley if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point; 3920fdc7489SMatthew Knepley } 3930fdc7489SMatthew Knepley 3940fdc7489SMatthew Knepley tetgenio::polygon *poly = in.facetlist[idx].polygonlist; 3950fdc7489SMatthew Knepley poly->numberofvertices = numVertices; 3960fdc7489SMatthew Knepley poly->vertexlist = new int[poly->numberofvertices]; 3970fdc7489SMatthew Knepley for (v = 0; v < numVertices; ++v) { 3980fdc7489SMatthew Knepley const PetscInt vIdx = points[v] - vStart; 3990fdc7489SMatthew Knepley poly->vertexlist[v] = vIdx; 4000fdc7489SMatthew Knepley } 4010fdc7489SMatthew Knepley 4029566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, f, &val)); 403af226901SMatthew G. Knepley if (val != defVal) in.facetmarkerlist[idx] = (int)val; 4040fdc7489SMatthew Knepley 4059566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points)); 4060fdc7489SMatthew Knepley } 4070fdc7489SMatthew Knepley } 4080fdc7489SMatthew Knepley 4099566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 4103a074057SBarry Smith in.numberofcorners = 4; 4113a074057SBarry Smith in.numberoftetrahedra = cEnd - cStart; 4123a074057SBarry Smith in.tetrahedronvolumelist = (double *)maxVolumes; 4133a074057SBarry Smith if (in.numberoftetrahedra > 0) { 4143a074057SBarry Smith in.tetrahedronlist = new int[in.numberoftetrahedra * in.numberofcorners]; 4153a074057SBarry Smith for (c = cStart; c < cEnd; ++c) { 4163a074057SBarry Smith const PetscInt idx = c - cStart; 417*f22e26b7SPierre Jolivet PetscInt *closure = nullptr; 4183a074057SBarry Smith PetscInt closureSize; 4193a074057SBarry Smith 4209566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 42163a3b9bcSJacob 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); 4220fdc7489SMatthew Knepley for (v = 0; v < 4; ++v) in.tetrahedronlist[idx * in.numberofcorners + v] = closure[(v + closureSize - 4) * 2] - vStart; 4239566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 4243a074057SBarry Smith } 4253a074057SBarry Smith } 4260fdc7489SMatthew Knepley 427dd400576SPatrick Sanan if (rank == 0) { 4283a074057SBarry Smith char args[32]; 4293a074057SBarry Smith 4303a074057SBarry Smith /* Take away 'Q' for verbose output */ 431c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(args, "qezQra", sizeof(args))); 4323a074057SBarry Smith ::tetrahedralize(args, &in, &out); 4333a074057SBarry Smith } 4343a074057SBarry Smith 435*f22e26b7SPierre Jolivet in.tetrahedronvolumelist = nullptr; 4363a074057SBarry Smith { 4373a074057SBarry Smith const PetscInt numCorners = 4; 4383a074057SBarry Smith const PetscInt numCells = out.numberoftetrahedra; 4393a074057SBarry Smith const PetscInt numVertices = out.numberofpoints; 440*f22e26b7SPierre Jolivet PetscReal *meshCoords = nullptr; 441*f22e26b7SPierre Jolivet PetscInt *cells = nullptr; 4420fdc7489SMatthew Knepley PetscBool interpolate = isInterpolated == DMPLEX_INTERPOLATED_FULL ? PETSC_TRUE : PETSC_FALSE; 4433a074057SBarry Smith 444a4a685f2SJacob Faibussowitsch if (sizeof(PetscReal) == sizeof(out.pointlist[0])) { 445a4a685f2SJacob Faibussowitsch meshCoords = (PetscReal *)out.pointlist; 446a4a685f2SJacob Faibussowitsch } else { 447a4a685f2SJacob Faibussowitsch PetscInt i; 448a4a685f2SJacob Faibussowitsch 449a4a685f2SJacob Faibussowitsch meshCoords = new PetscReal[dim * numVertices]; 4500fdc7489SMatthew Knepley for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal)out.pointlist[i]; 451a4a685f2SJacob Faibussowitsch } 452a4a685f2SJacob Faibussowitsch if (sizeof(PetscInt) == sizeof(out.tetrahedronlist[0])) { 453a4a685f2SJacob Faibussowitsch cells = (PetscInt *)out.tetrahedronlist; 454a4a685f2SJacob Faibussowitsch } else { 455a4a685f2SJacob Faibussowitsch PetscInt i; 456a4a685f2SJacob Faibussowitsch 457a4a685f2SJacob Faibussowitsch cells = new PetscInt[numCells * numCorners]; 4580fdc7489SMatthew Knepley for (i = 0; i < numCells * numCorners; ++i) cells[i] = (PetscInt)out.tetrahedronlist[i]; 459a4a685f2SJacob Faibussowitsch } 460a4a685f2SJacob Faibussowitsch 4619566063dSJacob Faibussowitsch PetscCall(DMPlexInvertCells_Tetgen(numCells, numCorners, cells)); 4629566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined)); 463ad540459SPierre Jolivet if (sizeof(PetscReal) != sizeof(out.pointlist[0])) delete[] meshCoords; 464ad540459SPierre Jolivet if (sizeof(PetscInt) != sizeof(out.tetrahedronlist[0])) delete[] cells; 4650fdc7489SMatthew Knepley 4663a074057SBarry Smith /* Set labels */ 4679566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dmRefined)); 4683a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 46948a46eb9SPierre Jolivet if (out.pointmarkerlist[v]) PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, v + numCells, out.pointmarkerlist[v])); 4703a074057SBarry Smith } 4713a074057SBarry Smith if (interpolate) { 4720fdc7489SMatthew Knepley PetscInt e, f; 4733a074057SBarry Smith 4740fdc7489SMatthew Knepley for (e = 0; e < out.numberofedges; ++e) { 4753a074057SBarry Smith if (out.edgemarkerlist[e]) { 4763a074057SBarry Smith const PetscInt vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells}; 4773a074057SBarry Smith const PetscInt *edges; 4783a074057SBarry Smith PetscInt numEdges; 4793a074057SBarry Smith 4809566063dSJacob Faibussowitsch PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 48163a3b9bcSJacob Faibussowitsch PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges); 4829566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, edges[0], out.edgemarkerlist[e])); 4839566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 4843a074057SBarry Smith } 4853a074057SBarry Smith } 4860fdc7489SMatthew Knepley for (f = 0; f < out.numberoftrifaces; ++f) { 4873a074057SBarry Smith if (out.trifacemarkerlist[f]) { 4883a074057SBarry Smith const PetscInt vertices[3] = {out.trifacelist[f * 3 + 0] + numCells, out.trifacelist[f * 3 + 1] + numCells, out.trifacelist[f * 3 + 2] + numCells}; 4893a074057SBarry Smith const PetscInt *faces; 4903a074057SBarry Smith PetscInt numFaces; 4913a074057SBarry Smith 4929566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces)); 49363a3b9bcSJacob Faibussowitsch PetscCheck(numFaces == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces); 4949566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, faces[0], out.trifacemarkerlist[f])); 4959566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces)); 4963a074057SBarry Smith } 4973a074057SBarry Smith } 4983a074057SBarry Smith } 4990fdc7489SMatthew Knepley 5009566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj)); 5019318fe57SMatthew G. Knepley if (modelObj) { 5020fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS 5030fdc7489SMatthew Knepley DMLabel bodyLabel; 5040fdc7489SMatthew Knepley PetscInt cStart, cEnd, c, eStart, eEnd, fStart, fEnd; 505c1cad2e7SMatthew G. Knepley PetscBool islite = PETSC_FALSE; 5060fdc7489SMatthew Knepley ego *bodies; 5070fdc7489SMatthew Knepley ego model, geom; 5080fdc7489SMatthew Knepley int Nb, oclass, mtype, *senses; 5090fdc7489SMatthew Knepley 5100fdc7489SMatthew Knepley /* Get Attached EGADS Model from Original DMPlex */ 5119566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj)); 512c1cad2e7SMatthew G. Knepley if (modelObj) { 5139566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(modelObj, (void **)&model)); 514*f22e26b7SPierre Jolivet PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, nullptr, &Nb, &bodies, &senses)); 5150fdc7489SMatthew Knepley /* Transfer EGADS Model to Volumetric Mesh */ 5169566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)*dmRefined, "EGADS Model", (PetscObject)modelObj)); 517c1cad2e7SMatthew G. Knepley } else { 5189566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSLite Model", (PetscObject *)&modelObj)); 519c1cad2e7SMatthew G. Knepley if (modelObj) { 5209566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(modelObj, (void **)&model)); 521*f22e26b7SPierre Jolivet PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, nullptr, &Nb, &bodies, &senses)); 522c1cad2e7SMatthew G. Knepley /* Transfer EGADS Model to Volumetric Mesh */ 5239566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)*dmRefined, "EGADSLite Model", (PetscObject)modelObj)); 524c1cad2e7SMatthew G. Knepley islite = PETSC_TRUE; 525c1cad2e7SMatthew G. Knepley } 526c1cad2e7SMatthew G. Knepley } 527c1cad2e7SMatthew G. Knepley if (!modelObj) goto skip_egads; 5280fdc7489SMatthew Knepley 5290fdc7489SMatthew Knepley /* Set Cell Labels */ 5309566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dmRefined, "EGADS Body ID", &bodyLabel)); 5319566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dmRefined, 0, &cStart, &cEnd)); 5329566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dmRefined, 1, &fStart, &fEnd)); 5339566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(*dmRefined, 1, &eStart, &eEnd)); 5340fdc7489SMatthew Knepley 5350fdc7489SMatthew Knepley for (c = cStart; c < cEnd; ++c) { 5360fdc7489SMatthew Knepley PetscReal centroid[3] = {0., 0., 0.}; 5370fdc7489SMatthew Knepley PetscInt b; 5380fdc7489SMatthew Knepley 53935cb6cd3SPierre Jolivet /* Determine what body the cell's centroid is located in */ 5400fdc7489SMatthew Knepley if (!interpolate) { 5410fdc7489SMatthew Knepley PetscSection coordSection; 5420fdc7489SMatthew Knepley Vec coordinates; 543*f22e26b7SPierre Jolivet PetscScalar *coords = nullptr; 5440fdc7489SMatthew Knepley PetscInt coordSize, s, d; 5450fdc7489SMatthew Knepley 5469566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(*dmRefined, &coordinates)); 5479566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dmRefined, &coordSection)); 5489566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords)); 5499371c9d4SSatish Balay for (s = 0; s < coordSize; ++s) 5509371c9d4SSatish Balay for (d = 0; d < dim; ++d) centroid[d] += coords[s * dim + d]; 5519566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords)); 552*f22e26b7SPierre Jolivet } else PetscCall(DMPlexComputeCellGeometryFVM(*dmRefined, c, nullptr, centroid, nullptr)); 5530fdc7489SMatthew Knepley for (b = 0; b < Nb; ++b) { 5549371c9d4SSatish Balay if (islite) { 5559371c9d4SSatish Balay if (EGlite_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 5569371c9d4SSatish Balay } else { 5579371c9d4SSatish Balay if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 5589371c9d4SSatish Balay } 5590fdc7489SMatthew Knepley } 5600fdc7489SMatthew Knepley if (b < Nb) { 5610fdc7489SMatthew Knepley PetscInt cval = b, eVal, fVal; 562*f22e26b7SPierre Jolivet PetscInt *closure = nullptr, Ncl, cl; 5630fdc7489SMatthew Knepley 5649566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, c, cval)); 5659566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure)); 5660fdc7489SMatthew Knepley for (cl = 0; cl < Ncl; cl += 2) { 5670fdc7489SMatthew Knepley const PetscInt p = closure[cl]; 5680fdc7489SMatthew Knepley 5690fdc7489SMatthew Knepley if (p >= eStart && p < eEnd) { 5709566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, p, &eVal)); 5719566063dSJacob Faibussowitsch if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 5720fdc7489SMatthew Knepley } 5730fdc7489SMatthew Knepley if (p >= fStart && p < fEnd) { 5749566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, p, &fVal)); 5759566063dSJacob Faibussowitsch if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 5760fdc7489SMatthew Knepley } 5770fdc7489SMatthew Knepley } 5789566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure)); 5790fdc7489SMatthew Knepley } 5800fdc7489SMatthew Knepley } 581c1cad2e7SMatthew G. Knepley skip_egads:; 5820fdc7489SMatthew Knepley #endif 5839318fe57SMatthew G. Knepley } 5849566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE)); 5853a074057SBarry Smith } 5863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5873a074057SBarry Smith } 588