13a074057SBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 23a074057SBarry Smith 30fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS 40fdc7489SMatthew Knepley #include <egads.h> 50fdc7489SMatthew Knepley #endif 60fdc7489SMatthew Knepley 73a074057SBarry Smith #include <ctetgen.h> 83a074057SBarry Smith 93a074057SBarry Smith /* This is to fix the tetrahedron orientation from TetGen */ 109371c9d4SSatish Balay static PetscErrorCode DMPlexInvertCells_CTetgen(PetscInt numCells, PetscInt numCorners, PetscInt cells[]) { 113a074057SBarry Smith PetscInt bound = numCells * numCorners, coff; 123a074057SBarry Smith 133a074057SBarry Smith PetscFunctionBegin; 149371c9d4SSatish Balay #define SWAP(a, b) \ 159371c9d4SSatish Balay do { \ 169371c9d4SSatish Balay PetscInt tmp = (a); \ 179371c9d4SSatish Balay (a) = (b); \ 189371c9d4SSatish Balay (b) = tmp; \ 199371c9d4SSatish Balay } while (0) 2096ca5757SLisandro Dalcin for (coff = 0; coff < bound; coff += numCorners) SWAP(cells[coff], cells[coff + 1]); 2196ca5757SLisandro Dalcin #undef SWAP 223a074057SBarry Smith PetscFunctionReturn(0); 233a074057SBarry Smith } 243a074057SBarry Smith 259371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexGenerate_CTetgen(DM boundary, PetscBool interpolate, DM *dm) { 263a074057SBarry Smith MPI_Comm comm; 273a074057SBarry Smith const PetscInt dim = 3; 283a074057SBarry Smith PLC *in, *out; 290fdc7489SMatthew Knepley DMUniversalLabel universal; 30af226901SMatthew G. Knepley PetscInt vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, defVal, verbose = 0; 310fdc7489SMatthew Knepley DMPlexInterpolatedFlag isInterpolated; 323a074057SBarry Smith PetscMPIInt rank; 333a074057SBarry Smith 343a074057SBarry Smith PetscFunctionBegin; 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)boundary)->prefix, "-ctetgen_verbose", &verbose, NULL)); 369566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)boundary, &comm)); 379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 389566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(boundary, &isInterpolated)); 399566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelCreate(boundary, &universal)); 40af226901SMatthew G. Knepley PetscCall(DMLabelGetDefaultValue(universal->label, &defVal)); 410fdc7489SMatthew Knepley 429566063dSJacob Faibussowitsch PetscCall(PLCCreate(&in)); 439566063dSJacob Faibussowitsch PetscCall(PLCCreate(&out)); 443a074057SBarry Smith 459566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd)); 463a074057SBarry Smith in->numberofpoints = vEnd - vStart; 473a074057SBarry Smith if (in->numberofpoints > 0) { 483a074057SBarry Smith PetscSection coordSection; 493a074057SBarry Smith Vec coordinates; 500fdc7489SMatthew Knepley const PetscScalar *array; 513a074057SBarry Smith 529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->numberofpoints * dim, &in->pointlist)); 5341e9d8b5SMatthew G. Knepley PetscCall(PetscCalloc1(in->numberofpoints, &in->pointmarkerlist)); 549566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(boundary, &coordinates)); 559566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(boundary, &coordSection)); 569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &array)); 573a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 583a074057SBarry Smith const PetscInt idx = v - vStart; 590fdc7489SMatthew Knepley PetscInt off, d, m; 603a074057SBarry Smith 619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 620fdc7489SMatthew Knepley for (d = 0; d < dim; ++d) in->pointlist[idx * dim + d] = PetscRealPart(array[off + d]); 639566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, v, &m)); 64af226901SMatthew G. Knepley if (m != defVal) in->pointmarkerlist[idx] = (int)m; 653a074057SBarry Smith } 669566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &array)); 673a074057SBarry Smith } 683a074057SBarry Smith 699566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(boundary, 1, &eStart, &eEnd)); 700fdc7489SMatthew Knepley in->numberofedges = eEnd - eStart; 710fdc7489SMatthew Knepley if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberofedges > 0) { 729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->numberofedges * 2, &in->edgelist)); 739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->numberofedges, &in->edgemarkerlist)); 740fdc7489SMatthew Knepley for (e = eStart; e < eEnd; ++e) { 750fdc7489SMatthew Knepley const PetscInt idx = e - eStart; 760fdc7489SMatthew Knepley const PetscInt *cone; 770fdc7489SMatthew Knepley PetscInt coneSize, val; 780fdc7489SMatthew Knepley 799566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(boundary, e, &coneSize)); 809566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(boundary, e, &cone)); 810fdc7489SMatthew Knepley in->edgelist[idx * 2] = cone[0] - vStart; 820fdc7489SMatthew Knepley in->edgelist[idx * 2 + 1] = cone[1] - vStart; 830fdc7489SMatthew Knepley 849566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, e, &val)); 85af226901SMatthew G. Knepley if (val != defVal) in->edgemarkerlist[idx] = (int)val; 860fdc7489SMatthew Knepley } 870fdc7489SMatthew Knepley } 880fdc7489SMatthew Knepley 899566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd)); 903a074057SBarry Smith in->numberoffacets = fEnd - fStart; 913a074057SBarry Smith if (in->numberoffacets > 0) { 929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->numberoffacets, &in->facetlist)); 939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->numberoffacets, &in->facetmarkerlist)); 943a074057SBarry Smith for (f = fStart; f < fEnd; ++f) { 953a074057SBarry Smith const PetscInt idx = f - fStart; 963a074057SBarry Smith PetscInt *points = NULL, numPoints, p, numVertices = 0, v, m = -1; 973a074057SBarry Smith polygon *poly; 983a074057SBarry Smith 993a074057SBarry Smith in->facetlist[idx].numberofpolygons = 1; 1009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->facetlist[idx].numberofpolygons, &in->facetlist[idx].polygonlist)); 1013a074057SBarry Smith in->facetlist[idx].numberofholes = 0; 1023a074057SBarry Smith in->facetlist[idx].holelist = NULL; 1033a074057SBarry Smith 1049566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points)); 1053a074057SBarry Smith for (p = 0; p < numPoints * 2; p += 2) { 1063a074057SBarry Smith const PetscInt point = points[p]; 1073a074057SBarry Smith if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point; 1083a074057SBarry Smith } 1093a074057SBarry Smith 1103a074057SBarry Smith poly = in->facetlist[idx].polygonlist; 1113a074057SBarry Smith poly->numberofvertices = numVertices; 1129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(poly->numberofvertices, &poly->vertexlist)); 1133a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 1143a074057SBarry Smith const PetscInt vIdx = points[v] - vStart; 1153a074057SBarry Smith poly->vertexlist[v] = vIdx; 1163a074057SBarry Smith } 1179566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, f, &m)); 118af226901SMatthew G. Knepley if (m != defVal) in->facetmarkerlist[idx] = (int)m; 1199566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points)); 1203a074057SBarry Smith } 1213a074057SBarry Smith } 122dd400576SPatrick Sanan if (rank == 0) { 1233a074057SBarry Smith TetGenOpts t; 1243a074057SBarry Smith 1259566063dSJacob Faibussowitsch PetscCall(TetGenOptsInitialize(&t)); 1263a074057SBarry Smith t.in = boundary; /* Should go away */ 1273a074057SBarry Smith t.plc = 1; 1283a074057SBarry Smith t.quality = 1; 1293a074057SBarry Smith t.edgesout = 1; 1303a074057SBarry Smith t.zeroindex = 1; 1313a074057SBarry Smith t.quiet = 1; 1323a074057SBarry Smith t.verbose = verbose; 1330fdc7489SMatthew Knepley #if 0 1340fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS 1350fdc7489SMatthew Knepley /* Need to add in more TetGen code */ 1360fdc7489SMatthew Knepley t.nobisect = 1; /* Add Y to preserve Surface Mesh for EGADS */ 1370fdc7489SMatthew Knepley #endif 1380fdc7489SMatthew Knepley #endif 1390fdc7489SMatthew Knepley 1409566063dSJacob Faibussowitsch PetscCall(TetGenCheckOpts(&t)); 1419566063dSJacob Faibussowitsch PetscCall(TetGenTetrahedralize(&t, in, out)); 1423a074057SBarry Smith } 1433a074057SBarry Smith { 1443a074057SBarry Smith const PetscInt numCorners = 4; 1453a074057SBarry Smith const PetscInt numCells = out->numberoftetrahedra; 1463a074057SBarry Smith const PetscInt numVertices = out->numberofpoints; 1470fdc7489SMatthew Knepley PetscReal *meshCoords = NULL; 1480fdc7489SMatthew Knepley PetscInt *cells = NULL; 1493a074057SBarry Smith 150a4a685f2SJacob Faibussowitsch if (sizeof(PetscReal) == sizeof(out->pointlist[0])) { 151a4a685f2SJacob Faibussowitsch meshCoords = (PetscReal *)out->pointlist; 152a4a685f2SJacob Faibussowitsch } else { 1533a074057SBarry Smith PetscInt i; 1543a074057SBarry Smith 1559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * numVertices, &meshCoords)); 1560fdc7489SMatthew Knepley for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal)out->pointlist[i]; 1573a074057SBarry Smith } 158a4a685f2SJacob Faibussowitsch if (sizeof(PetscInt) == sizeof(out->tetrahedronlist[0])) { 159a4a685f2SJacob Faibussowitsch cells = (PetscInt *)out->tetrahedronlist; 160a4a685f2SJacob Faibussowitsch } else { 161a4a685f2SJacob Faibussowitsch PetscInt i; 162a4a685f2SJacob Faibussowitsch 1639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells * numCorners, &cells)); 1640fdc7489SMatthew Knepley for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out->tetrahedronlist[i]; 165a4a685f2SJacob Faibussowitsch } 1663a074057SBarry Smith 1679566063dSJacob Faibussowitsch PetscCall(DMPlexInvertCells_CTetgen(numCells, numCorners, cells)); 1689566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm)); 169*48a46eb9SPierre Jolivet if (sizeof(PetscReal) != sizeof(out->pointlist[0])) PetscCall(PetscFree(meshCoords)); 170*48a46eb9SPierre Jolivet if (sizeof(PetscInt) != sizeof(out->tetrahedronlist[0])) PetscCall(PetscFree(cells)); 1710fdc7489SMatthew Knepley 1723a074057SBarry Smith /* Set labels */ 1739566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dm)); 1743a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 175*48a46eb9SPierre Jolivet if (out->pointmarkerlist[v]) PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, v + numCells, out->pointmarkerlist[v])); 1763a074057SBarry Smith } 1773a074057SBarry Smith if (interpolate) { 1783a074057SBarry Smith PetscInt e; 1793a074057SBarry Smith 1803a074057SBarry Smith for (e = 0; e < out->numberofedges; e++) { 1813a074057SBarry Smith if (out->edgemarkerlist[e]) { 1823a074057SBarry Smith const PetscInt vertices[2] = {out->edgelist[e * 2 + 0] + numCells, out->edgelist[e * 2 + 1] + numCells}; 1833a074057SBarry Smith const PetscInt *edges; 1843a074057SBarry Smith PetscInt numEdges; 1853a074057SBarry Smith 1869566063dSJacob Faibussowitsch PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges)); 18763a3b9bcSJacob Faibussowitsch PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges); 1889566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, edges[0], out->edgemarkerlist[e])); 1899566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges)); 1903a074057SBarry Smith } 1913a074057SBarry Smith } 1923a074057SBarry Smith for (f = 0; f < out->numberoftrifaces; f++) { 1933a074057SBarry Smith if (out->trifacemarkerlist[f]) { 1943a074057SBarry Smith const PetscInt vertices[3] = {out->trifacelist[f * 3 + 0] + numCells, out->trifacelist[f * 3 + 1] + numCells, out->trifacelist[f * 3 + 2] + numCells}; 1953a074057SBarry Smith const PetscInt *faces; 1963a074057SBarry Smith PetscInt numFaces; 1973a074057SBarry Smith 1989566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces)); 19963a3b9bcSJacob Faibussowitsch PetscCheck(numFaces == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces); 2009566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, faces[0], out->trifacemarkerlist[f])); 2019566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces)); 2023a074057SBarry Smith } 2033a074057SBarry Smith } 2043a074057SBarry Smith } 2050fdc7489SMatthew Knepley 2060fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS 2070fdc7489SMatthew Knepley { 2080fdc7489SMatthew Knepley DMLabel bodyLabel; 2090fdc7489SMatthew Knepley PetscContainer modelObj; 2100fdc7489SMatthew Knepley PetscInt cStart, cEnd, c, eStart, eEnd, fStart, fEnd; 2110fdc7489SMatthew Knepley ego *bodies; 2120fdc7489SMatthew Knepley ego model, geom; 2130fdc7489SMatthew Knepley int Nb, oclass, mtype, *senses; 2140fdc7489SMatthew Knepley 2150fdc7489SMatthew Knepley /* Get Attached EGADS Model from Original DMPlex */ 2169566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)boundary, "EGADS Model", (PetscObject *)&modelObj)); 2170fdc7489SMatthew Knepley if (modelObj) { 2189566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(modelObj, (void **)&model)); 2199566063dSJacob Faibussowitsch PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 2200fdc7489SMatthew Knepley /* Transfer EGADS Model to Volumetric Mesh */ 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)*dm, "EGADS Model", (PetscObject)modelObj)); 2220fdc7489SMatthew Knepley 2230fdc7489SMatthew Knepley /* Set Cell Labels */ 2249566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "EGADS Body ID", &bodyLabel)); 2259566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 2269566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd)); 2279566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(*dm, 1, &eStart, &eEnd)); 2280fdc7489SMatthew Knepley 2290fdc7489SMatthew Knepley for (c = cStart; c < cEnd; ++c) { 2300fdc7489SMatthew Knepley PetscReal centroid[3] = {0., 0., 0.}; 2310fdc7489SMatthew Knepley PetscInt b; 2320fdc7489SMatthew Knepley 2330fdc7489SMatthew Knepley /* Deterimine what body the cell's centroid is located in */ 2340fdc7489SMatthew Knepley if (!interpolate) { 2350fdc7489SMatthew Knepley PetscSection coordSection; 2360fdc7489SMatthew Knepley Vec coordinates; 2370fdc7489SMatthew Knepley PetscScalar *coords = NULL; 2380fdc7489SMatthew Knepley PetscInt coordSize, s, d; 2390fdc7489SMatthew Knepley 2409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(*dm, &coordinates)); 2419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 2429566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(*dm, coordSection, coordinates, c, &coordSize, &coords)); 2439371c9d4SSatish Balay for (s = 0; s < coordSize; ++s) 2449371c9d4SSatish Balay for (d = 0; d < dim; ++d) centroid[d] += coords[s * dim + d]; 2459566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(*dm, coordSection, coordinates, c, &coordSize, &coords)); 2461baa6e33SBarry Smith } else PetscCall(DMPlexComputeCellGeometryFVM(*dm, c, NULL, centroid, NULL)); 2470fdc7489SMatthew Knepley for (b = 0; b < Nb; ++b) { 2480fdc7489SMatthew Knepley if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 2490fdc7489SMatthew Knepley } 2500fdc7489SMatthew Knepley if (b < Nb) { 2510fdc7489SMatthew Knepley PetscInt cval = b, eVal, fVal; 2520fdc7489SMatthew Knepley PetscInt *closure = NULL, Ncl, cl; 2530fdc7489SMatthew Knepley 2549566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, c, cval)); 2559566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure)); 256c1cad2e7SMatthew G. Knepley for (cl = 0; cl < Ncl; ++cl) { 257c1cad2e7SMatthew G. Knepley const PetscInt p = closure[cl * 2]; 2580fdc7489SMatthew Knepley 2590fdc7489SMatthew Knepley if (p >= eStart && p < eEnd) { 2609566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, p, &eVal)); 2619566063dSJacob Faibussowitsch if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 2620fdc7489SMatthew Knepley } 2630fdc7489SMatthew Knepley if (p >= fStart && p < fEnd) { 2649566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, p, &fVal)); 2659566063dSJacob Faibussowitsch if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 2660fdc7489SMatthew Knepley } 2670fdc7489SMatthew Knepley } 2689566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure)); 2690fdc7489SMatthew Knepley } 2700fdc7489SMatthew Knepley } 2710fdc7489SMatthew Knepley } 2720fdc7489SMatthew Knepley } 2730fdc7489SMatthew Knepley #endif 2749566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE)); 2753a074057SBarry Smith } 2763a074057SBarry Smith 2779566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelDestroy(&universal)); 2789566063dSJacob Faibussowitsch PetscCall(PLCDestroy(&in)); 2799566063dSJacob Faibussowitsch PetscCall(PLCDestroy(&out)); 2803a074057SBarry Smith PetscFunctionReturn(0); 2813a074057SBarry Smith } 2823a074057SBarry Smith 2839371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined) { 2843a074057SBarry Smith MPI_Comm comm; 2853a074057SBarry Smith const PetscInt dim = 3; 2863a074057SBarry Smith PLC *in, *out; 2870fdc7489SMatthew Knepley DMUniversalLabel universal; 288af226901SMatthew G. Knepley PetscInt vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, cStart, cEnd, c, defVal, verbose = 0; 2890fdc7489SMatthew Knepley DMPlexInterpolatedFlag isInterpolated; 2903a074057SBarry Smith PetscMPIInt rank; 2913a074057SBarry Smith 2923a074057SBarry Smith PetscFunctionBegin; 2939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)dm)->prefix, "-ctetgen_verbose", &verbose, NULL)); 2949566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2969566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(dm, &isInterpolated)); 2979566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelCreate(dm, &universal)); 298af226901SMatthew G. Knepley PetscCall(DMLabelGetDefaultValue(universal->label, &defVal)); 2990fdc7489SMatthew Knepley 3009566063dSJacob Faibussowitsch PetscCall(PLCCreate(&in)); 3019566063dSJacob Faibussowitsch PetscCall(PLCCreate(&out)); 3023a074057SBarry Smith 3039566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 3043a074057SBarry Smith in->numberofpoints = vEnd - vStart; 3053a074057SBarry Smith if (in->numberofpoints > 0) { 3063a074057SBarry Smith PetscSection coordSection; 3073a074057SBarry Smith Vec coordinates; 3083a074057SBarry Smith PetscScalar *array; 3093a074057SBarry Smith 3109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->numberofpoints * dim, &in->pointlist)); 31141e9d8b5SMatthew G. Knepley PetscCall(PetscCalloc1(in->numberofpoints, &in->pointmarkerlist)); 3129566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 3139566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 3149566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &array)); 3153a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 3163a074057SBarry Smith const PetscInt idx = v - vStart; 3170fdc7489SMatthew Knepley PetscInt off, d, m; 3183a074057SBarry Smith 3199566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 3200fdc7489SMatthew Knepley for (d = 0; d < dim; ++d) in->pointlist[idx * dim + d] = PetscRealPart(array[off + d]); 3219566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, v, &m)); 322af226901SMatthew G. Knepley if (m != defVal) in->pointmarkerlist[idx] = (int)m; 3233a074057SBarry Smith } 3249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &array)); 3253a074057SBarry Smith } 3263a074057SBarry Smith 3279566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd)); 3280fdc7489SMatthew Knepley in->numberofedges = eEnd - eStart; 3290fdc7489SMatthew Knepley if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberofedges > 0) { 3309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->numberofedges * 2, &in->edgelist)); 3319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->numberofedges, &in->edgemarkerlist)); 3320fdc7489SMatthew Knepley for (e = eStart; e < eEnd; ++e) { 3330fdc7489SMatthew Knepley const PetscInt idx = e - eStart; 3340fdc7489SMatthew Knepley const PetscInt *cone; 3350fdc7489SMatthew Knepley PetscInt coneSize, val; 3360fdc7489SMatthew Knepley 3379566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, e, &coneSize)); 3389566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, e, &cone)); 3390fdc7489SMatthew Knepley in->edgelist[idx * 2] = cone[0] - vStart; 3400fdc7489SMatthew Knepley in->edgelist[idx * 2 + 1] = cone[1] - vStart; 3410fdc7489SMatthew Knepley 3429566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, e, &val)); 343af226901SMatthew G. Knepley if (val != defVal) in->edgemarkerlist[idx] = (int)val; 3440fdc7489SMatthew Knepley } 3450fdc7489SMatthew Knepley } 3460fdc7489SMatthew Knepley 3479566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 3480fdc7489SMatthew Knepley in->numberoftrifaces = 0; 3490fdc7489SMatthew Knepley for (f = fStart; f < fEnd; ++f) { 3500fdc7489SMatthew Knepley PetscInt supportSize; 3510fdc7489SMatthew Knepley 3529566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, f, &supportSize)); 3530fdc7489SMatthew Knepley if (supportSize == 1) ++in->numberoftrifaces; 3540fdc7489SMatthew Knepley } 3550fdc7489SMatthew Knepley if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberoftrifaces > 0) { 3560fdc7489SMatthew Knepley PetscInt tf = 0; 3570fdc7489SMatthew Knepley 3589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->numberoftrifaces * 3, &in->trifacelist)); 3599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->numberoftrifaces, &in->trifacemarkerlist)); 3600fdc7489SMatthew Knepley for (f = fStart; f < fEnd; ++f) { 3610fdc7489SMatthew Knepley PetscInt *points = NULL; 3620fdc7489SMatthew Knepley PetscInt supportSize, numPoints, p, Nv = 0, val; 3630fdc7489SMatthew Knepley 3649566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, f, &supportSize)); 3650fdc7489SMatthew Knepley if (supportSize != 1) continue; 3669566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points)); 3670fdc7489SMatthew Knepley for (p = 0; p < numPoints * 2; p += 2) { 3680fdc7489SMatthew Knepley const PetscInt point = points[p]; 3690fdc7489SMatthew Knepley if ((point >= vStart) && (point < vEnd)) in->trifacelist[tf * 3 + Nv++] = point - vStart; 3700fdc7489SMatthew Knepley } 3719566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points)); 37263a3b9bcSJacob Faibussowitsch PetscCheck(Nv == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " has %" PetscInt_FMT " vertices, not 3", f, Nv); 3739566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(universal->label, f, &val)); 374af226901SMatthew G. Knepley if (val != defVal) in->trifacemarkerlist[tf] = (int)val; 3750fdc7489SMatthew Knepley ++tf; 3760fdc7489SMatthew Knepley } 3770fdc7489SMatthew Knepley } 3780fdc7489SMatthew Knepley 3799566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 3803a074057SBarry Smith in->numberofcorners = 4; 3813a074057SBarry Smith in->numberoftetrahedra = cEnd - cStart; 3823a074057SBarry Smith in->tetrahedronvolumelist = maxVolumes; 3833a074057SBarry Smith if (in->numberoftetrahedra > 0) { 3849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in->numberoftetrahedra * in->numberofcorners, &in->tetrahedronlist)); 3853a074057SBarry Smith for (c = cStart; c < cEnd; ++c) { 3863a074057SBarry Smith const PetscInt idx = c - cStart; 3873a074057SBarry Smith PetscInt *closure = NULL; 3883a074057SBarry Smith PetscInt closureSize; 3893a074057SBarry Smith 3909566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 3915f80ce2aSJacob 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); 3920fdc7489SMatthew Knepley for (v = 0; v < 4; ++v) in->tetrahedronlist[idx * in->numberofcorners + v] = closure[(v + closureSize - 4) * 2] - vStart; 3939566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 3943a074057SBarry Smith } 3953a074057SBarry Smith } 3960fdc7489SMatthew Knepley 397dd400576SPatrick Sanan if (rank == 0) { 3983a074057SBarry Smith TetGenOpts t; 3993a074057SBarry Smith 4009566063dSJacob Faibussowitsch PetscCall(TetGenOptsInitialize(&t)); 4013a074057SBarry Smith 4023a074057SBarry Smith t.in = dm; /* Should go away */ 4033a074057SBarry Smith t.refine = 1; 4043a074057SBarry Smith t.varvolume = 1; 4053a074057SBarry Smith t.quality = 1; 4063a074057SBarry Smith t.edgesout = 1; 4073a074057SBarry Smith t.zeroindex = 1; 4083a074057SBarry Smith t.quiet = 1; 4093a074057SBarry Smith t.verbose = verbose; /* Change this */ 4103a074057SBarry Smith 4119566063dSJacob Faibussowitsch PetscCall(TetGenCheckOpts(&t)); 4129566063dSJacob Faibussowitsch PetscCall(TetGenTetrahedralize(&t, in, out)); 4133a074057SBarry Smith } 4140fdc7489SMatthew Knepley 4153a074057SBarry Smith in->tetrahedronvolumelist = NULL; 4163a074057SBarry Smith { 4173a074057SBarry Smith const PetscInt numCorners = 4; 4183a074057SBarry Smith const PetscInt numCells = out->numberoftetrahedra; 4193a074057SBarry Smith const PetscInt numVertices = out->numberofpoints; 4200fdc7489SMatthew Knepley PetscReal *meshCoords = NULL; 4210fdc7489SMatthew Knepley PetscInt *cells = NULL; 4220fdc7489SMatthew Knepley PetscBool interpolate = isInterpolated == DMPLEX_INTERPOLATED_FULL ? PETSC_TRUE : PETSC_FALSE; 4233a074057SBarry Smith 424a4a685f2SJacob Faibussowitsch if (sizeof(PetscReal) == sizeof(out->pointlist[0])) { 425a4a685f2SJacob Faibussowitsch meshCoords = (PetscReal *)out->pointlist; 426a4a685f2SJacob Faibussowitsch } else { 4273a074057SBarry Smith PetscInt i; 4283a074057SBarry Smith 4299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * numVertices, &meshCoords)); 4300fdc7489SMatthew Knepley for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal)out->pointlist[i]; 4313a074057SBarry Smith } 432a4a685f2SJacob Faibussowitsch if (sizeof(PetscInt) == sizeof(out->tetrahedronlist[0])) { 433a4a685f2SJacob Faibussowitsch cells = (PetscInt *)out->tetrahedronlist; 434a4a685f2SJacob Faibussowitsch } else { 435a4a685f2SJacob Faibussowitsch PetscInt i; 436a4a685f2SJacob Faibussowitsch 4379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells * numCorners, &cells)); 4380fdc7489SMatthew Knepley for (i = 0; i < numCells * numCorners; ++i) cells[i] = (PetscInt)out->tetrahedronlist[i]; 439a4a685f2SJacob Faibussowitsch } 4403a074057SBarry Smith 4419566063dSJacob Faibussowitsch PetscCall(DMPlexInvertCells_CTetgen(numCells, numCorners, cells)); 4429566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined)); 4439566063dSJacob Faibussowitsch if (sizeof(PetscReal) != sizeof(out->pointlist[0])) PetscCall(PetscFree(meshCoords)); 4449566063dSJacob Faibussowitsch if (sizeof(PetscInt) != sizeof(out->tetrahedronlist[0])) PetscCall(PetscFree(cells)); 4450fdc7489SMatthew Knepley 4463a074057SBarry Smith /* Set labels */ 4479566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dmRefined)); 4483a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 449*48a46eb9SPierre Jolivet if (out->pointmarkerlist[v]) PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, v + numCells, out->pointmarkerlist[v])); 4503a074057SBarry Smith } 4513a074057SBarry Smith if (interpolate) { 4523a074057SBarry Smith PetscInt e, f; 4533a074057SBarry Smith 4543a074057SBarry Smith for (e = 0; e < out->numberofedges; e++) { 4553a074057SBarry Smith if (out->edgemarkerlist[e]) { 4563a074057SBarry Smith const PetscInt vertices[2] = {out->edgelist[e * 2 + 0] + numCells, out->edgelist[e * 2 + 1] + numCells}; 4573a074057SBarry Smith const PetscInt *edges; 4583a074057SBarry Smith PetscInt numEdges; 4593a074057SBarry Smith 4609566063dSJacob Faibussowitsch PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 46163a3b9bcSJacob Faibussowitsch PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges); 4629566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, edges[0], out->edgemarkerlist[e])); 4639566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 4643a074057SBarry Smith } 4653a074057SBarry Smith } 4663a074057SBarry Smith for (f = 0; f < out->numberoftrifaces; f++) { 4673a074057SBarry Smith if (out->trifacemarkerlist[f]) { 4683a074057SBarry Smith const PetscInt vertices[3] = {out->trifacelist[f * 3 + 0] + numCells, out->trifacelist[f * 3 + 1] + numCells, out->trifacelist[f * 3 + 2] + numCells}; 4693a074057SBarry Smith const PetscInt *faces; 4703a074057SBarry Smith PetscInt numFaces; 4713a074057SBarry Smith 4729566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces)); 47363a3b9bcSJacob Faibussowitsch PetscCheck(numFaces == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces); 4749566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, faces[0], out->trifacemarkerlist[f])); 4759566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces)); 4763a074057SBarry Smith } 4773a074057SBarry Smith } 4783a074057SBarry Smith } 4790fdc7489SMatthew Knepley 4800fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS 4810fdc7489SMatthew Knepley { 4820fdc7489SMatthew Knepley DMLabel bodyLabel; 4830fdc7489SMatthew Knepley PetscContainer modelObj; 4840fdc7489SMatthew Knepley PetscInt cStart, cEnd, c, eStart, eEnd, fStart, fEnd; 4850fdc7489SMatthew Knepley ego *bodies; 4860fdc7489SMatthew Knepley ego model, geom; 4870fdc7489SMatthew Knepley int Nb, oclass, mtype, *senses; 4880fdc7489SMatthew Knepley 4890fdc7489SMatthew Knepley /* Get Attached EGADS Model from Original DMPlex */ 4909566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj)); 4910fdc7489SMatthew Knepley if (modelObj) { 4929566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(modelObj, (void **)&model)); 4939566063dSJacob Faibussowitsch PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 4940fdc7489SMatthew Knepley /* Transfer EGADS Model to Volumetric Mesh */ 4959566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)*dmRefined, "EGADS Model", (PetscObject)modelObj)); 4960fdc7489SMatthew Knepley 4970fdc7489SMatthew Knepley /* Set Cell Labels */ 4989566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dmRefined, "EGADS Body ID", &bodyLabel)); 4999566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dmRefined, 0, &cStart, &cEnd)); 5009566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dmRefined, 1, &fStart, &fEnd)); 5019566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(*dmRefined, 1, &eStart, &eEnd)); 5020fdc7489SMatthew Knepley 5030fdc7489SMatthew Knepley for (c = cStart; c < cEnd; ++c) { 5040fdc7489SMatthew Knepley PetscReal centroid[3] = {0., 0., 0.}; 5050fdc7489SMatthew Knepley PetscInt b; 5060fdc7489SMatthew Knepley 5070fdc7489SMatthew Knepley /* Deterimine what body the cell's centroid is located in */ 5080fdc7489SMatthew Knepley if (!interpolate) { 5090fdc7489SMatthew Knepley PetscSection coordSection; 5100fdc7489SMatthew Knepley Vec coordinates; 5110fdc7489SMatthew Knepley PetscScalar *coords = NULL; 5120fdc7489SMatthew Knepley PetscInt coordSize, s, d; 5130fdc7489SMatthew Knepley 5149566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(*dmRefined, &coordinates)); 5159566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dmRefined, &coordSection)); 5169566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords)); 5179371c9d4SSatish Balay for (s = 0; s < coordSize; ++s) 5189371c9d4SSatish Balay for (d = 0; d < dim; ++d) centroid[d] += coords[s * dim + d]; 5199566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords)); 5201baa6e33SBarry Smith } else PetscCall(DMPlexComputeCellGeometryFVM(*dmRefined, c, NULL, centroid, NULL)); 5210fdc7489SMatthew Knepley for (b = 0; b < Nb; ++b) { 5220fdc7489SMatthew Knepley if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 5230fdc7489SMatthew Knepley } 5240fdc7489SMatthew Knepley if (b < Nb) { 5250fdc7489SMatthew Knepley PetscInt cval = b, eVal, fVal; 5260fdc7489SMatthew Knepley PetscInt *closure = NULL, Ncl, cl; 5270fdc7489SMatthew Knepley 5289566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, c, cval)); 5299566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure)); 5300fdc7489SMatthew Knepley for (cl = 0; cl < Ncl; cl += 2) { 5310fdc7489SMatthew Knepley const PetscInt p = closure[cl]; 5320fdc7489SMatthew Knepley 5330fdc7489SMatthew Knepley if (p >= eStart && p < eEnd) { 5349566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, p, &eVal)); 5359566063dSJacob Faibussowitsch if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 5360fdc7489SMatthew Knepley } 5370fdc7489SMatthew Knepley if (p >= fStart && p < fEnd) { 5389566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, p, &fVal)); 5399566063dSJacob Faibussowitsch if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 5400fdc7489SMatthew Knepley } 5410fdc7489SMatthew Knepley } 5429566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure)); 5430fdc7489SMatthew Knepley } 5440fdc7489SMatthew Knepley } 5450fdc7489SMatthew Knepley } 5460fdc7489SMatthew Knepley } 5470fdc7489SMatthew Knepley #endif 5489566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE)); 5493a074057SBarry Smith } 5509566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelDestroy(&universal)); 5519566063dSJacob Faibussowitsch PetscCall(PLCDestroy(&in)); 5529566063dSJacob Faibussowitsch PetscCall(PLCDestroy(&out)); 5533a074057SBarry Smith PetscFunctionReturn(0); 5543a074057SBarry Smith } 555