xref: /petsc/src/dm/impls/plex/generators/ctetgen/ctetgenerate.c (revision a4a685f2b656679e5de30bfff53536422329f553)
13a074057SBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
23a074057SBarry Smith 
33a074057SBarry Smith #include <ctetgen.h>
43a074057SBarry Smith 
53a074057SBarry Smith /* This is to fix the tetrahedron orientation from TetGen */
6*a4a685f2SJacob Faibussowitsch static PetscErrorCode DMPlexInvertCells_CTetgen(PetscInt numCells, PetscInt numCorners, PetscInt cells[])
73a074057SBarry Smith {
83a074057SBarry Smith   PetscInt bound = numCells*numCorners, coff;
93a074057SBarry Smith 
103a074057SBarry Smith   PetscFunctionBegin;
11*a4a685f2SJacob Faibussowitsch #define SWAP(a,b) do { PetscInt tmp = (a); (a) = (b); (b) = tmp; } while(0)
1296ca5757SLisandro Dalcin   for (coff = 0; coff < bound; coff += numCorners) SWAP(cells[coff],cells[coff+1]);
1396ca5757SLisandro Dalcin #undef SWAP
143a074057SBarry Smith   PetscFunctionReturn(0);
153a074057SBarry Smith }
163a074057SBarry Smith 
173a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerate_CTetgen(DM boundary, PetscBool interpolate, DM *dm)
183a074057SBarry Smith {
193a074057SBarry Smith   MPI_Comm       comm;
203a074057SBarry Smith   const PetscInt dim       = 3;
213a074057SBarry Smith   const char    *labelName = "marker";
223a074057SBarry Smith   PLC           *in, *out;
233a074057SBarry Smith   DMLabel        label;
243a074057SBarry Smith   PetscInt       verbose = 0, vStart, vEnd, v, fStart, fEnd, f;
253a074057SBarry Smith   PetscMPIInt    rank;
263a074057SBarry Smith   PetscErrorCode ierr;
273a074057SBarry Smith 
283a074057SBarry Smith   PetscFunctionBegin;
293a074057SBarry Smith   ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr);
303a074057SBarry Smith   ierr = PetscOptionsGetInt(NULL,((PetscObject) boundary)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr);
313a074057SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
323a074057SBarry Smith   ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr);
333a074057SBarry Smith   ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr);
343a074057SBarry Smith   ierr = PLCCreate(&in);CHKERRQ(ierr);
353a074057SBarry Smith   ierr = PLCCreate(&out);CHKERRQ(ierr);
363a074057SBarry Smith 
373a074057SBarry Smith   in->numberofpoints = vEnd - vStart;
383a074057SBarry Smith   if (in->numberofpoints > 0) {
393a074057SBarry Smith     PetscSection coordSection;
403a074057SBarry Smith     Vec          coordinates;
413a074057SBarry Smith     PetscScalar *array;
423a074057SBarry Smith 
433a074057SBarry Smith     ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr);
443a074057SBarry Smith     ierr = PetscMalloc1(in->numberofpoints,       &in->pointmarkerlist);CHKERRQ(ierr);
453a074057SBarry Smith     ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr);
463a074057SBarry Smith     ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr);
473a074057SBarry Smith     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
483a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
493a074057SBarry Smith       const PetscInt idx = v - vStart;
503a074057SBarry Smith       PetscInt       off, d, m = -1;
513a074057SBarry Smith 
523a074057SBarry Smith       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
533a074057SBarry Smith       for (d = 0; d < dim; ++d) {
543a074057SBarry Smith         in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
553a074057SBarry Smith       }
563a074057SBarry Smith       if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);}
573a074057SBarry Smith 
583a074057SBarry Smith       in->pointmarkerlist[idx] = (int) m;
593a074057SBarry Smith     }
603a074057SBarry Smith     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
613a074057SBarry Smith   }
623a074057SBarry Smith   ierr  = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr);
633a074057SBarry Smith 
643a074057SBarry Smith   in->numberoffacets = fEnd - fStart;
653a074057SBarry Smith   if (in->numberoffacets > 0) {
663a074057SBarry Smith     ierr = PetscMalloc1(in->numberoffacets, &in->facetlist);CHKERRQ(ierr);
673a074057SBarry Smith     ierr = PetscMalloc1(in->numberoffacets,   &in->facetmarkerlist);CHKERRQ(ierr);
683a074057SBarry Smith     for (f = fStart; f < fEnd; ++f) {
693a074057SBarry Smith       const PetscInt idx     = f - fStart;
703a074057SBarry Smith       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, m = -1;
713a074057SBarry Smith       polygon       *poly;
723a074057SBarry Smith 
733a074057SBarry Smith       in->facetlist[idx].numberofpolygons = 1;
743a074057SBarry Smith 
753a074057SBarry Smith       ierr = PetscMalloc1(in->facetlist[idx].numberofpolygons, &in->facetlist[idx].polygonlist);CHKERRQ(ierr);
763a074057SBarry Smith 
773a074057SBarry Smith       in->facetlist[idx].numberofholes    = 0;
783a074057SBarry Smith       in->facetlist[idx].holelist         = NULL;
793a074057SBarry Smith 
803a074057SBarry Smith       ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
813a074057SBarry Smith       for (p = 0; p < numPoints*2; p += 2) {
823a074057SBarry Smith         const PetscInt point = points[p];
833a074057SBarry Smith         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
843a074057SBarry Smith       }
853a074057SBarry Smith 
863a074057SBarry Smith       poly                   = in->facetlist[idx].polygonlist;
873a074057SBarry Smith       poly->numberofvertices = numVertices;
883a074057SBarry Smith       ierr                   = PetscMalloc1(poly->numberofvertices, &poly->vertexlist);CHKERRQ(ierr);
893a074057SBarry Smith       for (v = 0; v < numVertices; ++v) {
903a074057SBarry Smith         const PetscInt vIdx = points[v] - vStart;
913a074057SBarry Smith         poly->vertexlist[v] = vIdx;
923a074057SBarry Smith       }
933a074057SBarry Smith       if (label) {ierr = DMLabelGetValue(label, f, &m);CHKERRQ(ierr);}
943a074057SBarry Smith       in->facetmarkerlist[idx] = (int) m;
953a074057SBarry Smith       ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
963a074057SBarry Smith     }
973a074057SBarry Smith   }
983a074057SBarry Smith   if (!rank) {
993a074057SBarry Smith     TetGenOpts t;
1003a074057SBarry Smith 
1013a074057SBarry Smith     ierr        = TetGenOptsInitialize(&t);CHKERRQ(ierr);
1023a074057SBarry Smith     t.in        = boundary; /* Should go away */
1033a074057SBarry Smith     t.plc       = 1;
1043a074057SBarry Smith     t.quality   = 1;
1053a074057SBarry Smith     t.edgesout  = 1;
1063a074057SBarry Smith     t.zeroindex = 1;
1073a074057SBarry Smith     t.quiet     = 1;
1083a074057SBarry Smith     t.verbose   = verbose;
1093a074057SBarry Smith     ierr        = TetGenCheckOpts(&t);CHKERRQ(ierr);
1103a074057SBarry Smith     ierr        = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr);
1113a074057SBarry Smith   }
1123a074057SBarry Smith   {
1133a074057SBarry Smith     DMLabel        glabel      = NULL;
1143a074057SBarry Smith     const PetscInt numCorners  = 4;
1153a074057SBarry Smith     const PetscInt numCells    = out->numberoftetrahedra;
1163a074057SBarry Smith     const PetscInt numVertices = out->numberofpoints;
117*a4a685f2SJacob Faibussowitsch     PetscReal      *meshCoords;
118*a4a685f2SJacob Faibussowitsch     PetscInt       *cells;
1193a074057SBarry Smith 
120*a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) == sizeof (out->pointlist[0])) {
121*a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *) out->pointlist;
122*a4a685f2SJacob Faibussowitsch     } else {
1233a074057SBarry Smith       PetscInt i;
1243a074057SBarry Smith 
125*a4a685f2SJacob Faibussowitsch       ierr = PetscMalloc1(dim * numVertices,&meshCoords);CHKERRQ(ierr);
126*a4a685f2SJacob Faibussowitsch       for (i = 0; i < dim * numVertices; i++) {
1273a074057SBarry Smith         meshCoords[i] = (PetscReal) out->pointlist[i];
1283a074057SBarry Smith       }
1293a074057SBarry Smith     }
130*a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) == sizeof (out->tetrahedronlist[0])) {
131*a4a685f2SJacob Faibussowitsch       cells = (PetscInt *) out->tetrahedronlist;
132*a4a685f2SJacob Faibussowitsch     } else {
133*a4a685f2SJacob Faibussowitsch       PetscInt i;
134*a4a685f2SJacob Faibussowitsch 
135*a4a685f2SJacob Faibussowitsch       ierr = PetscMalloc1(numCells * numCorners, &cells);CHKERRQ(ierr);
136*a4a685f2SJacob Faibussowitsch       for (i = 0; i < numCells * numCorners; i++) {
137*a4a685f2SJacob Faibussowitsch         cells[i] = (PetscInt) out->tetrahedronlist[i];
138*a4a685f2SJacob Faibussowitsch       }
139*a4a685f2SJacob Faibussowitsch     }
1403a074057SBarry Smith 
14196ca5757SLisandro Dalcin     ierr = DMPlexInvertCells_CTetgen(numCells, numCorners, cells);CHKERRQ(ierr);
142*a4a685f2SJacob Faibussowitsch     ierr = DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr);
143*a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) != sizeof (out->pointlist[0])) {
1443a074057SBarry Smith       ierr = PetscFree(meshCoords);CHKERRQ(ierr);
1453a074057SBarry Smith     }
146*a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) != sizeof (out->tetrahedronlist[0])) {
147*a4a685f2SJacob Faibussowitsch       ierr = PetscFree(cells);CHKERRQ(ierr);
148*a4a685f2SJacob Faibussowitsch     }
149*a4a685f2SJacob Faibussowitsch     if (label) {
150*a4a685f2SJacob Faibussowitsch       ierr = DMCreateLabel(*dm, labelName);CHKERRQ(ierr);
151*a4a685f2SJacob Faibussowitsch       ierr = DMGetLabel(*dm, labelName, &glabel);CHKERRQ(ierr);
152*a4a685f2SJacob Faibussowitsch     }
1533a074057SBarry Smith     /* Set labels */
1543a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
1553a074057SBarry Smith       if (out->pointmarkerlist[v]) {
1563a074057SBarry Smith         if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);}
1573a074057SBarry Smith       }
1583a074057SBarry Smith     }
1593a074057SBarry Smith     if (interpolate) {
1603a074057SBarry Smith       PetscInt e;
1613a074057SBarry Smith 
1623a074057SBarry Smith       for (e = 0; e < out->numberofedges; e++) {
1633a074057SBarry Smith         if (out->edgemarkerlist[e]) {
1643a074057SBarry Smith           const PetscInt  vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells};
1653a074057SBarry Smith           const PetscInt *edges;
1663a074057SBarry Smith           PetscInt        numEdges;
1673a074057SBarry Smith 
1683a074057SBarry Smith           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1693a074057SBarry Smith           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
1703a074057SBarry Smith           if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);}
1713a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1723a074057SBarry Smith         }
1733a074057SBarry Smith       }
1743a074057SBarry Smith       for (f = 0; f < out->numberoftrifaces; f++) {
1753a074057SBarry Smith         if (out->trifacemarkerlist[f]) {
1763a074057SBarry Smith           const PetscInt  vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells};
1773a074057SBarry Smith           const PetscInt *faces;
1783a074057SBarry Smith           PetscInt        numFaces;
1793a074057SBarry Smith 
1803a074057SBarry Smith           ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1813a074057SBarry Smith           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
1823a074057SBarry Smith           if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);}
1833a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1843a074057SBarry Smith         }
1853a074057SBarry Smith       }
1863a074057SBarry Smith     }
1873a074057SBarry Smith     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
1883a074057SBarry Smith   }
1893a074057SBarry Smith 
1903a074057SBarry Smith   ierr = PLCDestroy(&in);CHKERRQ(ierr);
1913a074057SBarry Smith   ierr = PLCDestroy(&out);CHKERRQ(ierr);
1923a074057SBarry Smith   PetscFunctionReturn(0);
1933a074057SBarry Smith }
1943a074057SBarry Smith 
1953a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined)
1963a074057SBarry Smith {
1973a074057SBarry Smith   MPI_Comm       comm;
1983a074057SBarry Smith   const PetscInt dim       = 3;
1993a074057SBarry Smith   const char    *labelName = "marker";
2003a074057SBarry Smith   PLC           *in, *out;
2013a074057SBarry Smith   DMLabel        label;
2023a074057SBarry Smith   PetscInt       verbose = 0, vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal;
2033a074057SBarry Smith   PetscMPIInt    rank;
2043a074057SBarry Smith   PetscErrorCode ierr;
2053a074057SBarry Smith 
2063a074057SBarry Smith   PetscFunctionBegin;
2073a074057SBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
2083a074057SBarry Smith   ierr = PetscOptionsGetInt(NULL,((PetscObject) dm)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr);
2093a074057SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2103a074057SBarry Smith   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2113a074057SBarry Smith   ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
2123a074057SBarry Smith   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
2133a074057SBarry Smith   ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr);
2143a074057SBarry Smith   ierr = PLCCreate(&in);CHKERRQ(ierr);
2153a074057SBarry Smith   ierr = PLCCreate(&out);CHKERRQ(ierr);
2163a074057SBarry Smith 
2173a074057SBarry Smith   in->numberofpoints = vEnd - vStart;
2183a074057SBarry Smith   if (in->numberofpoints > 0) {
2193a074057SBarry Smith     PetscSection coordSection;
2203a074057SBarry Smith     Vec          coordinates;
2213a074057SBarry Smith     PetscScalar *array;
2223a074057SBarry Smith 
2233a074057SBarry Smith     ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr);
2243a074057SBarry Smith     ierr = PetscMalloc1(in->numberofpoints, &in->pointmarkerlist);CHKERRQ(ierr);
2253a074057SBarry Smith     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2263a074057SBarry Smith     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2273a074057SBarry Smith     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
2283a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
2293a074057SBarry Smith       const PetscInt idx = v - vStart;
2303a074057SBarry Smith       PetscInt       off, d, m = -1;
2313a074057SBarry Smith 
2323a074057SBarry Smith       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
2333a074057SBarry Smith       for (d = 0; d < dim; ++d) {
2343a074057SBarry Smith         in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
2353a074057SBarry Smith       }
2363a074057SBarry Smith       if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);}
2373a074057SBarry Smith 
2383a074057SBarry Smith       in->pointmarkerlist[idx] = (int) m;
2393a074057SBarry Smith     }
2403a074057SBarry Smith     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
2413a074057SBarry Smith   }
2423a074057SBarry Smith   ierr  = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2433a074057SBarry Smith 
2443a074057SBarry Smith   in->numberofcorners       = 4;
2453a074057SBarry Smith   in->numberoftetrahedra    = cEnd - cStart;
2463a074057SBarry Smith   in->tetrahedronvolumelist = maxVolumes;
2473a074057SBarry Smith   if (in->numberoftetrahedra > 0) {
2483a074057SBarry Smith     ierr = PetscMalloc1(in->numberoftetrahedra*in->numberofcorners, &in->tetrahedronlist);CHKERRQ(ierr);
2493a074057SBarry Smith     for (c = cStart; c < cEnd; ++c) {
2503a074057SBarry Smith       const PetscInt idx      = c - cStart;
2513a074057SBarry Smith       PetscInt      *closure = NULL;
2523a074057SBarry Smith       PetscInt       closureSize;
2533a074057SBarry Smith 
2543a074057SBarry Smith       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
2553a074057SBarry 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);
2563a074057SBarry Smith       for (v = 0; v < 4; ++v) {
2573a074057SBarry Smith         in->tetrahedronlist[idx*in->numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
2583a074057SBarry Smith       }
2593a074057SBarry Smith       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
2603a074057SBarry Smith     }
2613a074057SBarry Smith   }
2623a074057SBarry Smith   if (!rank) {
2633a074057SBarry Smith     TetGenOpts t;
2643a074057SBarry Smith 
2653a074057SBarry Smith     ierr = TetGenOptsInitialize(&t);CHKERRQ(ierr);
2663a074057SBarry Smith 
2673a074057SBarry Smith     t.in        = dm; /* Should go away */
2683a074057SBarry Smith     t.refine    = 1;
2693a074057SBarry Smith     t.varvolume = 1;
2703a074057SBarry Smith     t.quality   = 1;
2713a074057SBarry Smith     t.edgesout  = 1;
2723a074057SBarry Smith     t.zeroindex = 1;
2733a074057SBarry Smith     t.quiet     = 1;
2743a074057SBarry Smith     t.verbose   = verbose; /* Change this */
2753a074057SBarry Smith 
2763a074057SBarry Smith     ierr = TetGenCheckOpts(&t);CHKERRQ(ierr);
2773a074057SBarry Smith     ierr = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr);
2783a074057SBarry Smith   }
2793a074057SBarry Smith   in->tetrahedronvolumelist = NULL;
2803a074057SBarry Smith   {
2813a074057SBarry Smith     DMLabel        rlabel      = NULL;
2823a074057SBarry Smith     const PetscInt numCorners  = 4;
2833a074057SBarry Smith     const PetscInt numCells    = out->numberoftetrahedra;
2843a074057SBarry Smith     const PetscInt numVertices = out->numberofpoints;
285*a4a685f2SJacob Faibussowitsch     PetscReal      *meshCoords;
286*a4a685f2SJacob Faibussowitsch     PetscInt       *cells;
2873a074057SBarry Smith     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
2883a074057SBarry Smith 
289*a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) == sizeof (out->pointlist[0])) {
290*a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *) out->pointlist;
291*a4a685f2SJacob Faibussowitsch     } else {
2923a074057SBarry Smith       PetscInt i;
2933a074057SBarry Smith 
294*a4a685f2SJacob Faibussowitsch       ierr = PetscMalloc1(dim * numVertices,&meshCoords);CHKERRQ(ierr);
295*a4a685f2SJacob Faibussowitsch       for (i = 0; i < dim * numVertices; i++) {
2963a074057SBarry Smith         meshCoords[i] = (PetscReal) out->pointlist[i];
2973a074057SBarry Smith       }
2983a074057SBarry Smith     }
299*a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) == sizeof (out->tetrahedronlist[0])) {
300*a4a685f2SJacob Faibussowitsch       cells = (PetscInt *) out->tetrahedronlist;
301*a4a685f2SJacob Faibussowitsch     } else {
302*a4a685f2SJacob Faibussowitsch       PetscInt i;
303*a4a685f2SJacob Faibussowitsch 
304*a4a685f2SJacob Faibussowitsch       ierr = PetscMalloc1(numCells * numCorners, &cells);CHKERRQ(ierr);
305*a4a685f2SJacob Faibussowitsch       for (i = 0; i < numCells * numCorners; i++) {
306*a4a685f2SJacob Faibussowitsch         cells[i] = (PetscInt) out->tetrahedronlist[i];
307*a4a685f2SJacob Faibussowitsch       }
308*a4a685f2SJacob Faibussowitsch     }
3093a074057SBarry Smith 
31096ca5757SLisandro Dalcin     ierr = DMPlexInvertCells_CTetgen(numCells, numCorners, cells);CHKERRQ(ierr);
311*a4a685f2SJacob Faibussowitsch     ierr = DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
312*a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) != sizeof (out->pointlist[0])) {
3133a074057SBarry Smith       ierr = PetscFree(meshCoords);CHKERRQ(ierr);
3143a074057SBarry Smith     }
315*a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) != sizeof (out->tetrahedronlist[0])) {
316*a4a685f2SJacob Faibussowitsch       ierr = PetscFree(cells);CHKERRQ(ierr);
317*a4a685f2SJacob Faibussowitsch     }
318*a4a685f2SJacob Faibussowitsch     if (label) {
319*a4a685f2SJacob Faibussowitsch       ierr = DMCreateLabel(*dmRefined, labelName);CHKERRQ(ierr);
320*a4a685f2SJacob Faibussowitsch       ierr = DMGetLabel(*dmRefined, labelName, &rlabel);CHKERRQ(ierr);
321*a4a685f2SJacob Faibussowitsch     }
3223a074057SBarry Smith     /* Set labels */
3233a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
3243a074057SBarry Smith       if (out->pointmarkerlist[v]) {
3253a074057SBarry Smith         if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);}
3263a074057SBarry Smith       }
3273a074057SBarry Smith     }
3283a074057SBarry Smith     if (interpolate) {
3293a074057SBarry Smith       PetscInt e, f;
3303a074057SBarry Smith 
3313a074057SBarry Smith       for (e = 0; e < out->numberofedges; e++) {
3323a074057SBarry Smith         if (out->edgemarkerlist[e]) {
3333a074057SBarry Smith           const PetscInt  vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells};
3343a074057SBarry Smith           const PetscInt *edges;
3353a074057SBarry Smith           PetscInt        numEdges;
3363a074057SBarry Smith 
3373a074057SBarry Smith           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
3383a074057SBarry Smith           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
3393a074057SBarry Smith           if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);}
3403a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
3413a074057SBarry Smith         }
3423a074057SBarry Smith       }
3433a074057SBarry Smith       for (f = 0; f < out->numberoftrifaces; f++) {
3443a074057SBarry Smith         if (out->trifacemarkerlist[f]) {
3453a074057SBarry Smith           const PetscInt  vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells};
3463a074057SBarry Smith           const PetscInt *faces;
3473a074057SBarry Smith           PetscInt        numFaces;
3483a074057SBarry Smith 
3493a074057SBarry Smith           ierr = DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
3503a074057SBarry Smith           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
3513a074057SBarry Smith           if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);}
3523a074057SBarry Smith           ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
3533a074057SBarry Smith         }
3543a074057SBarry Smith       }
3553a074057SBarry Smith     }
3563a074057SBarry Smith     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
3573a074057SBarry Smith   }
3583a074057SBarry Smith   ierr = PLCDestroy(&in);CHKERRQ(ierr);
3593a074057SBarry Smith   ierr = PLCDestroy(&out);CHKERRQ(ierr);
3603a074057SBarry Smith   PetscFunctionReturn(0);
3613a074057SBarry Smith }
362