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 */ 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; 260fdc7489SMatthew Knepley DMUniversalLabel universal; 270fdc7489SMatthew Knepley PetscInt vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, verbose = 0; 280fdc7489SMatthew Knepley DMPlexInterpolatedFlag isInterpolated; 293a074057SBarry Smith PetscMPIInt rank; 303a074057SBarry Smith 313a074057SBarry Smith PetscFunctionBegin; 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,((PetscObject) boundary)->prefix, "-ctetgen_verbose", &verbose, NULL)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)boundary,&comm)); 34*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsInterpolatedCollective(boundary, &isInterpolated)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelCreate(boundary, &universal)); 370fdc7489SMatthew Knepley 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(PLCCreate(&in)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(PLCCreate(&out)); 403a074057SBarry Smith 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd)); 423a074057SBarry Smith in->numberofpoints = vEnd - vStart; 433a074057SBarry Smith if (in->numberofpoints > 0) { 443a074057SBarry Smith PetscSection coordSection; 453a074057SBarry Smith Vec coordinates; 460fdc7489SMatthew Knepley const PetscScalar *array; 473a074057SBarry Smith 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberofpoints*dim, &in->pointlist)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberofpoints, &in->pointmarkerlist)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(boundary, &coordinates)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(boundary, &coordSection)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(coordinates, &array)); 533a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 543a074057SBarry Smith const PetscInt idx = v - vStart; 550fdc7489SMatthew Knepley PetscInt off, d, m; 563a074057SBarry Smith 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(coordSection, v, &off)); 580fdc7489SMatthew Knepley for (d = 0; d < dim; ++d) in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(universal->label, v, &m)); 603a074057SBarry Smith in->pointmarkerlist[idx] = (int) m; 613a074057SBarry Smith } 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(coordinates, &array)); 633a074057SBarry Smith } 643a074057SBarry Smith 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(boundary, 1, &eStart, &eEnd)); 660fdc7489SMatthew Knepley in->numberofedges = eEnd - eStart; 670fdc7489SMatthew Knepley if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberofedges > 0) { 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberofedges*2, &in->edgelist)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberofedges, &in->edgemarkerlist)); 700fdc7489SMatthew Knepley for (e = eStart; e < eEnd; ++e) { 710fdc7489SMatthew Knepley const PetscInt idx = e - eStart; 720fdc7489SMatthew Knepley const PetscInt *cone; 730fdc7489SMatthew Knepley PetscInt coneSize, val; 740fdc7489SMatthew Knepley 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetConeSize(boundary, e, &coneSize)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(boundary, e, &cone)); 770fdc7489SMatthew Knepley in->edgelist[idx*2] = cone[0] - vStart; 780fdc7489SMatthew Knepley in->edgelist[idx*2 + 1] = cone[1] - vStart; 790fdc7489SMatthew Knepley 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(universal->label, e, &val)); 810fdc7489SMatthew Knepley in->edgemarkerlist[idx] = (int) val; 820fdc7489SMatthew Knepley } 830fdc7489SMatthew Knepley } 840fdc7489SMatthew Knepley 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd)); 863a074057SBarry Smith in->numberoffacets = fEnd - fStart; 873a074057SBarry Smith if (in->numberoffacets > 0) { 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberoffacets, &in->facetlist)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberoffacets, &in->facetmarkerlist)); 903a074057SBarry Smith for (f = fStart; f < fEnd; ++f) { 913a074057SBarry Smith const PetscInt idx = f - fStart; 923a074057SBarry Smith PetscInt *points = NULL, numPoints, p, numVertices = 0, v, m = -1; 933a074057SBarry Smith polygon *poly; 943a074057SBarry Smith 953a074057SBarry Smith in->facetlist[idx].numberofpolygons = 1; 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->facetlist[idx].numberofpolygons, &in->facetlist[idx].polygonlist)); 973a074057SBarry Smith in->facetlist[idx].numberofholes = 0; 983a074057SBarry Smith in->facetlist[idx].holelist = NULL; 993a074057SBarry Smith 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points)); 1013a074057SBarry Smith for (p = 0; p < numPoints*2; p += 2) { 1023a074057SBarry Smith const PetscInt point = points[p]; 1033a074057SBarry Smith if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point; 1043a074057SBarry Smith } 1053a074057SBarry Smith 1063a074057SBarry Smith poly = in->facetlist[idx].polygonlist; 1073a074057SBarry Smith poly->numberofvertices = numVertices; 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(poly->numberofvertices, &poly->vertexlist)); 1093a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 1103a074057SBarry Smith const PetscInt vIdx = points[v] - vStart; 1113a074057SBarry Smith poly->vertexlist[v] = vIdx; 1123a074057SBarry Smith } 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(universal->label, f, &m)); 1143a074057SBarry Smith in->facetmarkerlist[idx] = (int) m; 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points)); 1163a074057SBarry Smith } 1173a074057SBarry Smith } 118dd400576SPatrick Sanan if (rank == 0) { 1193a074057SBarry Smith TetGenOpts t; 1203a074057SBarry Smith 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(TetGenOptsInitialize(&t)); 1223a074057SBarry Smith t.in = boundary; /* Should go away */ 1233a074057SBarry Smith t.plc = 1; 1243a074057SBarry Smith t.quality = 1; 1253a074057SBarry Smith t.edgesout = 1; 1263a074057SBarry Smith t.zeroindex = 1; 1273a074057SBarry Smith t.quiet = 1; 1283a074057SBarry Smith t.verbose = verbose; 1290fdc7489SMatthew Knepley #if 0 1300fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS 1310fdc7489SMatthew Knepley /* Need to add in more TetGen code */ 1320fdc7489SMatthew Knepley t.nobisect = 1; /* Add Y to preserve Surface Mesh for EGADS */ 1330fdc7489SMatthew Knepley #endif 1340fdc7489SMatthew Knepley #endif 1350fdc7489SMatthew Knepley 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(TetGenCheckOpts(&t)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(TetGenTetrahedralize(&t, in, out)); 1383a074057SBarry Smith } 1393a074057SBarry Smith { 1403a074057SBarry Smith const PetscInt numCorners = 4; 1413a074057SBarry Smith const PetscInt numCells = out->numberoftetrahedra; 1423a074057SBarry Smith const PetscInt numVertices = out->numberofpoints; 1430fdc7489SMatthew Knepley PetscReal *meshCoords = NULL; 1440fdc7489SMatthew Knepley PetscInt *cells = NULL; 1453a074057SBarry Smith 146a4a685f2SJacob Faibussowitsch if (sizeof (PetscReal) == sizeof (out->pointlist[0])) { 147a4a685f2SJacob Faibussowitsch meshCoords = (PetscReal *) out->pointlist; 148a4a685f2SJacob Faibussowitsch } else { 1493a074057SBarry Smith PetscInt i; 1503a074057SBarry Smith 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(dim * numVertices, &meshCoords)); 1520fdc7489SMatthew Knepley for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out->pointlist[i]; 1533a074057SBarry Smith } 154a4a685f2SJacob Faibussowitsch if (sizeof (PetscInt) == sizeof (out->tetrahedronlist[0])) { 155a4a685f2SJacob Faibussowitsch cells = (PetscInt *) out->tetrahedronlist; 156a4a685f2SJacob Faibussowitsch } else { 157a4a685f2SJacob Faibussowitsch PetscInt i; 158a4a685f2SJacob Faibussowitsch 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(numCells * numCorners, &cells)); 1600fdc7489SMatthew Knepley for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt) out->tetrahedronlist[i]; 161a4a685f2SJacob Faibussowitsch } 1623a074057SBarry Smith 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInvertCells_CTetgen(numCells, numCorners, cells)); 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm)); 165a4a685f2SJacob Faibussowitsch if (sizeof (PetscReal) != sizeof (out->pointlist[0])) { 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(meshCoords)); 1673a074057SBarry Smith } 168a4a685f2SJacob Faibussowitsch if (sizeof (PetscInt) != sizeof (out->tetrahedronlist[0])) { 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(cells)); 170a4a685f2SJacob Faibussowitsch } 1710fdc7489SMatthew Knepley 1723a074057SBarry Smith /* Set labels */ 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dm)); 1743a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 1753a074057SBarry Smith if (out->pointmarkerlist[v]) { 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, v+numCells, out->pointmarkerlist[v])); 1773a074057SBarry Smith } 1783a074057SBarry Smith } 1793a074057SBarry Smith if (interpolate) { 1803a074057SBarry Smith PetscInt e; 1813a074057SBarry Smith 1823a074057SBarry Smith for (e = 0; e < out->numberofedges; e++) { 1833a074057SBarry Smith if (out->edgemarkerlist[e]) { 1843a074057SBarry Smith const PetscInt vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells}; 1853a074057SBarry Smith const PetscInt *edges; 1863a074057SBarry Smith PetscInt numEdges; 1873a074057SBarry Smith 188*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges)); 189*5f80ce2aSJacob Faibussowitsch PetscCheck(numEdges == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 190*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, edges[0], out->edgemarkerlist[e])); 191*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges)); 1923a074057SBarry Smith } 1933a074057SBarry Smith } 1943a074057SBarry Smith for (f = 0; f < out->numberoftrifaces; f++) { 1953a074057SBarry Smith if (out->trifacemarkerlist[f]) { 1963a074057SBarry Smith const PetscInt vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells}; 1973a074057SBarry Smith const PetscInt *faces; 1983a074057SBarry Smith PetscInt numFaces; 1993a074057SBarry Smith 200*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces)); 201*5f80ce2aSJacob Faibussowitsch PetscCheck(numFaces == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces); 202*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, faces[0], out->trifacemarkerlist[f])); 203*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces)); 2043a074057SBarry Smith } 2053a074057SBarry Smith } 2063a074057SBarry Smith } 2070fdc7489SMatthew Knepley 2080fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS 2090fdc7489SMatthew Knepley { 2100fdc7489SMatthew Knepley DMLabel bodyLabel; 2110fdc7489SMatthew Knepley PetscContainer modelObj; 2120fdc7489SMatthew Knepley PetscInt cStart, cEnd, c, eStart, eEnd, fStart, fEnd; 2130fdc7489SMatthew Knepley ego *bodies; 2140fdc7489SMatthew Knepley ego model, geom; 2150fdc7489SMatthew Knepley int Nb, oclass, mtype, *senses; 2160fdc7489SMatthew Knepley 2170fdc7489SMatthew Knepley /* Get Attached EGADS Model from Original DMPlex */ 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectQuery((PetscObject) boundary, "EGADS Model", (PetscObject *) &modelObj)); 2190fdc7489SMatthew Knepley if (modelObj) { 220*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerGetPointer(modelObj, (void **) &model)); 221*5f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 2220fdc7489SMatthew Knepley /* Transfer EGADS Model to Volumetric Mesh */ 223*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose((PetscObject) *dm, "EGADS Model", (PetscObject) modelObj)); 2240fdc7489SMatthew Knepley 2250fdc7489SMatthew Knepley /* Set Cell Labels */ 226*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(*dm, "EGADS Body ID", &bodyLabel)); 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 228*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd)); 229*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(*dm, 1, &eStart, &eEnd)); 2300fdc7489SMatthew Knepley 2310fdc7489SMatthew Knepley for (c = cStart; c < cEnd; ++c) { 2320fdc7489SMatthew Knepley PetscReal centroid[3] = {0., 0., 0.}; 2330fdc7489SMatthew Knepley PetscInt b; 2340fdc7489SMatthew Knepley 2350fdc7489SMatthew Knepley /* Deterimine what body the cell's centroid is located in */ 2360fdc7489SMatthew Knepley if (!interpolate) { 2370fdc7489SMatthew Knepley PetscSection coordSection; 2380fdc7489SMatthew Knepley Vec coordinates; 2390fdc7489SMatthew Knepley PetscScalar *coords = NULL; 2400fdc7489SMatthew Knepley PetscInt coordSize, s, d; 2410fdc7489SMatthew Knepley 242*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(*dm, &coordinates)); 243*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(*dm, &coordSection)); 244*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(*dm, coordSection, coordinates, c, &coordSize, &coords)); 2450fdc7489SMatthew Knepley for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d]; 246*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(*dm, coordSection, coordinates, c, &coordSize, &coords)); 2470fdc7489SMatthew Knepley } else { 248*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeCellGeometryFVM(*dm, c, NULL, centroid, NULL)); 2490fdc7489SMatthew Knepley } 2500fdc7489SMatthew Knepley for (b = 0; b < Nb; ++b) { 2510fdc7489SMatthew Knepley if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 2520fdc7489SMatthew Knepley } 2530fdc7489SMatthew Knepley if (b < Nb) { 2540fdc7489SMatthew Knepley PetscInt cval = b, eVal, fVal; 2550fdc7489SMatthew Knepley PetscInt *closure = NULL, Ncl, cl; 2560fdc7489SMatthew Knepley 257*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, c, cval)); 258*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure)); 259c1cad2e7SMatthew G. Knepley for (cl = 0; cl < Ncl; ++cl) { 260c1cad2e7SMatthew G. Knepley const PetscInt p = closure[cl*2]; 2610fdc7489SMatthew Knepley 2620fdc7489SMatthew Knepley if (p >= eStart && p < eEnd) { 263*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(bodyLabel, p, &eVal)); 264*5f80ce2aSJacob Faibussowitsch if (eVal < 0) CHKERRQ(DMLabelSetValue(bodyLabel, p, cval)); 2650fdc7489SMatthew Knepley } 2660fdc7489SMatthew Knepley if (p >= fStart && p < fEnd) { 267*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(bodyLabel, p, &fVal)); 268*5f80ce2aSJacob Faibussowitsch if (fVal < 0) CHKERRQ(DMLabelSetValue(bodyLabel, p, cval)); 2690fdc7489SMatthew Knepley } 2700fdc7489SMatthew Knepley } 271*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure)); 2720fdc7489SMatthew Knepley } 2730fdc7489SMatthew Knepley } 2740fdc7489SMatthew Knepley } 2750fdc7489SMatthew Knepley } 2760fdc7489SMatthew Knepley #endif 277*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetRefinementUniform(*dm, PETSC_FALSE)); 2783a074057SBarry Smith } 2793a074057SBarry Smith 280*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelDestroy(&universal)); 281*5f80ce2aSJacob Faibussowitsch CHKERRQ(PLCDestroy(&in)); 282*5f80ce2aSJacob Faibussowitsch CHKERRQ(PLCDestroy(&out)); 2833a074057SBarry Smith PetscFunctionReturn(0); 2843a074057SBarry Smith } 2853a074057SBarry Smith 2863a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined) 2873a074057SBarry Smith { 2883a074057SBarry Smith MPI_Comm comm; 2893a074057SBarry Smith const PetscInt dim = 3; 2903a074057SBarry Smith PLC *in, *out; 2910fdc7489SMatthew Knepley DMUniversalLabel universal; 2920fdc7489SMatthew Knepley PetscInt vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, cStart, cEnd, c, verbose = 0; 2930fdc7489SMatthew Knepley DMPlexInterpolatedFlag isInterpolated; 2943a074057SBarry Smith PetscMPIInt rank; 2953a074057SBarry Smith 2963a074057SBarry Smith PetscFunctionBegin; 297*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,((PetscObject) dm)->prefix, "-ctetgen_verbose", &verbose, NULL)); 298*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)dm,&comm)); 299*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 300*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsInterpolatedCollective(dm, &isInterpolated)); 301*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelCreate(dm, &universal)); 3020fdc7489SMatthew Knepley 303*5f80ce2aSJacob Faibussowitsch CHKERRQ(PLCCreate(&in)); 304*5f80ce2aSJacob Faibussowitsch CHKERRQ(PLCCreate(&out)); 3053a074057SBarry Smith 306*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 3073a074057SBarry Smith in->numberofpoints = vEnd - vStart; 3083a074057SBarry Smith if (in->numberofpoints > 0) { 3093a074057SBarry Smith PetscSection coordSection; 3103a074057SBarry Smith Vec coordinates; 3113a074057SBarry Smith PetscScalar *array; 3123a074057SBarry Smith 313*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberofpoints*dim, &in->pointlist)); 314*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberofpoints, &in->pointmarkerlist)); 315*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates)); 316*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(dm, &coordSection)); 317*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(coordinates, &array)); 3183a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 3193a074057SBarry Smith const PetscInt idx = v - vStart; 3200fdc7489SMatthew Knepley PetscInt off, d, m; 3213a074057SBarry Smith 322*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(coordSection, v, &off)); 3230fdc7489SMatthew Knepley for (d = 0; d < dim; ++d) in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 324*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(universal->label, v, &m)); 3253a074057SBarry Smith in->pointmarkerlist[idx] = (int) m; 3263a074057SBarry Smith } 327*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(coordinates, &array)); 3283a074057SBarry Smith } 3293a074057SBarry Smith 330*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd)); 3310fdc7489SMatthew Knepley in->numberofedges = eEnd - eStart; 3320fdc7489SMatthew Knepley if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberofedges > 0) { 333*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberofedges * 2, &in->edgelist)); 334*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberofedges, &in->edgemarkerlist)); 3350fdc7489SMatthew Knepley for (e = eStart; e < eEnd; ++e) { 3360fdc7489SMatthew Knepley const PetscInt idx = e - eStart; 3370fdc7489SMatthew Knepley const PetscInt *cone; 3380fdc7489SMatthew Knepley PetscInt coneSize, val; 3390fdc7489SMatthew Knepley 340*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetConeSize(dm, e, &coneSize)); 341*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(dm, e, &cone)); 3420fdc7489SMatthew Knepley in->edgelist[idx*2] = cone[0] - vStart; 3430fdc7489SMatthew Knepley in->edgelist[idx*2 + 1] = cone[1] - vStart; 3440fdc7489SMatthew Knepley 345*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(universal->label, e, &val)); 3460fdc7489SMatthew Knepley in->edgemarkerlist[idx] = (int) val; 3470fdc7489SMatthew Knepley } 3480fdc7489SMatthew Knepley } 3490fdc7489SMatthew Knepley 350*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 3510fdc7489SMatthew Knepley in->numberoftrifaces = 0; 3520fdc7489SMatthew Knepley for (f = fStart; f < fEnd; ++f) { 3530fdc7489SMatthew Knepley PetscInt supportSize; 3540fdc7489SMatthew Knepley 355*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetSupportSize(dm, f, &supportSize)); 3560fdc7489SMatthew Knepley if (supportSize == 1) ++in->numberoftrifaces; 3570fdc7489SMatthew Knepley } 3580fdc7489SMatthew Knepley if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberoftrifaces > 0) { 3590fdc7489SMatthew Knepley PetscInt tf = 0; 3600fdc7489SMatthew Knepley 361*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberoftrifaces*3, &in->trifacelist)); 362*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberoftrifaces, &in->trifacemarkerlist)); 3630fdc7489SMatthew Knepley for (f = fStart; f < fEnd; ++f) { 3640fdc7489SMatthew Knepley PetscInt *points = NULL; 3650fdc7489SMatthew Knepley PetscInt supportSize, numPoints, p, Nv = 0, val; 3660fdc7489SMatthew Knepley 367*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetSupportSize(dm, f, &supportSize)); 3680fdc7489SMatthew Knepley if (supportSize != 1) continue; 369*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points)); 3700fdc7489SMatthew Knepley for (p = 0; p < numPoints*2; p += 2) { 3710fdc7489SMatthew Knepley const PetscInt point = points[p]; 3720fdc7489SMatthew Knepley if ((point >= vStart) && (point < vEnd)) in->trifacelist[tf*3 + Nv++] = point - vStart; 3730fdc7489SMatthew Knepley } 374*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points)); 375*5f80ce2aSJacob Faibussowitsch PetscCheck(Nv == 3,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %D has %D vertices, not 3", f, Nv); 376*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(universal->label, f, &val)); 3770fdc7489SMatthew Knepley in->trifacemarkerlist[tf] = (int) val; 3780fdc7489SMatthew Knepley ++tf; 3790fdc7489SMatthew Knepley } 3800fdc7489SMatthew Knepley } 3810fdc7489SMatthew Knepley 382*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 3833a074057SBarry Smith in->numberofcorners = 4; 3843a074057SBarry Smith in->numberoftetrahedra = cEnd - cStart; 3853a074057SBarry Smith in->tetrahedronvolumelist = maxVolumes; 3863a074057SBarry Smith if (in->numberoftetrahedra > 0) { 387*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(in->numberoftetrahedra*in->numberofcorners, &in->tetrahedronlist)); 3883a074057SBarry Smith for (c = cStart; c < cEnd; ++c) { 3893a074057SBarry Smith const PetscInt idx = c - cStart; 3903a074057SBarry Smith PetscInt *closure = NULL; 3913a074057SBarry Smith PetscInt closureSize; 3923a074057SBarry Smith 393*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 394*5f80ce2aSJacob 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); 3950fdc7489SMatthew Knepley for (v = 0; v < 4; ++v) in->tetrahedronlist[idx*in->numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart; 396*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 3973a074057SBarry Smith } 3983a074057SBarry Smith } 3990fdc7489SMatthew Knepley 400dd400576SPatrick Sanan if (rank == 0) { 4013a074057SBarry Smith TetGenOpts t; 4023a074057SBarry Smith 403*5f80ce2aSJacob Faibussowitsch CHKERRQ(TetGenOptsInitialize(&t)); 4043a074057SBarry Smith 4053a074057SBarry Smith t.in = dm; /* Should go away */ 4063a074057SBarry Smith t.refine = 1; 4073a074057SBarry Smith t.varvolume = 1; 4083a074057SBarry Smith t.quality = 1; 4093a074057SBarry Smith t.edgesout = 1; 4103a074057SBarry Smith t.zeroindex = 1; 4113a074057SBarry Smith t.quiet = 1; 4123a074057SBarry Smith t.verbose = verbose; /* Change this */ 4133a074057SBarry Smith 414*5f80ce2aSJacob Faibussowitsch CHKERRQ(TetGenCheckOpts(&t)); 415*5f80ce2aSJacob Faibussowitsch CHKERRQ(TetGenTetrahedralize(&t, in, out)); 4163a074057SBarry Smith } 4170fdc7489SMatthew Knepley 4183a074057SBarry Smith in->tetrahedronvolumelist = NULL; 4193a074057SBarry Smith { 4203a074057SBarry Smith const PetscInt numCorners = 4; 4213a074057SBarry Smith const PetscInt numCells = out->numberoftetrahedra; 4223a074057SBarry Smith const PetscInt numVertices = out->numberofpoints; 4230fdc7489SMatthew Knepley PetscReal *meshCoords = NULL; 4240fdc7489SMatthew Knepley PetscInt *cells = NULL; 4250fdc7489SMatthew Knepley PetscBool interpolate = isInterpolated == DMPLEX_INTERPOLATED_FULL ? PETSC_TRUE : PETSC_FALSE; 4263a074057SBarry Smith 427a4a685f2SJacob Faibussowitsch if (sizeof (PetscReal) == sizeof (out->pointlist[0])) { 428a4a685f2SJacob Faibussowitsch meshCoords = (PetscReal *) out->pointlist; 429a4a685f2SJacob Faibussowitsch } else { 4303a074057SBarry Smith PetscInt i; 4313a074057SBarry Smith 432*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(dim * numVertices, &meshCoords)); 4330fdc7489SMatthew Knepley for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out->pointlist[i]; 4343a074057SBarry Smith } 435a4a685f2SJacob Faibussowitsch if (sizeof (PetscInt) == sizeof (out->tetrahedronlist[0])) { 436a4a685f2SJacob Faibussowitsch cells = (PetscInt *) out->tetrahedronlist; 437a4a685f2SJacob Faibussowitsch } else { 438a4a685f2SJacob Faibussowitsch PetscInt i; 439a4a685f2SJacob Faibussowitsch 440*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(numCells * numCorners, &cells)); 4410fdc7489SMatthew Knepley for (i = 0; i < numCells * numCorners; ++i) cells[i] = (PetscInt) out->tetrahedronlist[i]; 442a4a685f2SJacob Faibussowitsch } 4433a074057SBarry Smith 444*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInvertCells_CTetgen(numCells, numCorners, cells)); 445*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined)); 446*5f80ce2aSJacob Faibussowitsch if (sizeof (PetscReal) != sizeof (out->pointlist[0])) CHKERRQ(PetscFree(meshCoords)); 447*5f80ce2aSJacob Faibussowitsch if (sizeof (PetscInt) != sizeof (out->tetrahedronlist[0])) CHKERRQ(PetscFree(cells)); 4480fdc7489SMatthew Knepley 4493a074057SBarry Smith /* Set labels */ 450*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dmRefined)); 4513a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 4523a074057SBarry Smith if (out->pointmarkerlist[v]) { 453*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, v+numCells, out->pointmarkerlist[v])); 4543a074057SBarry Smith } 4553a074057SBarry Smith } 4563a074057SBarry Smith if (interpolate) { 4573a074057SBarry Smith PetscInt e, f; 4583a074057SBarry Smith 4593a074057SBarry Smith for (e = 0; e < out->numberofedges; e++) { 4603a074057SBarry Smith if (out->edgemarkerlist[e]) { 4613a074057SBarry Smith const PetscInt vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells}; 4623a074057SBarry Smith const PetscInt *edges; 4633a074057SBarry Smith PetscInt numEdges; 4643a074057SBarry Smith 465*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 466*5f80ce2aSJacob Faibussowitsch PetscCheck(numEdges == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 467*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, edges[0], out->edgemarkerlist[e])); 468*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 4693a074057SBarry Smith } 4703a074057SBarry Smith } 4713a074057SBarry Smith for (f = 0; f < out->numberoftrifaces; f++) { 4723a074057SBarry Smith if (out->trifacemarkerlist[f]) { 4733a074057SBarry Smith const PetscInt vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells}; 4743a074057SBarry Smith const PetscInt *faces; 4753a074057SBarry Smith PetscInt numFaces; 4763a074057SBarry Smith 477*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces)); 478*5f80ce2aSJacob Faibussowitsch PetscCheck(numFaces == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces); 479*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, faces[0], out->trifacemarkerlist[f])); 480*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces)); 4813a074057SBarry Smith } 4823a074057SBarry Smith } 4833a074057SBarry Smith } 4840fdc7489SMatthew Knepley 4850fdc7489SMatthew Knepley #ifdef PETSC_HAVE_EGADS 4860fdc7489SMatthew Knepley { 4870fdc7489SMatthew Knepley DMLabel bodyLabel; 4880fdc7489SMatthew Knepley PetscContainer modelObj; 4890fdc7489SMatthew Knepley PetscInt cStart, cEnd, c, eStart, eEnd, fStart, fEnd; 4900fdc7489SMatthew Knepley ego *bodies; 4910fdc7489SMatthew Knepley ego model, geom; 4920fdc7489SMatthew Knepley int Nb, oclass, mtype, *senses; 4930fdc7489SMatthew Knepley 4940fdc7489SMatthew Knepley /* Get Attached EGADS Model from Original DMPlex */ 495*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj)); 4960fdc7489SMatthew Knepley if (modelObj) { 497*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerGetPointer(modelObj, (void **) &model)); 498*5f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 4990fdc7489SMatthew Knepley /* Transfer EGADS Model to Volumetric Mesh */ 500*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose((PetscObject) *dmRefined, "EGADS Model", (PetscObject) modelObj)); 5010fdc7489SMatthew Knepley 5020fdc7489SMatthew Knepley /* Set Cell Labels */ 503*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(*dmRefined, "EGADS Body ID", &bodyLabel)); 504*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(*dmRefined, 0, &cStart, &cEnd)); 505*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(*dmRefined, 1, &fStart, &fEnd)); 506*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(*dmRefined, 1, &eStart, &eEnd)); 5070fdc7489SMatthew Knepley 5080fdc7489SMatthew Knepley for (c = cStart; c < cEnd; ++c) { 5090fdc7489SMatthew Knepley PetscReal centroid[3] = {0., 0., 0.}; 5100fdc7489SMatthew Knepley PetscInt b; 5110fdc7489SMatthew Knepley 5120fdc7489SMatthew Knepley /* Deterimine what body the cell's centroid is located in */ 5130fdc7489SMatthew Knepley if (!interpolate) { 5140fdc7489SMatthew Knepley PetscSection coordSection; 5150fdc7489SMatthew Knepley Vec coordinates; 5160fdc7489SMatthew Knepley PetscScalar *coords = NULL; 5170fdc7489SMatthew Knepley PetscInt coordSize, s, d; 5180fdc7489SMatthew Knepley 519*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(*dmRefined, &coordinates)); 520*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(*dmRefined, &coordSection)); 521*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords)); 5220fdc7489SMatthew Knepley for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d]; 523*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords)); 5240fdc7489SMatthew Knepley } else { 525*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeCellGeometryFVM(*dmRefined, c, NULL, centroid, NULL)); 5260fdc7489SMatthew Knepley } 5270fdc7489SMatthew Knepley for (b = 0; b < Nb; ++b) { 5280fdc7489SMatthew Knepley if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 5290fdc7489SMatthew Knepley } 5300fdc7489SMatthew Knepley if (b < Nb) { 5310fdc7489SMatthew Knepley PetscInt cval = b, eVal, fVal; 5320fdc7489SMatthew Knepley PetscInt *closure = NULL, Ncl, cl; 5330fdc7489SMatthew Knepley 534*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, c, cval)); 535*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure)); 5360fdc7489SMatthew Knepley for (cl = 0; cl < Ncl; cl += 2) { 5370fdc7489SMatthew Knepley const PetscInt p = closure[cl]; 5380fdc7489SMatthew Knepley 5390fdc7489SMatthew Knepley if (p >= eStart && p < eEnd) { 540*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(bodyLabel, p, &eVal)); 541*5f80ce2aSJacob Faibussowitsch if (eVal < 0) CHKERRQ(DMLabelSetValue(bodyLabel, p, cval)); 5420fdc7489SMatthew Knepley } 5430fdc7489SMatthew Knepley if (p >= fStart && p < fEnd) { 544*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(bodyLabel, p, &fVal)); 545*5f80ce2aSJacob Faibussowitsch if (fVal < 0) CHKERRQ(DMLabelSetValue(bodyLabel, p, cval)); 5460fdc7489SMatthew Knepley } 5470fdc7489SMatthew Knepley } 548*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure)); 5490fdc7489SMatthew Knepley } 5500fdc7489SMatthew Knepley } 5510fdc7489SMatthew Knepley } 5520fdc7489SMatthew Knepley } 5530fdc7489SMatthew Knepley #endif 554*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE)); 5553a074057SBarry Smith } 556*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMUniversalLabelDestroy(&universal)); 557*5f80ce2aSJacob Faibussowitsch CHKERRQ(PLCDestroy(&in)); 558*5f80ce2aSJacob Faibussowitsch CHKERRQ(PLCDestroy(&out)); 5593a074057SBarry Smith PetscFunctionReturn(0); 5603a074057SBarry Smith } 561