xref: /petsc/src/dm/impls/plex/generators/tetgen/tetgenerate.cxx (revision b83bf3d7549f6efc01ccb7c977ade1611394aee8)
13a074057SBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
23a074057SBarry Smith 
3*b83bf3d7SVaclav Hapla #if defined(PETSC_HAVE_TETGEN_TETLIBRARY_NEEDED)
4*b83bf3d7SVaclav Hapla #define TETLIBRARY
5*b83bf3d7SVaclav Hapla #endif
63a074057SBarry Smith #include <tetgen.h>
73a074057SBarry Smith 
83a074057SBarry Smith /* This is to fix the tetrahedron orientation from TetGen */
9a4a685f2SJacob Faibussowitsch static PetscErrorCode DMPlexInvertCells_Tetgen(PetscInt numCells, PetscInt numCorners, PetscInt cells[])
103a074057SBarry Smith {
113a074057SBarry Smith   PetscInt bound = numCells*numCorners, coff;
123a074057SBarry Smith 
133a074057SBarry Smith   PetscFunctionBegin;
14a4a685f2SJacob Faibussowitsch #define SWAP(a,b) do { PetscInt tmp = (a); (a) = (b); (b) = tmp; } while (0)
1596ca5757SLisandro Dalcin   for (coff = 0; coff < bound; coff += numCorners) SWAP(cells[coff],cells[coff+1]);
1696ca5757SLisandro Dalcin #undef SWAP
173a074057SBarry Smith   PetscFunctionReturn(0);
183a074057SBarry Smith }
193a074057SBarry Smith 
203a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerate_Tetgen(DM boundary, PetscBool interpolate, DM *dm)
213a074057SBarry Smith {
223a074057SBarry Smith   MPI_Comm       comm;
233a074057SBarry Smith   DM_Plex       *mesh      = (DM_Plex *) boundary->data;
243a074057SBarry Smith   const PetscInt dim       = 3;
253a074057SBarry Smith   const char    *labelName = "marker";
263a074057SBarry Smith   ::tetgenio     in;
273a074057SBarry Smith   ::tetgenio     out;
283a074057SBarry Smith   DMLabel        label;
293a074057SBarry Smith   PetscInt       vStart, vEnd, v, fStart, fEnd, f;
303a074057SBarry Smith   PetscMPIInt    rank;
313a074057SBarry Smith   PetscErrorCode ierr;
323a074057SBarry Smith 
333a074057SBarry Smith   PetscFunctionBegin;
343a074057SBarry Smith   ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr);
353a074057SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
363a074057SBarry Smith   ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr);
373a074057SBarry Smith   ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr);
383a074057SBarry Smith 
393a074057SBarry Smith   in.numberofpoints = vEnd - vStart;
403a074057SBarry Smith   if (in.numberofpoints > 0) {
413a074057SBarry Smith     PetscSection coordSection;
423a074057SBarry Smith     Vec          coordinates;
433a074057SBarry Smith     PetscScalar *array;
443a074057SBarry Smith 
453a074057SBarry Smith     in.pointlist       = new double[in.numberofpoints*dim];
463a074057SBarry Smith     in.pointmarkerlist = new int[in.numberofpoints];
473a074057SBarry Smith 
483a074057SBarry Smith     ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr);
493a074057SBarry Smith     ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr);
503a074057SBarry Smith     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
513a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
523a074057SBarry Smith       const PetscInt idx = v - vStart;
533a074057SBarry Smith       PetscInt       off, d;
543a074057SBarry Smith 
553a074057SBarry Smith       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
563a074057SBarry Smith       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
573a074057SBarry Smith       if (label) {
583a074057SBarry Smith         PetscInt val;
593a074057SBarry Smith 
603a074057SBarry Smith         ierr = DMLabelGetValue(label, v, &val);CHKERRQ(ierr);
613a074057SBarry Smith         in.pointmarkerlist[idx] = (int) val;
623a074057SBarry Smith       }
633a074057SBarry Smith     }
643a074057SBarry Smith     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
653a074057SBarry Smith   }
663a074057SBarry Smith   ierr  = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr);
673a074057SBarry Smith 
683a074057SBarry Smith   in.numberoffacets = fEnd - fStart;
693a074057SBarry Smith   if (in.numberoffacets > 0) {
703a074057SBarry Smith     in.facetlist       = new tetgenio::facet[in.numberoffacets];
713a074057SBarry Smith     in.facetmarkerlist = new int[in.numberoffacets];
723a074057SBarry Smith     for (f = fStart; f < fEnd; ++f) {
733a074057SBarry Smith       const PetscInt idx     = f - fStart;
743a074057SBarry Smith       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v;
753a074057SBarry Smith 
763a074057SBarry Smith       in.facetlist[idx].numberofpolygons = 1;
773a074057SBarry Smith       in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
783a074057SBarry Smith       in.facetlist[idx].numberofholes    = 0;
793a074057SBarry Smith       in.facetlist[idx].holelist         = NULL;
803a074057SBarry Smith 
813a074057SBarry Smith       ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
823a074057SBarry Smith       for (p = 0; p < numPoints*2; p += 2) {
833a074057SBarry Smith         const PetscInt point = points[p];
843a074057SBarry Smith         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
853a074057SBarry Smith       }
863a074057SBarry Smith 
873a074057SBarry Smith       tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
883a074057SBarry Smith       poly->numberofvertices = numVertices;
893a074057SBarry Smith       poly->vertexlist       = new int[poly->numberofvertices];
903a074057SBarry Smith       for (v = 0; v < numVertices; ++v) {
913a074057SBarry Smith         const PetscInt vIdx = points[v] - vStart;
923a074057SBarry Smith         poly->vertexlist[v] = vIdx;
933a074057SBarry Smith       }
943a074057SBarry Smith       if (label) {
953a074057SBarry Smith         PetscInt val;
963a074057SBarry Smith 
973a074057SBarry Smith         ierr = DMLabelGetValue(label, f, &val);CHKERRQ(ierr);
983a074057SBarry Smith         in.facetmarkerlist[idx] = (int) val;
993a074057SBarry Smith       }
1003a074057SBarry Smith       ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
1013a074057SBarry Smith     }
1023a074057SBarry Smith   }
1033a074057SBarry Smith   if (!rank) {
1043a074057SBarry Smith     char args[32];
1053a074057SBarry Smith 
1063a074057SBarry Smith     /* Take away 'Q' for verbose output */
1073a074057SBarry Smith     ierr = PetscStrcpy(args, "pqezQ");CHKERRQ(ierr);
1083a074057SBarry Smith     if (mesh->tetgenOpts) {::tetrahedralize(mesh->tetgenOpts, &in, &out);}
1093a074057SBarry Smith     else                  {::tetrahedralize(args, &in, &out);}
1103a074057SBarry Smith   }
1113a074057SBarry Smith   {
1123a074057SBarry Smith     DMLabel          glabel      = NULL;
1133a074057SBarry Smith     const PetscInt   numCorners  = 4;
1143a074057SBarry Smith     const PetscInt   numCells    = out.numberoftetrahedra;
1153a074057SBarry Smith     const PetscInt   numVertices = out.numberofpoints;
116a4a685f2SJacob Faibussowitsch     PetscReal        *meshCoords = NULL;
117a4a685f2SJacob Faibussowitsch     PetscInt         *cells      = NULL;
118a4a685f2SJacob Faibussowitsch 
119a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) == sizeof (out.pointlist[0])) {
120a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *) out.pointlist;
121a4a685f2SJacob Faibussowitsch     } else {
122a4a685f2SJacob Faibussowitsch       PetscInt i;
123a4a685f2SJacob Faibussowitsch 
124a4a685f2SJacob Faibussowitsch       meshCoords = new PetscReal[dim * numVertices];
125a4a685f2SJacob Faibussowitsch       for (i = 0; i < dim * numVertices; i++) {
126a4a685f2SJacob Faibussowitsch         meshCoords[i] = (PetscReal) out.pointlist[i];
127a4a685f2SJacob Faibussowitsch       }
128a4a685f2SJacob Faibussowitsch     }
129a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) == sizeof (out.tetrahedronlist[0])) {
130a4a685f2SJacob Faibussowitsch       cells = (PetscInt *) out.tetrahedronlist;
131a4a685f2SJacob Faibussowitsch     } else {
132a4a685f2SJacob Faibussowitsch       PetscInt i;
133a4a685f2SJacob Faibussowitsch 
134a4a685f2SJacob Faibussowitsch       cells = new PetscInt[numCells * numCorners];
135a4a685f2SJacob Faibussowitsch       for (i = 0; i < numCells * numCorners; i++) {
136a4a685f2SJacob Faibussowitsch         cells[i] = (PetscInt) out.tetrahedronlist[i];
137a4a685f2SJacob Faibussowitsch       }
138a4a685f2SJacob Faibussowitsch     }
1393a074057SBarry Smith 
14096ca5757SLisandro Dalcin     ierr = DMPlexInvertCells_Tetgen(numCells, numCorners, cells);CHKERRQ(ierr);
141a4a685f2SJacob Faibussowitsch     ierr = DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr);
14296ca5757SLisandro Dalcin     if (label) {ierr = DMCreateLabel(*dm, labelName);CHKERRQ(ierr); ierr = DMGetLabel(*dm, labelName, &glabel);CHKERRQ(ierr);}
1433a074057SBarry Smith     /* Set labels */
1443a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
1453a074057SBarry Smith       if (out.pointmarkerlist[v]) {
1463a074057SBarry Smith         if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);}
1473a074057SBarry Smith       }
1483a074057SBarry Smith     }
1493a074057SBarry Smith     if (interpolate) {
1503a074057SBarry Smith #if 0
1513a074057SBarry Smith       PetscInt e;
1523a074057SBarry Smith 
1533a074057SBarry Smith       /* This check is never actually executed for ctetgen (which never returns edgemarkers) and seems to be broken for
1543a074057SBarry Smith        * tetgen */
1553a074057SBarry Smith       for (e = 0; e < out.numberofedges; e++) {
1563a074057SBarry Smith         if (out.edgemarkerlist[e]) {
1573a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
1583a074057SBarry Smith           const PetscInt *edges;
1593a074057SBarry Smith           PetscInt        numEdges;
1603a074057SBarry Smith 
1613a074057SBarry Smith           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1623a074057SBarry Smith           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
1633a074057SBarry Smith           if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
1643a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1653a074057SBarry Smith         }
1663a074057SBarry Smith       }
1673a074057SBarry Smith #endif
1683a074057SBarry Smith       for (f = 0; f < out.numberoftrifaces; f++) {
1693a074057SBarry Smith         if (out.trifacemarkerlist[f]) {
1703a074057SBarry Smith           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
1713a074057SBarry Smith           const PetscInt *faces;
1723a074057SBarry Smith           PetscInt        numFaces;
1733a074057SBarry Smith 
1743a074057SBarry Smith           ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1753a074057SBarry Smith           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
1763a074057SBarry Smith           if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);}
1773a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1783a074057SBarry Smith         }
1793a074057SBarry Smith       }
1803a074057SBarry Smith     }
1813a074057SBarry Smith     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
1823a074057SBarry Smith   }
1833a074057SBarry Smith   PetscFunctionReturn(0);
1843a074057SBarry Smith }
1853a074057SBarry Smith 
1863a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexRefine_Tetgen(DM dm, double *maxVolumes, DM *dmRefined)
1873a074057SBarry Smith {
1883a074057SBarry Smith   MPI_Comm       comm;
1893a074057SBarry Smith   const PetscInt dim       = 3;
1903a074057SBarry Smith   const char    *labelName = "marker";
1913a074057SBarry Smith   ::tetgenio     in;
1923a074057SBarry Smith   ::tetgenio     out;
1933a074057SBarry Smith   DMLabel        label;
1943a074057SBarry Smith   PetscInt       vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal;
1953a074057SBarry Smith   PetscMPIInt    rank;
1963a074057SBarry Smith   PetscErrorCode ierr;
1973a074057SBarry Smith 
1983a074057SBarry Smith   PetscFunctionBegin;
1993a074057SBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
2003a074057SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2013a074057SBarry Smith   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2023a074057SBarry Smith   ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
2033a074057SBarry Smith   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
2043a074057SBarry Smith   ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr);
2053a074057SBarry Smith 
2063a074057SBarry Smith   in.numberofpoints = vEnd - vStart;
2073a074057SBarry Smith   if (in.numberofpoints > 0) {
2083a074057SBarry Smith     PetscSection coordSection;
2093a074057SBarry Smith     Vec          coordinates;
2103a074057SBarry Smith     PetscScalar *array;
2113a074057SBarry Smith 
2123a074057SBarry Smith     in.pointlist       = new double[in.numberofpoints*dim];
2133a074057SBarry Smith     in.pointmarkerlist = new int[in.numberofpoints];
2143a074057SBarry Smith 
2153a074057SBarry Smith     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2163a074057SBarry Smith     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2173a074057SBarry Smith     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
2183a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
2193a074057SBarry Smith       const PetscInt idx = v - vStart;
2203a074057SBarry Smith       PetscInt       off, d;
2213a074057SBarry Smith 
2223a074057SBarry Smith       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
2233a074057SBarry Smith       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
2243a074057SBarry Smith       if (label) {
2253a074057SBarry Smith         PetscInt val;
2263a074057SBarry Smith 
2273a074057SBarry Smith         ierr = DMLabelGetValue(label, v, &val);CHKERRQ(ierr);
2283a074057SBarry Smith         in.pointmarkerlist[idx] = (int) val;
2293a074057SBarry Smith       }
2303a074057SBarry Smith     }
2313a074057SBarry Smith     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
2323a074057SBarry Smith   }
2333a074057SBarry Smith   ierr  = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2343a074057SBarry Smith 
2353a074057SBarry Smith   in.numberofcorners       = 4;
2363a074057SBarry Smith   in.numberoftetrahedra    = cEnd - cStart;
2373a074057SBarry Smith   in.tetrahedronvolumelist = (double*) maxVolumes;
2383a074057SBarry Smith   if (in.numberoftetrahedra > 0) {
2393a074057SBarry Smith     in.tetrahedronlist = new int[in.numberoftetrahedra*in.numberofcorners];
2403a074057SBarry Smith     for (c = cStart; c < cEnd; ++c) {
2413a074057SBarry Smith       const PetscInt idx      = c - cStart;
2423a074057SBarry Smith       PetscInt      *closure = NULL;
2433a074057SBarry Smith       PetscInt       closureSize;
2443a074057SBarry Smith 
2453a074057SBarry Smith       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
2463a074057SBarry 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);
2473a074057SBarry Smith       for (v = 0; v < 4; ++v) {
2483a074057SBarry Smith         in.tetrahedronlist[idx*in.numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
2493a074057SBarry Smith       }
2503a074057SBarry Smith       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
2513a074057SBarry Smith     }
2523a074057SBarry Smith   }
2533a074057SBarry Smith   /* TODO: Put in boundary faces with markers */
2543a074057SBarry Smith   if (!rank) {
2553a074057SBarry Smith     char args[32];
2563a074057SBarry Smith 
2573a074057SBarry Smith #if 1
2583a074057SBarry Smith     /* Take away 'Q' for verbose output */
2593a074057SBarry Smith     ierr = PetscStrcpy(args, "qezQra");CHKERRQ(ierr);
2603a074057SBarry Smith #else
2613a074057SBarry Smith     ierr = PetscStrcpy(args, "qezraVVVV");CHKERRQ(ierr);
2623a074057SBarry Smith #endif
2633a074057SBarry Smith     ::tetrahedralize(args, &in, &out);
2643a074057SBarry Smith   }
2653a074057SBarry Smith   in.tetrahedronvolumelist = NULL;
2663a074057SBarry Smith 
2673a074057SBarry Smith   {
2683a074057SBarry Smith     DMLabel          rlabel      = NULL;
2693a074057SBarry Smith     const PetscInt   numCorners  = 4;
2703a074057SBarry Smith     const PetscInt   numCells    = out.numberoftetrahedra;
2713a074057SBarry Smith     const PetscInt   numVertices = out.numberofpoints;
272a4a685f2SJacob Faibussowitsch     PetscReal        *meshCoords = NULL;
273a4a685f2SJacob Faibussowitsch     PetscInt         *cells      = NULL;
2743a074057SBarry Smith     PetscBool        interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
2753a074057SBarry Smith 
276a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) == sizeof (out.pointlist[0])) {
277a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *) out.pointlist;
278a4a685f2SJacob Faibussowitsch     } else {
279a4a685f2SJacob Faibussowitsch       PetscInt i;
280a4a685f2SJacob Faibussowitsch 
281a4a685f2SJacob Faibussowitsch       meshCoords = new PetscReal[dim * numVertices];
282a4a685f2SJacob Faibussowitsch       for (i = 0; i < dim * numVertices; i++) {
283a4a685f2SJacob Faibussowitsch         meshCoords[i] = (PetscReal) out.pointlist[i];
284a4a685f2SJacob Faibussowitsch       }
285a4a685f2SJacob Faibussowitsch     }
286a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) == sizeof (out.tetrahedronlist[0])) {
287a4a685f2SJacob Faibussowitsch       cells = (PetscInt *) out.tetrahedronlist;
288a4a685f2SJacob Faibussowitsch     } else {
289a4a685f2SJacob Faibussowitsch       PetscInt i;
290a4a685f2SJacob Faibussowitsch 
291a4a685f2SJacob Faibussowitsch       cells = new PetscInt[numCells * numCorners];
292a4a685f2SJacob Faibussowitsch       for (i = 0; i < numCells * numCorners; i++) {
293a4a685f2SJacob Faibussowitsch         cells[i] = (PetscInt) out.tetrahedronlist[i];
294a4a685f2SJacob Faibussowitsch       }
295a4a685f2SJacob Faibussowitsch     }
296a4a685f2SJacob Faibussowitsch 
29796ca5757SLisandro Dalcin     ierr = DMPlexInvertCells_Tetgen(numCells, numCorners, cells);CHKERRQ(ierr);
298a4a685f2SJacob Faibussowitsch     ierr = DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
299a4a685f2SJacob Faibussowitsch     if (label) {
300a4a685f2SJacob Faibussowitsch       ierr = DMCreateLabel(*dmRefined, labelName);CHKERRQ(ierr);
301a4a685f2SJacob Faibussowitsch       ierr = DMGetLabel(*dmRefined, labelName, &rlabel);CHKERRQ(ierr);
302a4a685f2SJacob Faibussowitsch     }
3033a074057SBarry Smith     /* Set labels */
3043a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
3053a074057SBarry Smith       if (out.pointmarkerlist[v]) {
3063a074057SBarry Smith         if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);}
3073a074057SBarry Smith       }
3083a074057SBarry Smith     }
3093a074057SBarry Smith     if (interpolate) {
3103a074057SBarry Smith       PetscInt f;
3113a074057SBarry Smith #if 0
3123a074057SBarry Smith       PetscInt e;
3133a074057SBarry Smith 
3143a074057SBarry Smith       for (e = 0; e < out.numberofedges; e++) {
3153a074057SBarry Smith         if (out.edgemarkerlist[e]) {
3163a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
3173a074057SBarry Smith           const PetscInt *edges;
3183a074057SBarry Smith           PetscInt        numEdges;
3193a074057SBarry Smith 
3203a074057SBarry Smith           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
3213a074057SBarry Smith           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
3223a074057SBarry Smith           if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
3233a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
3243a074057SBarry Smith         }
3253a074057SBarry Smith       }
3263a074057SBarry Smith #endif
3273a074057SBarry Smith       for (f = 0; f < out.numberoftrifaces; f++) {
3283a074057SBarry Smith         if (out.trifacemarkerlist[f]) {
3293a074057SBarry Smith           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
3303a074057SBarry Smith           const PetscInt *faces;
3313a074057SBarry Smith           PetscInt        numFaces;
3323a074057SBarry Smith 
3333a074057SBarry Smith           ierr = DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
3343a074057SBarry Smith           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
3353a074057SBarry Smith           if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);}
3363a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
3373a074057SBarry Smith         }
3383a074057SBarry Smith       }
3393a074057SBarry Smith     }
3403a074057SBarry Smith     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
3413a074057SBarry Smith   }
3423a074057SBarry Smith   PetscFunctionReturn(0);
3433a074057SBarry Smith }
344