1*3a074057SBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2*3a074057SBarry Smith 3*3a074057SBarry Smith #include <triangle.h> 4*3a074057SBarry Smith 5*3a074057SBarry Smith static PetscErrorCode InitInput_Triangle(struct triangulateio *inputCtx) 6*3a074057SBarry Smith { 7*3a074057SBarry Smith PetscFunctionBegin; 8*3a074057SBarry Smith inputCtx->numberofpoints = 0; 9*3a074057SBarry Smith inputCtx->numberofpointattributes = 0; 10*3a074057SBarry Smith inputCtx->pointlist = NULL; 11*3a074057SBarry Smith inputCtx->pointattributelist = NULL; 12*3a074057SBarry Smith inputCtx->pointmarkerlist = NULL; 13*3a074057SBarry Smith inputCtx->numberofsegments = 0; 14*3a074057SBarry Smith inputCtx->segmentlist = NULL; 15*3a074057SBarry Smith inputCtx->segmentmarkerlist = NULL; 16*3a074057SBarry Smith inputCtx->numberoftriangleattributes = 0; 17*3a074057SBarry Smith inputCtx->trianglelist = NULL; 18*3a074057SBarry Smith inputCtx->numberofholes = 0; 19*3a074057SBarry Smith inputCtx->holelist = NULL; 20*3a074057SBarry Smith inputCtx->numberofregions = 0; 21*3a074057SBarry Smith inputCtx->regionlist = NULL; 22*3a074057SBarry Smith PetscFunctionReturn(0); 23*3a074057SBarry Smith } 24*3a074057SBarry Smith 25*3a074057SBarry Smith static PetscErrorCode InitOutput_Triangle(struct triangulateio *outputCtx) 26*3a074057SBarry Smith { 27*3a074057SBarry Smith PetscFunctionBegin; 28*3a074057SBarry Smith outputCtx->numberofpoints = 0; 29*3a074057SBarry Smith outputCtx->pointlist = NULL; 30*3a074057SBarry Smith outputCtx->pointattributelist = NULL; 31*3a074057SBarry Smith outputCtx->pointmarkerlist = NULL; 32*3a074057SBarry Smith outputCtx->numberoftriangles = 0; 33*3a074057SBarry Smith outputCtx->trianglelist = NULL; 34*3a074057SBarry Smith outputCtx->triangleattributelist = NULL; 35*3a074057SBarry Smith outputCtx->neighborlist = NULL; 36*3a074057SBarry Smith outputCtx->segmentlist = NULL; 37*3a074057SBarry Smith outputCtx->segmentmarkerlist = NULL; 38*3a074057SBarry Smith outputCtx->numberofedges = 0; 39*3a074057SBarry Smith outputCtx->edgelist = NULL; 40*3a074057SBarry Smith outputCtx->edgemarkerlist = NULL; 41*3a074057SBarry Smith PetscFunctionReturn(0); 42*3a074057SBarry Smith } 43*3a074057SBarry Smith 44*3a074057SBarry Smith static PetscErrorCode FiniOutput_Triangle(struct triangulateio *outputCtx) 45*3a074057SBarry Smith { 46*3a074057SBarry Smith PetscFunctionBegin; 47*3a074057SBarry Smith free(outputCtx->pointlist); 48*3a074057SBarry Smith free(outputCtx->pointmarkerlist); 49*3a074057SBarry Smith free(outputCtx->segmentlist); 50*3a074057SBarry Smith free(outputCtx->segmentmarkerlist); 51*3a074057SBarry Smith free(outputCtx->edgelist); 52*3a074057SBarry Smith free(outputCtx->edgemarkerlist); 53*3a074057SBarry Smith free(outputCtx->trianglelist); 54*3a074057SBarry Smith free(outputCtx->neighborlist); 55*3a074057SBarry Smith PetscFunctionReturn(0); 56*3a074057SBarry Smith } 57*3a074057SBarry Smith 58*3a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerate_Triangle(DM boundary, PetscBool interpolate, DM *dm) 59*3a074057SBarry Smith { 60*3a074057SBarry Smith MPI_Comm comm; 61*3a074057SBarry Smith DM_Plex *mesh = (DM_Plex *) boundary->data; 62*3a074057SBarry Smith PetscInt dim = 2; 63*3a074057SBarry Smith const PetscBool createConvexHull = PETSC_FALSE; 64*3a074057SBarry Smith const PetscBool constrained = PETSC_FALSE; 65*3a074057SBarry Smith const char *labelName = "marker"; 66*3a074057SBarry Smith const char *labelName2 = "Face Sets"; 67*3a074057SBarry Smith struct triangulateio in; 68*3a074057SBarry Smith struct triangulateio out; 69*3a074057SBarry Smith DMLabel label, label2; 70*3a074057SBarry Smith PetscInt vStart, vEnd, v, eStart, eEnd, e; 71*3a074057SBarry Smith PetscMPIInt rank; 72*3a074057SBarry Smith PetscErrorCode ierr; 73*3a074057SBarry Smith 74*3a074057SBarry Smith PetscFunctionBegin; 75*3a074057SBarry Smith ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr); 76*3a074057SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 77*3a074057SBarry Smith ierr = InitInput_Triangle(&in);CHKERRQ(ierr); 78*3a074057SBarry Smith ierr = InitOutput_Triangle(&out);CHKERRQ(ierr); 79*3a074057SBarry Smith ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr); 80*3a074057SBarry Smith ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr); 81*3a074057SBarry Smith ierr = DMGetLabel(boundary, labelName2, &label2);CHKERRQ(ierr); 82*3a074057SBarry Smith 83*3a074057SBarry Smith in.numberofpoints = vEnd - vStart; 84*3a074057SBarry Smith if (in.numberofpoints > 0) { 85*3a074057SBarry Smith PetscSection coordSection; 86*3a074057SBarry Smith Vec coordinates; 87*3a074057SBarry Smith PetscScalar *array; 88*3a074057SBarry Smith 89*3a074057SBarry Smith ierr = PetscMalloc1(in.numberofpoints*dim, &in.pointlist);CHKERRQ(ierr); 90*3a074057SBarry Smith ierr = PetscMalloc1(in.numberofpoints, &in.pointmarkerlist);CHKERRQ(ierr); 91*3a074057SBarry Smith ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr); 92*3a074057SBarry Smith ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr); 93*3a074057SBarry Smith ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 94*3a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 95*3a074057SBarry Smith const PetscInt idx = v - vStart; 96*3a074057SBarry Smith PetscInt off, d; 97*3a074057SBarry Smith 98*3a074057SBarry Smith ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 99*3a074057SBarry Smith for (d = 0; d < dim; ++d) { 100*3a074057SBarry Smith in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 101*3a074057SBarry Smith } 102*3a074057SBarry Smith if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);} 103*3a074057SBarry Smith } 104*3a074057SBarry Smith ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 105*3a074057SBarry Smith } 106*3a074057SBarry Smith ierr = DMPlexGetHeightStratum(boundary, 0, &eStart, &eEnd);CHKERRQ(ierr); 107*3a074057SBarry Smith in.numberofsegments = eEnd - eStart; 108*3a074057SBarry Smith if (in.numberofsegments > 0) { 109*3a074057SBarry Smith ierr = PetscMalloc1(in.numberofsegments*2, &in.segmentlist);CHKERRQ(ierr); 110*3a074057SBarry Smith ierr = PetscMalloc1(in.numberofsegments, &in.segmentmarkerlist);CHKERRQ(ierr); 111*3a074057SBarry Smith for (e = eStart; e < eEnd; ++e) { 112*3a074057SBarry Smith const PetscInt idx = e - eStart; 113*3a074057SBarry Smith const PetscInt *cone; 114*3a074057SBarry Smith 115*3a074057SBarry Smith ierr = DMPlexGetCone(boundary, e, &cone);CHKERRQ(ierr); 116*3a074057SBarry Smith 117*3a074057SBarry Smith in.segmentlist[idx*2+0] = cone[0] - vStart; 118*3a074057SBarry Smith in.segmentlist[idx*2+1] = cone[1] - vStart; 119*3a074057SBarry Smith 120*3a074057SBarry Smith if (label) {ierr = DMLabelGetValue(label, e, &in.segmentmarkerlist[idx]);CHKERRQ(ierr);} 121*3a074057SBarry Smith } 122*3a074057SBarry Smith } 123*3a074057SBarry Smith #if 0 /* Do not currently support holes */ 124*3a074057SBarry Smith PetscReal *holeCoords; 125*3a074057SBarry Smith PetscInt h, d; 126*3a074057SBarry Smith 127*3a074057SBarry Smith ierr = DMPlexGetHoles(boundary, &in.numberofholes, &holeCords);CHKERRQ(ierr); 128*3a074057SBarry Smith if (in.numberofholes > 0) { 129*3a074057SBarry Smith ierr = PetscMalloc1(in.numberofholes*dim, &in.holelist);CHKERRQ(ierr); 130*3a074057SBarry Smith for (h = 0; h < in.numberofholes; ++h) { 131*3a074057SBarry Smith for (d = 0; d < dim; ++d) { 132*3a074057SBarry Smith in.holelist[h*dim+d] = holeCoords[h*dim+d]; 133*3a074057SBarry Smith } 134*3a074057SBarry Smith } 135*3a074057SBarry Smith } 136*3a074057SBarry Smith #endif 137*3a074057SBarry Smith if (!rank) { 138*3a074057SBarry Smith char args[32]; 139*3a074057SBarry Smith 140*3a074057SBarry Smith /* Take away 'Q' for verbose output */ 141*3a074057SBarry Smith ierr = PetscStrcpy(args, "pqezQ");CHKERRQ(ierr); 142*3a074057SBarry Smith if (createConvexHull) {ierr = PetscStrcat(args, "c");CHKERRQ(ierr);} 143*3a074057SBarry Smith if (constrained) {ierr = PetscStrcpy(args, "zepDQ");CHKERRQ(ierr);} 144*3a074057SBarry Smith if (mesh->triangleOpts) {triangulate(mesh->triangleOpts, &in, &out, NULL);} 145*3a074057SBarry Smith else {triangulate(args, &in, &out, NULL);} 146*3a074057SBarry Smith } 147*3a074057SBarry Smith ierr = PetscFree(in.pointlist);CHKERRQ(ierr); 148*3a074057SBarry Smith ierr = PetscFree(in.pointmarkerlist);CHKERRQ(ierr); 149*3a074057SBarry Smith ierr = PetscFree(in.segmentlist);CHKERRQ(ierr); 150*3a074057SBarry Smith ierr = PetscFree(in.segmentmarkerlist);CHKERRQ(ierr); 151*3a074057SBarry Smith ierr = PetscFree(in.holelist);CHKERRQ(ierr); 152*3a074057SBarry Smith 153*3a074057SBarry Smith { 154*3a074057SBarry Smith DMLabel glabel = NULL; 155*3a074057SBarry Smith DMLabel glabel2 = NULL; 156*3a074057SBarry Smith const PetscInt numCorners = 3; 157*3a074057SBarry Smith const PetscInt numCells = out.numberoftriangles; 158*3a074057SBarry Smith const PetscInt numVertices = out.numberofpoints; 159*3a074057SBarry Smith const int *cells = out.trianglelist; 160*3a074057SBarry Smith const double *meshCoords = out.pointlist; 161*3a074057SBarry Smith 162*3a074057SBarry Smith ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr); 163*3a074057SBarry Smith if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);} 164*3a074057SBarry Smith if (label2) {ierr = DMCreateLabel(*dm, labelName2); ierr = DMGetLabel(*dm, labelName2, &glabel2);} 165*3a074057SBarry Smith /* Set labels */ 166*3a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 167*3a074057SBarry Smith if (out.pointmarkerlist[v]) { 168*3a074057SBarry Smith if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);} 169*3a074057SBarry Smith } 170*3a074057SBarry Smith } 171*3a074057SBarry Smith if (interpolate) { 172*3a074057SBarry Smith for (e = 0; e < out.numberofedges; e++) { 173*3a074057SBarry Smith if (out.edgemarkerlist[e]) { 174*3a074057SBarry Smith const PetscInt vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells}; 175*3a074057SBarry Smith const PetscInt *edges; 176*3a074057SBarry Smith PetscInt numEdges; 177*3a074057SBarry Smith 178*3a074057SBarry Smith ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 179*3a074057SBarry Smith if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 180*3a074057SBarry Smith if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 181*3a074057SBarry Smith if (glabel2) {ierr = DMLabelSetValue(glabel2, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 182*3a074057SBarry Smith ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 183*3a074057SBarry Smith } 184*3a074057SBarry Smith } 185*3a074057SBarry Smith } 186*3a074057SBarry Smith ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr); 187*3a074057SBarry Smith } 188*3a074057SBarry Smith #if 0 /* Do not currently support holes */ 189*3a074057SBarry Smith ierr = DMPlexCopyHoles(*dm, boundary);CHKERRQ(ierr); 190*3a074057SBarry Smith #endif 191*3a074057SBarry Smith ierr = FiniOutput_Triangle(&out);CHKERRQ(ierr); 192*3a074057SBarry Smith PetscFunctionReturn(0); 193*3a074057SBarry Smith } 194*3a074057SBarry Smith 195*3a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexRefine_Triangle(DM dm, double *maxVolumes, DM *dmRefined) 196*3a074057SBarry Smith { 197*3a074057SBarry Smith MPI_Comm comm; 198*3a074057SBarry Smith PetscInt dim = 2; 199*3a074057SBarry Smith const char *labelName = "marker"; 200*3a074057SBarry Smith struct triangulateio in; 201*3a074057SBarry Smith struct triangulateio out; 202*3a074057SBarry Smith DMLabel label; 203*3a074057SBarry Smith PetscInt vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal; 204*3a074057SBarry Smith PetscMPIInt rank; 205*3a074057SBarry Smith PetscErrorCode ierr; 206*3a074057SBarry Smith 207*3a074057SBarry Smith PetscFunctionBegin; 208*3a074057SBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 209*3a074057SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 210*3a074057SBarry Smith ierr = InitInput_Triangle(&in);CHKERRQ(ierr); 211*3a074057SBarry Smith ierr = InitOutput_Triangle(&out);CHKERRQ(ierr); 212*3a074057SBarry Smith ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 213*3a074057SBarry Smith ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 214*3a074057SBarry Smith ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 215*3a074057SBarry Smith ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr); 216*3a074057SBarry Smith 217*3a074057SBarry Smith in.numberofpoints = vEnd - vStart; 218*3a074057SBarry Smith if (in.numberofpoints > 0) { 219*3a074057SBarry Smith PetscSection coordSection; 220*3a074057SBarry Smith Vec coordinates; 221*3a074057SBarry Smith PetscScalar *array; 222*3a074057SBarry Smith 223*3a074057SBarry Smith ierr = PetscMalloc1(in.numberofpoints*dim, &in.pointlist);CHKERRQ(ierr); 224*3a074057SBarry Smith ierr = PetscMalloc1(in.numberofpoints, &in.pointmarkerlist);CHKERRQ(ierr); 225*3a074057SBarry Smith ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 226*3a074057SBarry Smith ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 227*3a074057SBarry Smith ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 228*3a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 229*3a074057SBarry Smith const PetscInt idx = v - vStart; 230*3a074057SBarry Smith PetscInt off, d; 231*3a074057SBarry Smith 232*3a074057SBarry Smith ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 233*3a074057SBarry Smith for (d = 0; d < dim; ++d) { 234*3a074057SBarry Smith in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 235*3a074057SBarry Smith } 236*3a074057SBarry Smith if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);} 237*3a074057SBarry Smith } 238*3a074057SBarry Smith ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 239*3a074057SBarry Smith } 240*3a074057SBarry Smith ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 241*3a074057SBarry Smith 242*3a074057SBarry Smith in.numberofcorners = 3; 243*3a074057SBarry Smith in.numberoftriangles = cEnd - cStart; 244*3a074057SBarry Smith 245*3a074057SBarry Smith in.trianglearealist = (double*) maxVolumes; 246*3a074057SBarry Smith if (in.numberoftriangles > 0) { 247*3a074057SBarry Smith ierr = PetscMalloc1(in.numberoftriangles*in.numberofcorners, &in.trianglelist);CHKERRQ(ierr); 248*3a074057SBarry Smith for (c = cStart; c < cEnd; ++c) { 249*3a074057SBarry Smith const PetscInt idx = c - cStart; 250*3a074057SBarry Smith PetscInt *closure = NULL; 251*3a074057SBarry Smith PetscInt closureSize; 252*3a074057SBarry Smith 253*3a074057SBarry Smith ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 254*3a074057SBarry Smith if ((closureSize != 4) && (closureSize != 7)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a triangle, %D vertices in closure", closureSize); 255*3a074057SBarry Smith for (v = 0; v < 3; ++v) { 256*3a074057SBarry Smith in.trianglelist[idx*in.numberofcorners + v] = closure[(v+closureSize-3)*2] - vStart; 257*3a074057SBarry Smith } 258*3a074057SBarry Smith ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 259*3a074057SBarry Smith } 260*3a074057SBarry Smith } 261*3a074057SBarry Smith /* TODO: Segment markers are missing on input */ 262*3a074057SBarry Smith #if 0 /* Do not currently support holes */ 263*3a074057SBarry Smith PetscReal *holeCoords; 264*3a074057SBarry Smith PetscInt h, d; 265*3a074057SBarry Smith 266*3a074057SBarry Smith ierr = DMPlexGetHoles(boundary, &in.numberofholes, &holeCords);CHKERRQ(ierr); 267*3a074057SBarry Smith if (in.numberofholes > 0) { 268*3a074057SBarry Smith ierr = PetscMalloc1(in.numberofholes*dim, &in.holelist);CHKERRQ(ierr); 269*3a074057SBarry Smith for (h = 0; h < in.numberofholes; ++h) { 270*3a074057SBarry Smith for (d = 0; d < dim; ++d) { 271*3a074057SBarry Smith in.holelist[h*dim+d] = holeCoords[h*dim+d]; 272*3a074057SBarry Smith } 273*3a074057SBarry Smith } 274*3a074057SBarry Smith } 275*3a074057SBarry Smith #endif 276*3a074057SBarry Smith if (!rank) { 277*3a074057SBarry Smith char args[32]; 278*3a074057SBarry Smith 279*3a074057SBarry Smith /* Take away 'Q' for verbose output */ 280*3a074057SBarry Smith ierr = PetscStrcpy(args, "pqezQra");CHKERRQ(ierr); 281*3a074057SBarry Smith triangulate(args, &in, &out, NULL); 282*3a074057SBarry Smith } 283*3a074057SBarry Smith ierr = PetscFree(in.pointlist);CHKERRQ(ierr); 284*3a074057SBarry Smith ierr = PetscFree(in.pointmarkerlist);CHKERRQ(ierr); 285*3a074057SBarry Smith ierr = PetscFree(in.segmentlist);CHKERRQ(ierr); 286*3a074057SBarry Smith ierr = PetscFree(in.segmentmarkerlist);CHKERRQ(ierr); 287*3a074057SBarry Smith ierr = PetscFree(in.trianglelist);CHKERRQ(ierr); 288*3a074057SBarry Smith 289*3a074057SBarry Smith { 290*3a074057SBarry Smith DMLabel rlabel = NULL; 291*3a074057SBarry Smith const PetscInt numCorners = 3; 292*3a074057SBarry Smith const PetscInt numCells = out.numberoftriangles; 293*3a074057SBarry Smith const PetscInt numVertices = out.numberofpoints; 294*3a074057SBarry Smith const int *cells = out.trianglelist; 295*3a074057SBarry Smith const double *meshCoords = out.pointlist; 296*3a074057SBarry Smith PetscBool interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE; 297*3a074057SBarry Smith 298*3a074057SBarry Smith ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr); 299*3a074057SBarry Smith if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);} 300*3a074057SBarry Smith /* Set labels */ 301*3a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 302*3a074057SBarry Smith if (out.pointmarkerlist[v]) { 303*3a074057SBarry Smith if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);} 304*3a074057SBarry Smith } 305*3a074057SBarry Smith } 306*3a074057SBarry Smith if (interpolate) { 307*3a074057SBarry Smith PetscInt e; 308*3a074057SBarry Smith 309*3a074057SBarry Smith for (e = 0; e < out.numberofedges; e++) { 310*3a074057SBarry Smith if (out.edgemarkerlist[e]) { 311*3a074057SBarry Smith const PetscInt vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells}; 312*3a074057SBarry Smith const PetscInt *edges; 313*3a074057SBarry Smith PetscInt numEdges; 314*3a074057SBarry Smith 315*3a074057SBarry Smith ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 316*3a074057SBarry Smith if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 317*3a074057SBarry Smith if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 318*3a074057SBarry Smith ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 319*3a074057SBarry Smith } 320*3a074057SBarry Smith } 321*3a074057SBarry Smith } 322*3a074057SBarry Smith ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr); 323*3a074057SBarry Smith } 324*3a074057SBarry Smith #if 0 /* Do not currently support holes */ 325*3a074057SBarry Smith ierr = DMPlexCopyHoles(*dm, boundary);CHKERRQ(ierr); 326*3a074057SBarry Smith #endif 327*3a074057SBarry Smith ierr = FiniOutput_Triangle(&out);CHKERRQ(ierr); 328*3a074057SBarry Smith PetscFunctionReturn(0); 329*3a074057SBarry Smith } 330