xref: /petsc/src/dm/impls/plex/generators/triangle/trigenerate.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
13a074057SBarry Smith #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
23a074057SBarry Smith 
33316697fSBarry Smith #if !defined(ANSI_DECLARATORS)
43316697fSBarry Smith #define ANSI_DECLARATORS
53316697fSBarry Smith #endif
63a074057SBarry Smith #include <triangle.h>
73a074057SBarry Smith 
89371c9d4SSatish Balay static PetscErrorCode InitInput_Triangle(struct triangulateio *inputCtx) {
93a074057SBarry Smith   PetscFunctionBegin;
103a074057SBarry Smith   inputCtx->numberofpoints             = 0;
113a074057SBarry Smith   inputCtx->numberofpointattributes    = 0;
123a074057SBarry Smith   inputCtx->pointlist                  = NULL;
133a074057SBarry Smith   inputCtx->pointattributelist         = NULL;
143a074057SBarry Smith   inputCtx->pointmarkerlist            = NULL;
153a074057SBarry Smith   inputCtx->numberofsegments           = 0;
163a074057SBarry Smith   inputCtx->segmentlist                = NULL;
173a074057SBarry Smith   inputCtx->segmentmarkerlist          = NULL;
183a074057SBarry Smith   inputCtx->numberoftriangleattributes = 0;
193a074057SBarry Smith   inputCtx->trianglelist               = NULL;
203a074057SBarry Smith   inputCtx->numberofholes              = 0;
213a074057SBarry Smith   inputCtx->holelist                   = NULL;
223a074057SBarry Smith   inputCtx->numberofregions            = 0;
233a074057SBarry Smith   inputCtx->regionlist                 = NULL;
243a074057SBarry Smith   PetscFunctionReturn(0);
253a074057SBarry Smith }
263a074057SBarry Smith 
279371c9d4SSatish Balay static PetscErrorCode InitOutput_Triangle(struct triangulateio *outputCtx) {
283a074057SBarry Smith   PetscFunctionBegin;
293a074057SBarry Smith   outputCtx->numberofpoints        = 0;
303a074057SBarry Smith   outputCtx->pointlist             = NULL;
313a074057SBarry Smith   outputCtx->pointattributelist    = NULL;
323a074057SBarry Smith   outputCtx->pointmarkerlist       = NULL;
333a074057SBarry Smith   outputCtx->numberoftriangles     = 0;
343a074057SBarry Smith   outputCtx->trianglelist          = NULL;
353a074057SBarry Smith   outputCtx->triangleattributelist = NULL;
363a074057SBarry Smith   outputCtx->neighborlist          = NULL;
373a074057SBarry Smith   outputCtx->segmentlist           = NULL;
383a074057SBarry Smith   outputCtx->segmentmarkerlist     = NULL;
393a074057SBarry Smith   outputCtx->numberofedges         = 0;
403a074057SBarry Smith   outputCtx->edgelist              = NULL;
413a074057SBarry Smith   outputCtx->edgemarkerlist        = NULL;
423a074057SBarry Smith   PetscFunctionReturn(0);
433a074057SBarry Smith }
443a074057SBarry Smith 
459371c9d4SSatish Balay static PetscErrorCode FiniOutput_Triangle(struct triangulateio *outputCtx) {
463a074057SBarry Smith   PetscFunctionBegin;
473a074057SBarry Smith   free(outputCtx->pointlist);
483a074057SBarry Smith   free(outputCtx->pointmarkerlist);
493a074057SBarry Smith   free(outputCtx->segmentlist);
503a074057SBarry Smith   free(outputCtx->segmentmarkerlist);
513a074057SBarry Smith   free(outputCtx->edgelist);
523a074057SBarry Smith   free(outputCtx->edgemarkerlist);
533a074057SBarry Smith   free(outputCtx->trianglelist);
543a074057SBarry Smith   free(outputCtx->neighborlist);
553a074057SBarry Smith   PetscFunctionReturn(0);
563a074057SBarry Smith }
573a074057SBarry Smith 
589371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexGenerate_Triangle(DM boundary, PetscBool interpolate, DM *dm) {
593a074057SBarry Smith   MPI_Comm             comm;
603a074057SBarry Smith   DM_Plex             *mesh             = (DM_Plex *)boundary->data;
613a074057SBarry Smith   PetscInt             dim              = 2;
623a074057SBarry Smith   const PetscBool      createConvexHull = PETSC_FALSE;
633a074057SBarry Smith   const PetscBool      constrained      = PETSC_FALSE;
643a074057SBarry Smith   const char          *labelName        = "marker";
653a074057SBarry Smith   const char          *labelName2       = "Face Sets";
663a074057SBarry Smith   struct triangulateio in;
673a074057SBarry Smith   struct triangulateio out;
683a074057SBarry Smith   DMLabel              label, label2;
693a074057SBarry Smith   PetscInt             vStart, vEnd, v, eStart, eEnd, e;
703a074057SBarry Smith   PetscMPIInt          rank;
713a074057SBarry Smith 
723a074057SBarry Smith   PetscFunctionBegin;
739566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)boundary, &comm));
749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
759566063dSJacob Faibussowitsch   PetscCall(InitInput_Triangle(&in));
769566063dSJacob Faibussowitsch   PetscCall(InitOutput_Triangle(&out));
779566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd));
789566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(boundary, labelName, &label));
799566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(boundary, labelName2, &label2));
803a074057SBarry Smith 
813a074057SBarry Smith   in.numberofpoints = vEnd - vStart;
823a074057SBarry Smith   if (in.numberofpoints > 0) {
833a074057SBarry Smith     PetscSection coordSection;
843a074057SBarry Smith     Vec          coordinates;
853a074057SBarry Smith     PetscScalar *array;
863a074057SBarry Smith 
879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints * dim, &in.pointlist));
889566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints, &in.pointmarkerlist));
899566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(boundary, &coordinates));
909566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(boundary, &coordSection));
919566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &array));
923a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
933a074057SBarry Smith       const PetscInt idx = v - vStart;
94469e3fe5SMatthew G. Knepley       PetscInt       val, off, d;
953a074057SBarry Smith 
969566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
979371c9d4SSatish Balay       for (d = 0; d < dim; ++d) { in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]); }
98469e3fe5SMatthew G. Knepley       if (label) {
999566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, v, &val));
100469e3fe5SMatthew G. Knepley         in.pointmarkerlist[idx] = val;
101469e3fe5SMatthew G. Knepley       }
1023a074057SBarry Smith     }
1039566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &array));
1043a074057SBarry Smith   }
1059566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(boundary, 0, &eStart, &eEnd));
1063a074057SBarry Smith   in.numberofsegments = eEnd - eStart;
1073a074057SBarry Smith   if (in.numberofsegments > 0) {
1089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofsegments * 2, &in.segmentlist));
1099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofsegments, &in.segmentmarkerlist));
1103a074057SBarry Smith     for (e = eStart; e < eEnd; ++e) {
1113a074057SBarry Smith       const PetscInt  idx = e - eStart;
1123a074057SBarry Smith       const PetscInt *cone;
113469e3fe5SMatthew G. Knepley       PetscInt        val;
1143a074057SBarry Smith 
1159566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(boundary, e, &cone));
1163a074057SBarry Smith 
1173a074057SBarry Smith       in.segmentlist[idx * 2 + 0] = cone[0] - vStart;
1183a074057SBarry Smith       in.segmentlist[idx * 2 + 1] = cone[1] - vStart;
1193a074057SBarry Smith 
120469e3fe5SMatthew G. Knepley       if (label) {
1219566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, e, &val));
122469e3fe5SMatthew G. Knepley         in.segmentmarkerlist[idx] = val;
123469e3fe5SMatthew G. Knepley       }
1243a074057SBarry Smith     }
1253a074057SBarry Smith   }
1263a074057SBarry Smith #if 0 /* Do not currently support holes */
1273a074057SBarry Smith   PetscReal *holeCoords;
1283a074057SBarry Smith   PetscInt   h, d;
1293a074057SBarry Smith 
1309566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHoles(boundary, &in.numberofholes, &holeCords));
1313a074057SBarry Smith   if (in.numberofholes > 0) {
1329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofholes*dim, &in.holelist));
1333a074057SBarry Smith     for (h = 0; h < in.numberofholes; ++h) {
1343a074057SBarry Smith       for (d = 0; d < dim; ++d) {
1353a074057SBarry Smith         in.holelist[h*dim+d] = holeCoords[h*dim+d];
1363a074057SBarry Smith       }
1373a074057SBarry Smith     }
1383a074057SBarry Smith   }
1393a074057SBarry Smith #endif
140dd400576SPatrick Sanan   if (rank == 0) {
1413a074057SBarry Smith     char args[32];
1423a074057SBarry Smith 
1433a074057SBarry Smith     /* Take away 'Q' for verbose output */
1449566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(args, "pqezQ"));
1459566063dSJacob Faibussowitsch     if (createConvexHull) PetscCall(PetscStrcat(args, "c"));
1469566063dSJacob Faibussowitsch     if (constrained) PetscCall(PetscStrcpy(args, "zepDQ"));
1479371c9d4SSatish Balay     if (mesh->triangleOpts) {
1489371c9d4SSatish Balay       triangulate(mesh->triangleOpts, &in, &out, NULL);
1499371c9d4SSatish Balay     } else {
1509371c9d4SSatish Balay       triangulate(args, &in, &out, NULL);
1519371c9d4SSatish Balay     }
1523a074057SBarry Smith   }
1539566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointlist));
1549566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointmarkerlist));
1559566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentlist));
1569566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentmarkerlist));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.holelist));
1583a074057SBarry Smith 
1593a074057SBarry Smith   {
1603a074057SBarry Smith     DMLabel        glabel      = NULL;
1613a074057SBarry Smith     DMLabel        glabel2     = NULL;
1623a074057SBarry Smith     const PetscInt numCorners  = 3;
1633a074057SBarry Smith     const PetscInt numCells    = out.numberoftriangles;
1643a074057SBarry Smith     const PetscInt numVertices = out.numberofpoints;
165a4a685f2SJacob Faibussowitsch     PetscInt      *cells;
166a4a685f2SJacob Faibussowitsch     PetscReal     *meshCoords;
1673a074057SBarry Smith 
168a4a685f2SJacob Faibussowitsch     if (sizeof(PetscReal) == sizeof(out.pointlist[0])) {
169a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *)out.pointlist;
170a4a685f2SJacob Faibussowitsch     } else {
171a4a685f2SJacob Faibussowitsch       PetscInt i;
172a4a685f2SJacob Faibussowitsch 
1739566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dim * numVertices, &meshCoords));
1749371c9d4SSatish Balay       for (i = 0; i < dim * numVertices; i++) { meshCoords[i] = (PetscReal)out.pointlist[i]; }
175a4a685f2SJacob Faibussowitsch     }
176a4a685f2SJacob Faibussowitsch     if (sizeof(PetscInt) == sizeof(out.trianglelist[0])) {
177a4a685f2SJacob Faibussowitsch       cells = (PetscInt *)out.trianglelist;
178a4a685f2SJacob Faibussowitsch     } else {
179a4a685f2SJacob Faibussowitsch       PetscInt i;
180a4a685f2SJacob Faibussowitsch 
1819566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(numCells * numCorners, &cells));
1829371c9d4SSatish Balay       for (i = 0; i < numCells * numCorners; i++) { cells[i] = (PetscInt)out.trianglelist[i]; }
183a4a685f2SJacob Faibussowitsch     }
1849566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm));
185*48a46eb9SPierre Jolivet     if (sizeof(PetscReal) != sizeof(out.pointlist[0])) PetscCall(PetscFree(meshCoords));
186*48a46eb9SPierre Jolivet     if (sizeof(PetscInt) != sizeof(out.trianglelist[0])) PetscCall(PetscFree(cells));
187a4a685f2SJacob Faibussowitsch     if (label) {
1889566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(*dm, labelName));
1899566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dm, labelName, &glabel));
190a4a685f2SJacob Faibussowitsch     }
191a4a685f2SJacob Faibussowitsch     if (label2) {
1929566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(*dm, labelName2));
1939566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dm, labelName2, &glabel2));
194a4a685f2SJacob Faibussowitsch     }
1953a074057SBarry Smith     /* Set labels */
1963a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
1973a074057SBarry Smith       if (out.pointmarkerlist[v]) {
1989566063dSJacob Faibussowitsch         if (glabel) PetscCall(DMLabelSetValue(glabel, v + numCells, out.pointmarkerlist[v]));
1993a074057SBarry Smith       }
2003a074057SBarry Smith     }
2013a074057SBarry Smith     if (interpolate) {
2023a074057SBarry Smith       for (e = 0; e < out.numberofedges; e++) {
2033a074057SBarry Smith         if (out.edgemarkerlist[e]) {
2043a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells};
2053a074057SBarry Smith           const PetscInt *edges;
2063a074057SBarry Smith           PetscInt        numEdges;
2073a074057SBarry Smith 
2089566063dSJacob Faibussowitsch           PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges));
20963a3b9bcSJacob Faibussowitsch           PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
2109566063dSJacob Faibussowitsch           if (glabel) PetscCall(DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]));
2119566063dSJacob Faibussowitsch           if (glabel2) PetscCall(DMLabelSetValue(glabel2, edges[0], out.edgemarkerlist[e]));
2129566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges));
2133a074057SBarry Smith         }
2143a074057SBarry Smith       }
2153a074057SBarry Smith     }
2169566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
2173a074057SBarry Smith   }
2183a074057SBarry Smith #if 0 /* Do not currently support holes */
2199566063dSJacob Faibussowitsch   PetscCall(DMPlexCopyHoles(*dm, boundary));
2203a074057SBarry Smith #endif
2219566063dSJacob Faibussowitsch   PetscCall(FiniOutput_Triangle(&out));
2223a074057SBarry Smith   PetscFunctionReturn(0);
2233a074057SBarry Smith }
2243a074057SBarry Smith 
2259371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexRefine_Triangle(DM dm, PetscReal *inmaxVolumes, DM *dmRefined) {
2263a074057SBarry Smith   MPI_Comm             comm;
2273a074057SBarry Smith   PetscInt             dim       = 2;
2283a074057SBarry Smith   const char          *labelName = "marker";
2293a074057SBarry Smith   struct triangulateio in;
2303a074057SBarry Smith   struct triangulateio out;
2313a074057SBarry Smith   DMLabel              label;
2328457c761SMatthew G. Knepley   PetscInt             vStart, vEnd, v, gcStart, cStart, cEnd, c, depth, depthGlobal;
2333a074057SBarry Smith   PetscMPIInt          rank;
234800b850cSBarry Smith   double              *maxVolumes;
2353a074057SBarry Smith 
2363a074057SBarry Smith   PetscFunctionBegin;
2379566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2399566063dSJacob Faibussowitsch   PetscCall(InitInput_Triangle(&in));
2409566063dSJacob Faibussowitsch   PetscCall(InitOutput_Triangle(&out));
2419566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
2421c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm));
2439566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2449566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, labelName, &label));
2453a074057SBarry Smith 
2463a074057SBarry Smith   in.numberofpoints = vEnd - vStart;
2473a074057SBarry Smith   if (in.numberofpoints > 0) {
2483a074057SBarry Smith     PetscSection coordSection;
2493a074057SBarry Smith     Vec          coordinates;
2503a074057SBarry Smith     PetscScalar *array;
2513a074057SBarry Smith 
2529566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints * dim, &in.pointlist));
2539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints, &in.pointmarkerlist));
2549566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2559566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dm, &coordSection));
2569566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &array));
2573a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
2583a074057SBarry Smith       const PetscInt idx = v - vStart;
259469e3fe5SMatthew G. Knepley       PetscInt       off, d, val;
2603a074057SBarry Smith 
2619566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
2629371c9d4SSatish Balay       for (d = 0; d < dim; ++d) { in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]); }
263469e3fe5SMatthew G. Knepley       if (label) {
2649566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, v, &val));
265469e3fe5SMatthew G. Knepley         in.pointmarkerlist[idx] = val;
266469e3fe5SMatthew G. Knepley       }
2673a074057SBarry Smith     }
2689566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &array));
2693a074057SBarry Smith   }
2709566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2719566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGhostCellStratum(dm, &gcStart, NULL));
2728457c761SMatthew G. Knepley   if (gcStart >= 0) cEnd = gcStart;
2733a074057SBarry Smith 
2743a074057SBarry Smith   in.numberofcorners   = 3;
2753a074057SBarry Smith   in.numberoftriangles = cEnd - cStart;
2763a074057SBarry Smith 
277800b850cSBarry Smith #if !defined(PETSC_USE_REAL_DOUBLE)
2789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cEnd - cStart, &maxVolumes));
279800b850cSBarry Smith   for (c = 0; c < cEnd - cStart; ++c) maxVolumes[c] = (double)inmaxVolumes[c];
280800b850cSBarry Smith #else
281800b850cSBarry Smith   maxVolumes = inmaxVolumes;
282800b850cSBarry Smith #endif
283800b850cSBarry Smith 
2843a074057SBarry Smith   in.trianglearealist = (double *)maxVolumes;
2853a074057SBarry Smith   if (in.numberoftriangles > 0) {
2869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberoftriangles * in.numberofcorners, &in.trianglelist));
2873a074057SBarry Smith     for (c = cStart; c < cEnd; ++c) {
2883a074057SBarry Smith       const PetscInt idx     = c - cStart;
2893a074057SBarry Smith       PetscInt      *closure = NULL;
2903a074057SBarry Smith       PetscInt       closureSize;
2913a074057SBarry Smith 
2929566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
29363a3b9bcSJacob Faibussowitsch       PetscCheck(!(closureSize != 4) || !(closureSize != 7), comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a triangle, %" PetscInt_FMT " vertices in closure", closureSize);
2949371c9d4SSatish Balay       for (v = 0; v < 3; ++v) { in.trianglelist[idx * in.numberofcorners + v] = closure[(v + closureSize - 3) * 2] - vStart; }
2959566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
2963a074057SBarry Smith     }
2973a074057SBarry Smith   }
2983a074057SBarry Smith   /* TODO: Segment markers are missing on input */
2993a074057SBarry Smith #if 0 /* Do not currently support holes */
3003a074057SBarry Smith   PetscReal *holeCoords;
3013a074057SBarry Smith   PetscInt   h, d;
3023a074057SBarry Smith 
3039566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHoles(boundary, &in.numberofholes, &holeCords));
3043a074057SBarry Smith   if (in.numberofholes > 0) {
3059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofholes*dim, &in.holelist));
3063a074057SBarry Smith     for (h = 0; h < in.numberofholes; ++h) {
3073a074057SBarry Smith       for (d = 0; d < dim; ++d) {
3083a074057SBarry Smith         in.holelist[h*dim+d] = holeCoords[h*dim+d];
3093a074057SBarry Smith       }
3103a074057SBarry Smith     }
3113a074057SBarry Smith   }
3123a074057SBarry Smith #endif
313dd400576SPatrick Sanan   if (rank == 0) {
3143a074057SBarry Smith     char args[32];
3153a074057SBarry Smith 
3163a074057SBarry Smith     /* Take away 'Q' for verbose output */
3179566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(args, "pqezQra"));
3183a074057SBarry Smith     triangulate(args, &in, &out, NULL);
3193a074057SBarry Smith   }
3209566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointlist));
3219566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointmarkerlist));
3229566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentlist));
3239566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentmarkerlist));
3249566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.trianglelist));
3253a074057SBarry Smith 
3263a074057SBarry Smith   {
3273a074057SBarry Smith     DMLabel        rlabel      = NULL;
3283a074057SBarry Smith     const PetscInt numCorners  = 3;
3293a074057SBarry Smith     const PetscInt numCells    = out.numberoftriangles;
3303a074057SBarry Smith     const PetscInt numVertices = out.numberofpoints;
331a4a685f2SJacob Faibussowitsch     PetscInt      *cells;
332a4a685f2SJacob Faibussowitsch     PetscReal     *meshCoords;
3333a074057SBarry Smith     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
3343a074057SBarry Smith 
335a4a685f2SJacob Faibussowitsch     if (sizeof(PetscReal) == sizeof(out.pointlist[0])) {
336a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *)out.pointlist;
337a4a685f2SJacob Faibussowitsch     } else {
338a4a685f2SJacob Faibussowitsch       PetscInt i;
339a4a685f2SJacob Faibussowitsch 
3409566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dim * numVertices, &meshCoords));
3419371c9d4SSatish Balay       for (i = 0; i < dim * numVertices; i++) { meshCoords[i] = (PetscReal)out.pointlist[i]; }
342a4a685f2SJacob Faibussowitsch     }
343a4a685f2SJacob Faibussowitsch     if (sizeof(PetscInt) == sizeof(out.trianglelist[0])) {
344a4a685f2SJacob Faibussowitsch       cells = (PetscInt *)out.trianglelist;
345a4a685f2SJacob Faibussowitsch     } else {
346a4a685f2SJacob Faibussowitsch       PetscInt i;
347a4a685f2SJacob Faibussowitsch 
3489566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(numCells * numCorners, &cells));
3499371c9d4SSatish Balay       for (i = 0; i < numCells * numCorners; i++) { cells[i] = (PetscInt)out.trianglelist[i]; }
350a4a685f2SJacob Faibussowitsch     }
351a4a685f2SJacob Faibussowitsch 
3529566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined));
353a4a685f2SJacob Faibussowitsch     if (label) {
3549566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(*dmRefined, labelName));
3559566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dmRefined, labelName, &rlabel));
356a4a685f2SJacob Faibussowitsch     }
357*48a46eb9SPierre Jolivet     if (sizeof(PetscReal) != sizeof(out.pointlist[0])) PetscCall(PetscFree(meshCoords));
358*48a46eb9SPierre Jolivet     if (sizeof(PetscInt) != sizeof(out.trianglelist[0])) PetscCall(PetscFree(cells));
3593a074057SBarry Smith     /* Set labels */
3603a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
3613a074057SBarry Smith       if (out.pointmarkerlist[v]) {
3629566063dSJacob Faibussowitsch         if (rlabel) PetscCall(DMLabelSetValue(rlabel, v + numCells, out.pointmarkerlist[v]));
3633a074057SBarry Smith       }
3643a074057SBarry Smith     }
3653a074057SBarry Smith     if (interpolate) {
3663a074057SBarry Smith       PetscInt e;
3673a074057SBarry Smith 
3683a074057SBarry Smith       for (e = 0; e < out.numberofedges; e++) {
3693a074057SBarry Smith         if (out.edgemarkerlist[e]) {
3703a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells};
3713a074057SBarry Smith           const PetscInt *edges;
3723a074057SBarry Smith           PetscInt        numEdges;
3733a074057SBarry Smith 
3749566063dSJacob Faibussowitsch           PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges));
37563a3b9bcSJacob Faibussowitsch           PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
3769566063dSJacob Faibussowitsch           if (rlabel) PetscCall(DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]));
3779566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges));
3783a074057SBarry Smith         }
3793a074057SBarry Smith       }
3803a074057SBarry Smith     }
3819566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE));
3823a074057SBarry Smith   }
3833a074057SBarry Smith #if 0 /* Do not currently support holes */
3849566063dSJacob Faibussowitsch   PetscCall(DMPlexCopyHoles(*dm, boundary));
3853a074057SBarry Smith #endif
3869566063dSJacob Faibussowitsch   PetscCall(FiniOutput_Triangle(&out));
387800b850cSBarry Smith #if !defined(PETSC_USE_REAL_DOUBLE)
3889566063dSJacob Faibussowitsch   PetscCall(PetscFree(maxVolumes));
389800b850cSBarry Smith #endif
3903a074057SBarry Smith   PetscFunctionReturn(0);
3913a074057SBarry Smith }
392