1*3a074057SBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2*3a074057SBarry Smith 3*3a074057SBarry Smith #include <ctetgen.h> 4*3a074057SBarry Smith 5*3a074057SBarry Smith /* This is to fix the tetrahedron orientation from TetGen */ 6*3a074057SBarry Smith static PetscErrorCode DMPlexInvertCells_Internal(PetscInt dim, PetscInt numCells, PetscInt numCorners, int cells[]) 7*3a074057SBarry Smith { 8*3a074057SBarry Smith PetscInt bound = numCells*numCorners, coff; 9*3a074057SBarry Smith PetscErrorCode ierr; 10*3a074057SBarry Smith 11*3a074057SBarry Smith PetscFunctionBegin; 12*3a074057SBarry Smith for (coff = 0; coff < bound; coff += numCorners) { 13*3a074057SBarry Smith ierr = DMPlexInvertCell(dim, numCorners, &cells[coff]);CHKERRQ(ierr); 14*3a074057SBarry Smith } 15*3a074057SBarry Smith PetscFunctionReturn(0); 16*3a074057SBarry Smith } 17*3a074057SBarry Smith 18*3a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerate_CTetgen(DM boundary, PetscBool interpolate, DM *dm) 19*3a074057SBarry Smith { 20*3a074057SBarry Smith MPI_Comm comm; 21*3a074057SBarry Smith const PetscInt dim = 3; 22*3a074057SBarry Smith const char *labelName = "marker"; 23*3a074057SBarry Smith PLC *in, *out; 24*3a074057SBarry Smith DMLabel label; 25*3a074057SBarry Smith PetscInt verbose = 0, vStart, vEnd, v, fStart, fEnd, f; 26*3a074057SBarry Smith PetscMPIInt rank; 27*3a074057SBarry Smith PetscErrorCode ierr; 28*3a074057SBarry Smith 29*3a074057SBarry Smith PetscFunctionBegin; 30*3a074057SBarry Smith ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr); 31*3a074057SBarry Smith ierr = PetscOptionsGetInt(NULL,((PetscObject) boundary)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr); 32*3a074057SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 33*3a074057SBarry Smith ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr); 34*3a074057SBarry Smith ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr); 35*3a074057SBarry Smith ierr = PLCCreate(&in);CHKERRQ(ierr); 36*3a074057SBarry Smith ierr = PLCCreate(&out);CHKERRQ(ierr); 37*3a074057SBarry Smith 38*3a074057SBarry Smith in->numberofpoints = vEnd - vStart; 39*3a074057SBarry Smith if (in->numberofpoints > 0) { 40*3a074057SBarry Smith PetscSection coordSection; 41*3a074057SBarry Smith Vec coordinates; 42*3a074057SBarry Smith PetscScalar *array; 43*3a074057SBarry Smith 44*3a074057SBarry Smith ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr); 45*3a074057SBarry Smith ierr = PetscMalloc1(in->numberofpoints, &in->pointmarkerlist);CHKERRQ(ierr); 46*3a074057SBarry Smith ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr); 47*3a074057SBarry Smith ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr); 48*3a074057SBarry Smith ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 49*3a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 50*3a074057SBarry Smith const PetscInt idx = v - vStart; 51*3a074057SBarry Smith PetscInt off, d, m = -1; 52*3a074057SBarry Smith 53*3a074057SBarry Smith ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 54*3a074057SBarry Smith for (d = 0; d < dim; ++d) { 55*3a074057SBarry Smith in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 56*3a074057SBarry Smith } 57*3a074057SBarry Smith if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);} 58*3a074057SBarry Smith 59*3a074057SBarry Smith in->pointmarkerlist[idx] = (int) m; 60*3a074057SBarry Smith } 61*3a074057SBarry Smith ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 62*3a074057SBarry Smith } 63*3a074057SBarry Smith ierr = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr); 64*3a074057SBarry Smith 65*3a074057SBarry Smith in->numberoffacets = fEnd - fStart; 66*3a074057SBarry Smith if (in->numberoffacets > 0) { 67*3a074057SBarry Smith ierr = PetscMalloc1(in->numberoffacets, &in->facetlist);CHKERRQ(ierr); 68*3a074057SBarry Smith ierr = PetscMalloc1(in->numberoffacets, &in->facetmarkerlist);CHKERRQ(ierr); 69*3a074057SBarry Smith for (f = fStart; f < fEnd; ++f) { 70*3a074057SBarry Smith const PetscInt idx = f - fStart; 71*3a074057SBarry Smith PetscInt *points = NULL, numPoints, p, numVertices = 0, v, m = -1; 72*3a074057SBarry Smith polygon *poly; 73*3a074057SBarry Smith 74*3a074057SBarry Smith in->facetlist[idx].numberofpolygons = 1; 75*3a074057SBarry Smith 76*3a074057SBarry Smith ierr = PetscMalloc1(in->facetlist[idx].numberofpolygons, &in->facetlist[idx].polygonlist);CHKERRQ(ierr); 77*3a074057SBarry Smith 78*3a074057SBarry Smith in->facetlist[idx].numberofholes = 0; 79*3a074057SBarry Smith in->facetlist[idx].holelist = NULL; 80*3a074057SBarry Smith 81*3a074057SBarry Smith ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr); 82*3a074057SBarry Smith for (p = 0; p < numPoints*2; p += 2) { 83*3a074057SBarry Smith const PetscInt point = points[p]; 84*3a074057SBarry Smith if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point; 85*3a074057SBarry Smith } 86*3a074057SBarry Smith 87*3a074057SBarry Smith poly = in->facetlist[idx].polygonlist; 88*3a074057SBarry Smith poly->numberofvertices = numVertices; 89*3a074057SBarry Smith ierr = PetscMalloc1(poly->numberofvertices, &poly->vertexlist);CHKERRQ(ierr); 90*3a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 91*3a074057SBarry Smith const PetscInt vIdx = points[v] - vStart; 92*3a074057SBarry Smith poly->vertexlist[v] = vIdx; 93*3a074057SBarry Smith } 94*3a074057SBarry Smith if (label) {ierr = DMLabelGetValue(label, f, &m);CHKERRQ(ierr);} 95*3a074057SBarry Smith in->facetmarkerlist[idx] = (int) m; 96*3a074057SBarry Smith ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr); 97*3a074057SBarry Smith } 98*3a074057SBarry Smith } 99*3a074057SBarry Smith if (!rank) { 100*3a074057SBarry Smith TetGenOpts t; 101*3a074057SBarry Smith 102*3a074057SBarry Smith ierr = TetGenOptsInitialize(&t);CHKERRQ(ierr); 103*3a074057SBarry Smith t.in = boundary; /* Should go away */ 104*3a074057SBarry Smith t.plc = 1; 105*3a074057SBarry Smith t.quality = 1; 106*3a074057SBarry Smith t.edgesout = 1; 107*3a074057SBarry Smith t.zeroindex = 1; 108*3a074057SBarry Smith t.quiet = 1; 109*3a074057SBarry Smith t.verbose = verbose; 110*3a074057SBarry Smith ierr = TetGenCheckOpts(&t);CHKERRQ(ierr); 111*3a074057SBarry Smith ierr = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr); 112*3a074057SBarry Smith } 113*3a074057SBarry Smith { 114*3a074057SBarry Smith DMLabel glabel = NULL; 115*3a074057SBarry Smith const PetscInt numCorners = 4; 116*3a074057SBarry Smith const PetscInt numCells = out->numberoftetrahedra; 117*3a074057SBarry Smith const PetscInt numVertices = out->numberofpoints; 118*3a074057SBarry Smith double *meshCoords; 119*3a074057SBarry Smith int *cells = out->tetrahedronlist; 120*3a074057SBarry Smith 121*3a074057SBarry Smith if (sizeof (PetscReal) == sizeof (double)) { 122*3a074057SBarry Smith meshCoords = (double *) out->pointlist; 123*3a074057SBarry Smith } 124*3a074057SBarry Smith else { 125*3a074057SBarry Smith PetscInt i; 126*3a074057SBarry Smith 127*3a074057SBarry Smith ierr = PetscMalloc1(3 * numVertices,&meshCoords);CHKERRQ(ierr); 128*3a074057SBarry Smith for (i = 0; i < 3 * numVertices; i++) { 129*3a074057SBarry Smith meshCoords[i] = (PetscReal) out->pointlist[i]; 130*3a074057SBarry Smith } 131*3a074057SBarry Smith } 132*3a074057SBarry Smith 133*3a074057SBarry Smith ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr); 134*3a074057SBarry Smith ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr); 135*3a074057SBarry Smith if (sizeof (PetscReal) != sizeof (double)) { 136*3a074057SBarry Smith ierr = PetscFree(meshCoords);CHKERRQ(ierr); 137*3a074057SBarry Smith } 138*3a074057SBarry Smith if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);} 139*3a074057SBarry Smith /* Set labels */ 140*3a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 141*3a074057SBarry Smith if (out->pointmarkerlist[v]) { 142*3a074057SBarry Smith if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);} 143*3a074057SBarry Smith } 144*3a074057SBarry Smith } 145*3a074057SBarry Smith if (interpolate) { 146*3a074057SBarry Smith PetscInt e; 147*3a074057SBarry Smith 148*3a074057SBarry Smith for (e = 0; e < out->numberofedges; e++) { 149*3a074057SBarry Smith if (out->edgemarkerlist[e]) { 150*3a074057SBarry Smith const PetscInt vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells}; 151*3a074057SBarry Smith const PetscInt *edges; 152*3a074057SBarry Smith PetscInt numEdges; 153*3a074057SBarry Smith 154*3a074057SBarry Smith ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 155*3a074057SBarry Smith if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 156*3a074057SBarry Smith if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);} 157*3a074057SBarry Smith ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 158*3a074057SBarry Smith } 159*3a074057SBarry Smith } 160*3a074057SBarry Smith for (f = 0; f < out->numberoftrifaces; f++) { 161*3a074057SBarry Smith if (out->trifacemarkerlist[f]) { 162*3a074057SBarry Smith const PetscInt vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells}; 163*3a074057SBarry Smith const PetscInt *faces; 164*3a074057SBarry Smith PetscInt numFaces; 165*3a074057SBarry Smith 166*3a074057SBarry Smith ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 167*3a074057SBarry Smith if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces); 168*3a074057SBarry Smith if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);} 169*3a074057SBarry Smith ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 170*3a074057SBarry Smith } 171*3a074057SBarry Smith } 172*3a074057SBarry Smith } 173*3a074057SBarry Smith ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr); 174*3a074057SBarry Smith } 175*3a074057SBarry Smith 176*3a074057SBarry Smith ierr = PLCDestroy(&in);CHKERRQ(ierr); 177*3a074057SBarry Smith ierr = PLCDestroy(&out);CHKERRQ(ierr); 178*3a074057SBarry Smith PetscFunctionReturn(0); 179*3a074057SBarry Smith } 180*3a074057SBarry Smith 181*3a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined) 182*3a074057SBarry Smith { 183*3a074057SBarry Smith MPI_Comm comm; 184*3a074057SBarry Smith const PetscInt dim = 3; 185*3a074057SBarry Smith const char *labelName = "marker"; 186*3a074057SBarry Smith PLC *in, *out; 187*3a074057SBarry Smith DMLabel label; 188*3a074057SBarry Smith PetscInt verbose = 0, vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal; 189*3a074057SBarry Smith PetscMPIInt rank; 190*3a074057SBarry Smith PetscErrorCode ierr; 191*3a074057SBarry Smith 192*3a074057SBarry Smith PetscFunctionBegin; 193*3a074057SBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 194*3a074057SBarry Smith ierr = PetscOptionsGetInt(NULL,((PetscObject) dm)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr); 195*3a074057SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 196*3a074057SBarry Smith ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 197*3a074057SBarry Smith ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 198*3a074057SBarry Smith ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 199*3a074057SBarry Smith ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr); 200*3a074057SBarry Smith ierr = PLCCreate(&in);CHKERRQ(ierr); 201*3a074057SBarry Smith ierr = PLCCreate(&out);CHKERRQ(ierr); 202*3a074057SBarry Smith 203*3a074057SBarry Smith in->numberofpoints = vEnd - vStart; 204*3a074057SBarry Smith if (in->numberofpoints > 0) { 205*3a074057SBarry Smith PetscSection coordSection; 206*3a074057SBarry Smith Vec coordinates; 207*3a074057SBarry Smith PetscScalar *array; 208*3a074057SBarry Smith 209*3a074057SBarry Smith ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr); 210*3a074057SBarry Smith ierr = PetscMalloc1(in->numberofpoints, &in->pointmarkerlist);CHKERRQ(ierr); 211*3a074057SBarry Smith ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 212*3a074057SBarry Smith ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 213*3a074057SBarry Smith ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 214*3a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 215*3a074057SBarry Smith const PetscInt idx = v - vStart; 216*3a074057SBarry Smith PetscInt off, d, m = -1; 217*3a074057SBarry Smith 218*3a074057SBarry Smith ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 219*3a074057SBarry Smith for (d = 0; d < dim; ++d) { 220*3a074057SBarry Smith in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 221*3a074057SBarry Smith } 222*3a074057SBarry Smith if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);} 223*3a074057SBarry Smith 224*3a074057SBarry Smith in->pointmarkerlist[idx] = (int) m; 225*3a074057SBarry Smith } 226*3a074057SBarry Smith ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 227*3a074057SBarry Smith } 228*3a074057SBarry Smith ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 229*3a074057SBarry Smith 230*3a074057SBarry Smith in->numberofcorners = 4; 231*3a074057SBarry Smith in->numberoftetrahedra = cEnd - cStart; 232*3a074057SBarry Smith in->tetrahedronvolumelist = maxVolumes; 233*3a074057SBarry Smith if (in->numberoftetrahedra > 0) { 234*3a074057SBarry Smith ierr = PetscMalloc1(in->numberoftetrahedra*in->numberofcorners, &in->tetrahedronlist);CHKERRQ(ierr); 235*3a074057SBarry Smith for (c = cStart; c < cEnd; ++c) { 236*3a074057SBarry Smith const PetscInt idx = c - cStart; 237*3a074057SBarry Smith PetscInt *closure = NULL; 238*3a074057SBarry Smith PetscInt closureSize; 239*3a074057SBarry Smith 240*3a074057SBarry Smith ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 241*3a074057SBarry Smith if ((closureSize != 5) && (closureSize != 15)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %D vertices in closure", closureSize); 242*3a074057SBarry Smith for (v = 0; v < 4; ++v) { 243*3a074057SBarry Smith in->tetrahedronlist[idx*in->numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart; 244*3a074057SBarry Smith } 245*3a074057SBarry Smith ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 246*3a074057SBarry Smith } 247*3a074057SBarry Smith } 248*3a074057SBarry Smith if (!rank) { 249*3a074057SBarry Smith TetGenOpts t; 250*3a074057SBarry Smith 251*3a074057SBarry Smith ierr = TetGenOptsInitialize(&t);CHKERRQ(ierr); 252*3a074057SBarry Smith 253*3a074057SBarry Smith t.in = dm; /* Should go away */ 254*3a074057SBarry Smith t.refine = 1; 255*3a074057SBarry Smith t.varvolume = 1; 256*3a074057SBarry Smith t.quality = 1; 257*3a074057SBarry Smith t.edgesout = 1; 258*3a074057SBarry Smith t.zeroindex = 1; 259*3a074057SBarry Smith t.quiet = 1; 260*3a074057SBarry Smith t.verbose = verbose; /* Change this */ 261*3a074057SBarry Smith 262*3a074057SBarry Smith ierr = TetGenCheckOpts(&t);CHKERRQ(ierr); 263*3a074057SBarry Smith ierr = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr); 264*3a074057SBarry Smith } 265*3a074057SBarry Smith in->tetrahedronvolumelist = NULL; 266*3a074057SBarry Smith { 267*3a074057SBarry Smith DMLabel rlabel = NULL; 268*3a074057SBarry Smith const PetscInt numCorners = 4; 269*3a074057SBarry Smith const PetscInt numCells = out->numberoftetrahedra; 270*3a074057SBarry Smith const PetscInt numVertices = out->numberofpoints; 271*3a074057SBarry Smith double *meshCoords; 272*3a074057SBarry Smith int *cells = out->tetrahedronlist; 273*3a074057SBarry Smith PetscBool interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE; 274*3a074057SBarry Smith 275*3a074057SBarry Smith if (sizeof (PetscReal) == sizeof (double)) { 276*3a074057SBarry Smith meshCoords = (double *) out->pointlist; 277*3a074057SBarry Smith } 278*3a074057SBarry Smith else { 279*3a074057SBarry Smith PetscInt i; 280*3a074057SBarry Smith 281*3a074057SBarry Smith ierr = PetscMalloc1(3 * numVertices,&meshCoords);CHKERRQ(ierr); 282*3a074057SBarry Smith for (i = 0; i < 3 * numVertices; i++) { 283*3a074057SBarry Smith meshCoords[i] = (PetscReal) out->pointlist[i]; 284*3a074057SBarry Smith } 285*3a074057SBarry Smith } 286*3a074057SBarry Smith 287*3a074057SBarry Smith ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr); 288*3a074057SBarry Smith ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr); 289*3a074057SBarry Smith if (sizeof (PetscReal) != sizeof (double)) { 290*3a074057SBarry Smith ierr = PetscFree(meshCoords);CHKERRQ(ierr); 291*3a074057SBarry Smith } 292*3a074057SBarry Smith if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);} 293*3a074057SBarry Smith /* Set labels */ 294*3a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 295*3a074057SBarry Smith if (out->pointmarkerlist[v]) { 296*3a074057SBarry Smith if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);} 297*3a074057SBarry Smith } 298*3a074057SBarry Smith } 299*3a074057SBarry Smith if (interpolate) { 300*3a074057SBarry Smith PetscInt e, f; 301*3a074057SBarry Smith 302*3a074057SBarry Smith for (e = 0; e < out->numberofedges; e++) { 303*3a074057SBarry Smith if (out->edgemarkerlist[e]) { 304*3a074057SBarry Smith const PetscInt vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells}; 305*3a074057SBarry Smith const PetscInt *edges; 306*3a074057SBarry Smith PetscInt numEdges; 307*3a074057SBarry Smith 308*3a074057SBarry Smith ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 309*3a074057SBarry Smith if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 310*3a074057SBarry Smith if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);} 311*3a074057SBarry Smith ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 312*3a074057SBarry Smith } 313*3a074057SBarry Smith } 314*3a074057SBarry Smith for (f = 0; f < out->numberoftrifaces; f++) { 315*3a074057SBarry Smith if (out->trifacemarkerlist[f]) { 316*3a074057SBarry Smith const PetscInt vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells}; 317*3a074057SBarry Smith const PetscInt *faces; 318*3a074057SBarry Smith PetscInt numFaces; 319*3a074057SBarry Smith 320*3a074057SBarry Smith ierr = DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 321*3a074057SBarry Smith if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces); 322*3a074057SBarry Smith if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);} 323*3a074057SBarry Smith ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 324*3a074057SBarry Smith } 325*3a074057SBarry Smith } 326*3a074057SBarry Smith } 327*3a074057SBarry Smith ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr); 328*3a074057SBarry Smith } 329*3a074057SBarry Smith ierr = PLCDestroy(&in);CHKERRQ(ierr); 330*3a074057SBarry Smith ierr = PLCDestroy(&out);CHKERRQ(ierr); 331*3a074057SBarry Smith PetscFunctionReturn(0); 332*3a074057SBarry Smith } 333