13a074057SBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 23a074057SBarry Smith 3*3316697fSBarry Smith #if !defined(ANSI_DECLARATORS) 4*3316697fSBarry Smith #define ANSI_DECLARATORS 5*3316697fSBarry Smith #endif 63a074057SBarry Smith #include <triangle.h> 73a074057SBarry Smith 83a074057SBarry Smith static PetscErrorCode InitInput_Triangle(struct triangulateio *inputCtx) 93a074057SBarry Smith { 103a074057SBarry Smith PetscFunctionBegin; 113a074057SBarry Smith inputCtx->numberofpoints = 0; 123a074057SBarry Smith inputCtx->numberofpointattributes = 0; 133a074057SBarry Smith inputCtx->pointlist = NULL; 143a074057SBarry Smith inputCtx->pointattributelist = NULL; 153a074057SBarry Smith inputCtx->pointmarkerlist = NULL; 163a074057SBarry Smith inputCtx->numberofsegments = 0; 173a074057SBarry Smith inputCtx->segmentlist = NULL; 183a074057SBarry Smith inputCtx->segmentmarkerlist = NULL; 193a074057SBarry Smith inputCtx->numberoftriangleattributes = 0; 203a074057SBarry Smith inputCtx->trianglelist = NULL; 213a074057SBarry Smith inputCtx->numberofholes = 0; 223a074057SBarry Smith inputCtx->holelist = NULL; 233a074057SBarry Smith inputCtx->numberofregions = 0; 243a074057SBarry Smith inputCtx->regionlist = NULL; 253a074057SBarry Smith PetscFunctionReturn(0); 263a074057SBarry Smith } 273a074057SBarry Smith 283a074057SBarry Smith static PetscErrorCode InitOutput_Triangle(struct triangulateio *outputCtx) 293a074057SBarry Smith { 303a074057SBarry Smith PetscFunctionBegin; 313a074057SBarry Smith outputCtx->numberofpoints = 0; 323a074057SBarry Smith outputCtx->pointlist = NULL; 333a074057SBarry Smith outputCtx->pointattributelist = NULL; 343a074057SBarry Smith outputCtx->pointmarkerlist = NULL; 353a074057SBarry Smith outputCtx->numberoftriangles = 0; 363a074057SBarry Smith outputCtx->trianglelist = NULL; 373a074057SBarry Smith outputCtx->triangleattributelist = NULL; 383a074057SBarry Smith outputCtx->neighborlist = NULL; 393a074057SBarry Smith outputCtx->segmentlist = NULL; 403a074057SBarry Smith outputCtx->segmentmarkerlist = NULL; 413a074057SBarry Smith outputCtx->numberofedges = 0; 423a074057SBarry Smith outputCtx->edgelist = NULL; 433a074057SBarry Smith outputCtx->edgemarkerlist = NULL; 443a074057SBarry Smith PetscFunctionReturn(0); 453a074057SBarry Smith } 463a074057SBarry Smith 473a074057SBarry Smith static PetscErrorCode FiniOutput_Triangle(struct triangulateio *outputCtx) 483a074057SBarry Smith { 493a074057SBarry Smith PetscFunctionBegin; 503a074057SBarry Smith free(outputCtx->pointlist); 513a074057SBarry Smith free(outputCtx->pointmarkerlist); 523a074057SBarry Smith free(outputCtx->segmentlist); 533a074057SBarry Smith free(outputCtx->segmentmarkerlist); 543a074057SBarry Smith free(outputCtx->edgelist); 553a074057SBarry Smith free(outputCtx->edgemarkerlist); 563a074057SBarry Smith free(outputCtx->trianglelist); 573a074057SBarry Smith free(outputCtx->neighborlist); 583a074057SBarry Smith PetscFunctionReturn(0); 593a074057SBarry Smith } 603a074057SBarry Smith 613a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerate_Triangle(DM boundary, PetscBool interpolate, DM *dm) 623a074057SBarry Smith { 633a074057SBarry Smith MPI_Comm comm; 643a074057SBarry Smith DM_Plex *mesh = (DM_Plex *) boundary->data; 653a074057SBarry Smith PetscInt dim = 2; 663a074057SBarry Smith const PetscBool createConvexHull = PETSC_FALSE; 673a074057SBarry Smith const PetscBool constrained = PETSC_FALSE; 683a074057SBarry Smith const char *labelName = "marker"; 693a074057SBarry Smith const char *labelName2 = "Face Sets"; 703a074057SBarry Smith struct triangulateio in; 713a074057SBarry Smith struct triangulateio out; 723a074057SBarry Smith DMLabel label, label2; 733a074057SBarry Smith PetscInt vStart, vEnd, v, eStart, eEnd, e; 743a074057SBarry Smith PetscMPIInt rank; 753a074057SBarry Smith PetscErrorCode ierr; 763a074057SBarry Smith 773a074057SBarry Smith PetscFunctionBegin; 783a074057SBarry Smith ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr); 793a074057SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 803a074057SBarry Smith ierr = InitInput_Triangle(&in);CHKERRQ(ierr); 813a074057SBarry Smith ierr = InitOutput_Triangle(&out);CHKERRQ(ierr); 823a074057SBarry Smith ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr); 833a074057SBarry Smith ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr); 843a074057SBarry Smith ierr = DMGetLabel(boundary, labelName2, &label2);CHKERRQ(ierr); 853a074057SBarry Smith 863a074057SBarry Smith in.numberofpoints = vEnd - vStart; 873a074057SBarry Smith if (in.numberofpoints > 0) { 883a074057SBarry Smith PetscSection coordSection; 893a074057SBarry Smith Vec coordinates; 903a074057SBarry Smith PetscScalar *array; 913a074057SBarry Smith 923a074057SBarry Smith ierr = PetscMalloc1(in.numberofpoints*dim, &in.pointlist);CHKERRQ(ierr); 933a074057SBarry Smith ierr = PetscMalloc1(in.numberofpoints, &in.pointmarkerlist);CHKERRQ(ierr); 943a074057SBarry Smith ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr); 953a074057SBarry Smith ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr); 963a074057SBarry Smith ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 973a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 983a074057SBarry Smith const PetscInt idx = v - vStart; 99469e3fe5SMatthew G. Knepley PetscInt val, off, d; 1003a074057SBarry Smith 1013a074057SBarry Smith ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 1023a074057SBarry Smith for (d = 0; d < dim; ++d) { 1033a074057SBarry Smith in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 1043a074057SBarry Smith } 105469e3fe5SMatthew G. Knepley if (label) { 106469e3fe5SMatthew G. Knepley ierr = DMLabelGetValue(label, v, &val);CHKERRQ(ierr); 107469e3fe5SMatthew G. Knepley in.pointmarkerlist[idx] = val; 108469e3fe5SMatthew G. Knepley } 1093a074057SBarry Smith } 1103a074057SBarry Smith ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 1113a074057SBarry Smith } 1123a074057SBarry Smith ierr = DMPlexGetHeightStratum(boundary, 0, &eStart, &eEnd);CHKERRQ(ierr); 1133a074057SBarry Smith in.numberofsegments = eEnd - eStart; 1143a074057SBarry Smith if (in.numberofsegments > 0) { 1153a074057SBarry Smith ierr = PetscMalloc1(in.numberofsegments*2, &in.segmentlist);CHKERRQ(ierr); 1163a074057SBarry Smith ierr = PetscMalloc1(in.numberofsegments, &in.segmentmarkerlist);CHKERRQ(ierr); 1173a074057SBarry Smith for (e = eStart; e < eEnd; ++e) { 1183a074057SBarry Smith const PetscInt idx = e - eStart; 1193a074057SBarry Smith const PetscInt *cone; 120469e3fe5SMatthew G. Knepley PetscInt val; 1213a074057SBarry Smith 1223a074057SBarry Smith ierr = DMPlexGetCone(boundary, e, &cone);CHKERRQ(ierr); 1233a074057SBarry Smith 1243a074057SBarry Smith in.segmentlist[idx*2+0] = cone[0] - vStart; 1253a074057SBarry Smith in.segmentlist[idx*2+1] = cone[1] - vStart; 1263a074057SBarry Smith 127469e3fe5SMatthew G. Knepley if (label) { 128469e3fe5SMatthew G. Knepley ierr = DMLabelGetValue(label, e, &val);CHKERRQ(ierr); 129469e3fe5SMatthew G. Knepley in.segmentmarkerlist[idx] = val; 130469e3fe5SMatthew G. Knepley } 1313a074057SBarry Smith } 1323a074057SBarry Smith } 1333a074057SBarry Smith #if 0 /* Do not currently support holes */ 1343a074057SBarry Smith PetscReal *holeCoords; 1353a074057SBarry Smith PetscInt h, d; 1363a074057SBarry Smith 1373a074057SBarry Smith ierr = DMPlexGetHoles(boundary, &in.numberofholes, &holeCords);CHKERRQ(ierr); 1383a074057SBarry Smith if (in.numberofholes > 0) { 1393a074057SBarry Smith ierr = PetscMalloc1(in.numberofholes*dim, &in.holelist);CHKERRQ(ierr); 1403a074057SBarry Smith for (h = 0; h < in.numberofholes; ++h) { 1413a074057SBarry Smith for (d = 0; d < dim; ++d) { 1423a074057SBarry Smith in.holelist[h*dim+d] = holeCoords[h*dim+d]; 1433a074057SBarry Smith } 1443a074057SBarry Smith } 1453a074057SBarry Smith } 1463a074057SBarry Smith #endif 1473a074057SBarry Smith if (!rank) { 1483a074057SBarry Smith char args[32]; 1493a074057SBarry Smith 1503a074057SBarry Smith /* Take away 'Q' for verbose output */ 1513a074057SBarry Smith ierr = PetscStrcpy(args, "pqezQ");CHKERRQ(ierr); 1523a074057SBarry Smith if (createConvexHull) {ierr = PetscStrcat(args, "c");CHKERRQ(ierr);} 1533a074057SBarry Smith if (constrained) {ierr = PetscStrcpy(args, "zepDQ");CHKERRQ(ierr);} 1543a074057SBarry Smith if (mesh->triangleOpts) {triangulate(mesh->triangleOpts, &in, &out, NULL);} 1553a074057SBarry Smith else {triangulate(args, &in, &out, NULL);} 1563a074057SBarry Smith } 1573a074057SBarry Smith ierr = PetscFree(in.pointlist);CHKERRQ(ierr); 1583a074057SBarry Smith ierr = PetscFree(in.pointmarkerlist);CHKERRQ(ierr); 1593a074057SBarry Smith ierr = PetscFree(in.segmentlist);CHKERRQ(ierr); 1603a074057SBarry Smith ierr = PetscFree(in.segmentmarkerlist);CHKERRQ(ierr); 1613a074057SBarry Smith ierr = PetscFree(in.holelist);CHKERRQ(ierr); 1623a074057SBarry Smith 1633a074057SBarry Smith { 1643a074057SBarry Smith DMLabel glabel = NULL; 1653a074057SBarry Smith DMLabel glabel2 = NULL; 1663a074057SBarry Smith const PetscInt numCorners = 3; 1673a074057SBarry Smith const PetscInt numCells = out.numberoftriangles; 1683a074057SBarry Smith const PetscInt numVertices = out.numberofpoints; 1693a074057SBarry Smith const int *cells = out.trianglelist; 1703a074057SBarry Smith const double *meshCoords = out.pointlist; 1713a074057SBarry Smith 1723a074057SBarry Smith ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr); 1733a074057SBarry Smith if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);} 1743a074057SBarry Smith if (label2) {ierr = DMCreateLabel(*dm, labelName2); ierr = DMGetLabel(*dm, labelName2, &glabel2);} 1753a074057SBarry Smith /* Set labels */ 1763a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 1773a074057SBarry Smith if (out.pointmarkerlist[v]) { 1783a074057SBarry Smith if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);} 1793a074057SBarry Smith } 1803a074057SBarry Smith } 1813a074057SBarry Smith if (interpolate) { 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 1883a074057SBarry Smith ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 1893a074057SBarry Smith if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 1903a074057SBarry Smith if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 1913a074057SBarry Smith if (glabel2) {ierr = DMLabelSetValue(glabel2, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 1923a074057SBarry Smith ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 1933a074057SBarry Smith } 1943a074057SBarry Smith } 1953a074057SBarry Smith } 1963a074057SBarry Smith ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr); 1973a074057SBarry Smith } 1983a074057SBarry Smith #if 0 /* Do not currently support holes */ 1993a074057SBarry Smith ierr = DMPlexCopyHoles(*dm, boundary);CHKERRQ(ierr); 2003a074057SBarry Smith #endif 2013a074057SBarry Smith ierr = FiniOutput_Triangle(&out);CHKERRQ(ierr); 2023a074057SBarry Smith PetscFunctionReturn(0); 2033a074057SBarry Smith } 2043a074057SBarry Smith 205800b850cSBarry Smith PETSC_EXTERN PetscErrorCode DMPlexRefine_Triangle(DM dm, PetscReal *inmaxVolumes, DM *dmRefined) 2063a074057SBarry Smith { 2073a074057SBarry Smith MPI_Comm comm; 2083a074057SBarry Smith PetscInt dim = 2; 2093a074057SBarry Smith const char *labelName = "marker"; 2103a074057SBarry Smith struct triangulateio in; 2113a074057SBarry Smith struct triangulateio out; 2123a074057SBarry Smith DMLabel label; 2133a074057SBarry Smith PetscInt vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal; 2143a074057SBarry Smith PetscMPIInt rank; 2153a074057SBarry Smith PetscErrorCode ierr; 216800b850cSBarry Smith double *maxVolumes; 2173a074057SBarry Smith 2183a074057SBarry Smith PetscFunctionBegin; 2193a074057SBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2203a074057SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 2213a074057SBarry Smith ierr = InitInput_Triangle(&in);CHKERRQ(ierr); 2223a074057SBarry Smith ierr = InitOutput_Triangle(&out);CHKERRQ(ierr); 2233a074057SBarry Smith ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2243a074057SBarry Smith ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 2253a074057SBarry Smith ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 2263a074057SBarry Smith ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr); 2273a074057SBarry Smith 2283a074057SBarry Smith in.numberofpoints = vEnd - vStart; 2293a074057SBarry Smith if (in.numberofpoints > 0) { 2303a074057SBarry Smith PetscSection coordSection; 2313a074057SBarry Smith Vec coordinates; 2323a074057SBarry Smith PetscScalar *array; 2333a074057SBarry Smith 2343a074057SBarry Smith ierr = PetscMalloc1(in.numberofpoints*dim, &in.pointlist);CHKERRQ(ierr); 2353a074057SBarry Smith ierr = PetscMalloc1(in.numberofpoints, &in.pointmarkerlist);CHKERRQ(ierr); 2363a074057SBarry Smith ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2373a074057SBarry Smith ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2383a074057SBarry Smith ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 2393a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 2403a074057SBarry Smith const PetscInt idx = v - vStart; 241469e3fe5SMatthew G. Knepley PetscInt off, d, val; 2423a074057SBarry Smith 2433a074057SBarry Smith ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 2443a074057SBarry Smith for (d = 0; d < dim; ++d) { 2453a074057SBarry Smith in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 2463a074057SBarry Smith } 247469e3fe5SMatthew G. Knepley if (label) { 248469e3fe5SMatthew G. Knepley ierr = DMLabelGetValue(label, v, &val);CHKERRQ(ierr); 249469e3fe5SMatthew G. Knepley in.pointmarkerlist[idx] = val; 250469e3fe5SMatthew G. Knepley } 2513a074057SBarry Smith } 2523a074057SBarry Smith ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 2533a074057SBarry Smith } 2543a074057SBarry Smith ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2553a074057SBarry Smith 2563a074057SBarry Smith in.numberofcorners = 3; 2573a074057SBarry Smith in.numberoftriangles = cEnd - cStart; 2583a074057SBarry Smith 259800b850cSBarry Smith #if !defined(PETSC_USE_REAL_DOUBLE) 260800b850cSBarry Smith ierr = PetscMalloc1(cEnd - cStart,&maxVolumes);CHKERRQ(ierr); 261800b850cSBarry Smith for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = (double)inmaxVolumes[c]; 262800b850cSBarry Smith #else 263800b850cSBarry Smith maxVolumes = inmaxVolumes; 264800b850cSBarry Smith #endif 265800b850cSBarry Smith 2663a074057SBarry Smith in.trianglearealist = (double*) maxVolumes; 2673a074057SBarry Smith if (in.numberoftriangles > 0) { 2683a074057SBarry Smith ierr = PetscMalloc1(in.numberoftriangles*in.numberofcorners, &in.trianglelist);CHKERRQ(ierr); 2693a074057SBarry Smith for (c = cStart; c < cEnd; ++c) { 2703a074057SBarry Smith const PetscInt idx = c - cStart; 2713a074057SBarry Smith PetscInt *closure = NULL; 2723a074057SBarry Smith PetscInt closureSize; 2733a074057SBarry Smith 2743a074057SBarry Smith ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 2753a074057SBarry 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); 2763a074057SBarry Smith for (v = 0; v < 3; ++v) { 2773a074057SBarry Smith in.trianglelist[idx*in.numberofcorners + v] = closure[(v+closureSize-3)*2] - vStart; 2783a074057SBarry Smith } 2793a074057SBarry Smith ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 2803a074057SBarry Smith } 2813a074057SBarry Smith } 2823a074057SBarry Smith /* TODO: Segment markers are missing on input */ 2833a074057SBarry Smith #if 0 /* Do not currently support holes */ 2843a074057SBarry Smith PetscReal *holeCoords; 2853a074057SBarry Smith PetscInt h, d; 2863a074057SBarry Smith 2873a074057SBarry Smith ierr = DMPlexGetHoles(boundary, &in.numberofholes, &holeCords);CHKERRQ(ierr); 2883a074057SBarry Smith if (in.numberofholes > 0) { 2893a074057SBarry Smith ierr = PetscMalloc1(in.numberofholes*dim, &in.holelist);CHKERRQ(ierr); 2903a074057SBarry Smith for (h = 0; h < in.numberofholes; ++h) { 2913a074057SBarry Smith for (d = 0; d < dim; ++d) { 2923a074057SBarry Smith in.holelist[h*dim+d] = holeCoords[h*dim+d]; 2933a074057SBarry Smith } 2943a074057SBarry Smith } 2953a074057SBarry Smith } 2963a074057SBarry Smith #endif 2973a074057SBarry Smith if (!rank) { 2983a074057SBarry Smith char args[32]; 2993a074057SBarry Smith 3003a074057SBarry Smith /* Take away 'Q' for verbose output */ 3013a074057SBarry Smith ierr = PetscStrcpy(args, "pqezQra");CHKERRQ(ierr); 3023a074057SBarry Smith triangulate(args, &in, &out, NULL); 3033a074057SBarry Smith } 3043a074057SBarry Smith ierr = PetscFree(in.pointlist);CHKERRQ(ierr); 3053a074057SBarry Smith ierr = PetscFree(in.pointmarkerlist);CHKERRQ(ierr); 3063a074057SBarry Smith ierr = PetscFree(in.segmentlist);CHKERRQ(ierr); 3073a074057SBarry Smith ierr = PetscFree(in.segmentmarkerlist);CHKERRQ(ierr); 3083a074057SBarry Smith ierr = PetscFree(in.trianglelist);CHKERRQ(ierr); 3093a074057SBarry Smith 3103a074057SBarry Smith { 3113a074057SBarry Smith DMLabel rlabel = NULL; 3123a074057SBarry Smith const PetscInt numCorners = 3; 3133a074057SBarry Smith const PetscInt numCells = out.numberoftriangles; 3143a074057SBarry Smith const PetscInt numVertices = out.numberofpoints; 3153a074057SBarry Smith const int *cells = out.trianglelist; 3163a074057SBarry Smith const double *meshCoords = out.pointlist; 3173a074057SBarry Smith PetscBool interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE; 3183a074057SBarry Smith 3193a074057SBarry Smith ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr); 3203a074057SBarry Smith if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);} 3213a074057SBarry Smith /* Set labels */ 3223a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 3233a074057SBarry Smith if (out.pointmarkerlist[v]) { 3243a074057SBarry Smith if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);} 3253a074057SBarry Smith } 3263a074057SBarry Smith } 3273a074057SBarry Smith if (interpolate) { 3283a074057SBarry Smith PetscInt e; 3293a074057SBarry Smith 3303a074057SBarry Smith for (e = 0; e < out.numberofedges; e++) { 3313a074057SBarry Smith if (out.edgemarkerlist[e]) { 3323a074057SBarry Smith const PetscInt vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells}; 3333a074057SBarry Smith const PetscInt *edges; 3343a074057SBarry Smith PetscInt numEdges; 3353a074057SBarry Smith 3363a074057SBarry Smith ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 3373a074057SBarry Smith if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 3383a074057SBarry Smith if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 3393a074057SBarry Smith ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 3403a074057SBarry Smith } 3413a074057SBarry Smith } 3423a074057SBarry Smith } 3433a074057SBarry Smith ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr); 3443a074057SBarry Smith } 3453a074057SBarry Smith #if 0 /* Do not currently support holes */ 3463a074057SBarry Smith ierr = DMPlexCopyHoles(*dm, boundary);CHKERRQ(ierr); 3473a074057SBarry Smith #endif 3483a074057SBarry Smith ierr = FiniOutput_Triangle(&out);CHKERRQ(ierr); 349800b850cSBarry Smith #if !defined(PETSC_USE_REAL_DOUBLE) 350800b850cSBarry Smith ierr = PetscFree(maxVolumes);CHKERRQ(ierr); 351800b850cSBarry Smith #endif 3523a074057SBarry Smith PetscFunctionReturn(0); 3533a074057SBarry Smith } 354