xref: /petsc/src/dm/impls/plex/generators/triangle/trigenerate.c (revision 84e0bd68ca5b84e7bd4f8b5b1bb8b0c638b469dc)
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 
8d71ae5a4SJacob Faibussowitsch static PetscErrorCode InitInput_Triangle(struct triangulateio *inputCtx)
9d71ae5a4SJacob Faibussowitsch {
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;
253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
263a074057SBarry Smith }
273a074057SBarry Smith 
28d71ae5a4SJacob Faibussowitsch static PetscErrorCode InitOutput_Triangle(struct triangulateio *outputCtx)
29d71ae5a4SJacob Faibussowitsch {
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;
443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
453a074057SBarry Smith }
463a074057SBarry Smith 
47d71ae5a4SJacob Faibussowitsch static PetscErrorCode FiniOutput_Triangle(struct triangulateio *outputCtx)
48d71ae5a4SJacob Faibussowitsch {
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);
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
593a074057SBarry Smith }
603a074057SBarry Smith 
61d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexGenerate_Triangle(DM boundary, PetscBool interpolate, DM *dm)
62d71ae5a4SJacob Faibussowitsch {
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;
75*84e0bd68SStefano Zampini   PetscBool            flg;
76*84e0bd68SStefano Zampini   char                 opts[64];
773a074057SBarry Smith 
783a074057SBarry Smith   PetscFunctionBegin;
799566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)boundary, &comm));
809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
819566063dSJacob Faibussowitsch   PetscCall(InitInput_Triangle(&in));
829566063dSJacob Faibussowitsch   PetscCall(InitOutput_Triangle(&out));
839566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd));
849566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(boundary, labelName, &label));
859566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(boundary, labelName2, &label2));
86*84e0bd68SStefano Zampini   PetscCall(PetscOptionsGetString(((PetscObject)boundary)->options, ((PetscObject)boundary)->prefix, "-dm_plex_generate_triangle_opts", opts, sizeof(opts), &flg));
87*84e0bd68SStefano Zampini   if (flg) PetscCall(DMPlexTriangleSetOptions(boundary, opts));
883a074057SBarry Smith 
89784a7648SPierre Jolivet   PetscCall(PetscCIntCast(vEnd - vStart, &in.numberofpoints));
903a074057SBarry Smith   if (in.numberofpoints > 0) {
913a074057SBarry Smith     PetscSection coordSection;
923a074057SBarry Smith     Vec          coordinates;
933a074057SBarry Smith     PetscScalar *array;
943a074057SBarry Smith 
959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints * dim, &in.pointlist));
969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints, &in.pointmarkerlist));
979566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(boundary, &coordinates));
989566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(boundary, &coordSection));
999566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &array));
1003a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
1013a074057SBarry Smith       const PetscInt idx = v - vStart;
102469e3fe5SMatthew G. Knepley       PetscInt       val, off, d;
1033a074057SBarry Smith 
1049566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
105ad540459SPierre Jolivet       for (d = 0; d < dim; ++d) in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]);
106469e3fe5SMatthew G. Knepley       if (label) {
1079566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, v, &val));
108784a7648SPierre Jolivet         PetscCall(PetscCIntCast(val, &in.pointmarkerlist[idx]));
109469e3fe5SMatthew G. Knepley       }
1103a074057SBarry Smith     }
1119566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &array));
1123a074057SBarry Smith   }
1139566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(boundary, 0, &eStart, &eEnd));
114784a7648SPierre Jolivet   PetscCall(PetscCIntCast(eEnd - eStart, &in.numberofsegments));
1153a074057SBarry Smith   if (in.numberofsegments > 0) {
1169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofsegments * 2, &in.segmentlist));
1179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofsegments, &in.segmentmarkerlist));
1183a074057SBarry Smith     for (e = eStart; e < eEnd; ++e) {
1193a074057SBarry Smith       const PetscInt  idx = e - eStart;
1203a074057SBarry Smith       const PetscInt *cone;
121469e3fe5SMatthew G. Knepley       PetscInt        val;
1223a074057SBarry Smith 
1239566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(boundary, e, &cone));
1243a074057SBarry Smith 
125784a7648SPierre Jolivet       PetscCall(PetscCIntCast(cone[0] - vStart, &in.segmentlist[idx * 2 + 0]));
126784a7648SPierre Jolivet       PetscCall(PetscCIntCast(cone[1] - vStart, &in.segmentlist[idx * 2 + 1]));
1273a074057SBarry Smith 
128469e3fe5SMatthew G. Knepley       if (label) {
1299566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, e, &val));
130784a7648SPierre Jolivet         PetscCall(PetscCIntCast(val, &in.segmentmarkerlist[idx]));
131469e3fe5SMatthew G. Knepley       }
1323a074057SBarry Smith     }
1333a074057SBarry Smith   }
1343a074057SBarry Smith #if 0 /* Do not currently support holes */
1353a074057SBarry Smith   PetscReal *holeCoords;
1363a074057SBarry Smith   PetscInt   h, d;
1373a074057SBarry Smith 
1389566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHoles(boundary, &in.numberofholes, &holeCords));
1393a074057SBarry Smith   if (in.numberofholes > 0) {
1409566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofholes*dim, &in.holelist));
1413a074057SBarry Smith     for (h = 0; h < in.numberofholes; ++h) {
1423a074057SBarry Smith       for (d = 0; d < dim; ++d) {
1433a074057SBarry Smith         in.holelist[h*dim+d] = holeCoords[h*dim+d];
1443a074057SBarry Smith       }
1453a074057SBarry Smith     }
1463a074057SBarry Smith   }
1473a074057SBarry Smith #endif
148dd400576SPatrick Sanan   if (rank == 0) {
1493a074057SBarry Smith     char args[32];
1503a074057SBarry Smith 
1513a074057SBarry Smith     /* Take away 'Q' for verbose output */
152c6a7a370SJeremy L Thompson     PetscCall(PetscStrncpy(args, "pqezQ", sizeof(args)));
153c6a7a370SJeremy L Thompson     if (createConvexHull) PetscCall(PetscStrlcat(args, "c", sizeof(args)));
154c6a7a370SJeremy L Thompson     if (constrained) PetscCall(PetscStrncpy(args, "zepDQ", sizeof(args)));
1559371c9d4SSatish Balay     if (mesh->triangleOpts) {
1569371c9d4SSatish Balay       triangulate(mesh->triangleOpts, &in, &out, NULL);
1579371c9d4SSatish Balay     } else {
1589371c9d4SSatish Balay       triangulate(args, &in, &out, NULL);
1599371c9d4SSatish Balay     }
1603a074057SBarry Smith   }
1619566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointlist));
1629566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointmarkerlist));
1639566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentlist));
1649566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentmarkerlist));
1659566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.holelist));
1663a074057SBarry Smith 
1673a074057SBarry Smith   {
1683a074057SBarry Smith     DMLabel        glabel      = NULL;
1693a074057SBarry Smith     DMLabel        glabel2     = NULL;
1703a074057SBarry Smith     const PetscInt numCorners  = 3;
1713a074057SBarry Smith     const PetscInt numCells    = out.numberoftriangles;
1723a074057SBarry Smith     const PetscInt numVertices = out.numberofpoints;
173a4a685f2SJacob Faibussowitsch     PetscInt      *cells;
174a4a685f2SJacob Faibussowitsch     PetscReal     *meshCoords;
1753a074057SBarry Smith 
176a4a685f2SJacob Faibussowitsch     if (sizeof(PetscReal) == sizeof(out.pointlist[0])) {
177a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *)out.pointlist;
178a4a685f2SJacob Faibussowitsch     } else {
179a4a685f2SJacob Faibussowitsch       PetscInt i;
180a4a685f2SJacob Faibussowitsch 
1819566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dim * numVertices, &meshCoords));
182ad540459SPierre Jolivet       for (i = 0; i < dim * numVertices; i++) meshCoords[i] = (PetscReal)out.pointlist[i];
183a4a685f2SJacob Faibussowitsch     }
184a4a685f2SJacob Faibussowitsch     if (sizeof(PetscInt) == sizeof(out.trianglelist[0])) {
185a4a685f2SJacob Faibussowitsch       cells = (PetscInt *)out.trianglelist;
186a4a685f2SJacob Faibussowitsch     } else {
187a4a685f2SJacob Faibussowitsch       PetscInt i;
188a4a685f2SJacob Faibussowitsch 
1899566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(numCells * numCorners, &cells));
190ad540459SPierre Jolivet       for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out.trianglelist[i];
191a4a685f2SJacob Faibussowitsch     }
1929566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm));
19348a46eb9SPierre Jolivet     if (sizeof(PetscReal) != sizeof(out.pointlist[0])) PetscCall(PetscFree(meshCoords));
19448a46eb9SPierre Jolivet     if (sizeof(PetscInt) != sizeof(out.trianglelist[0])) PetscCall(PetscFree(cells));
195a4a685f2SJacob Faibussowitsch     if (label) {
1969566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(*dm, labelName));
1979566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dm, labelName, &glabel));
198a4a685f2SJacob Faibussowitsch     }
199a4a685f2SJacob Faibussowitsch     if (label2) {
2009566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(*dm, labelName2));
2019566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dm, labelName2, &glabel2));
202a4a685f2SJacob Faibussowitsch     }
2033a074057SBarry Smith     /* Set labels */
2043a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
2053a074057SBarry Smith       if (out.pointmarkerlist[v]) {
2069566063dSJacob Faibussowitsch         if (glabel) PetscCall(DMLabelSetValue(glabel, v + numCells, out.pointmarkerlist[v]));
2073a074057SBarry Smith       }
2083a074057SBarry Smith     }
2093a074057SBarry Smith     if (interpolate) {
2103a074057SBarry Smith       for (e = 0; e < out.numberofedges; e++) {
2113a074057SBarry Smith         if (out.edgemarkerlist[e]) {
2123a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells};
2133a074057SBarry Smith           const PetscInt *edges;
2143a074057SBarry Smith           PetscInt        numEdges;
2153a074057SBarry Smith 
2169566063dSJacob Faibussowitsch           PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges));
21763a3b9bcSJacob Faibussowitsch           PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
2189566063dSJacob Faibussowitsch           if (glabel) PetscCall(DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]));
2199566063dSJacob Faibussowitsch           if (glabel2) PetscCall(DMLabelSetValue(glabel2, edges[0], out.edgemarkerlist[e]));
2209566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges));
2213a074057SBarry Smith         }
2223a074057SBarry Smith       }
2233a074057SBarry Smith     }
2249566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
2253a074057SBarry Smith   }
2263a074057SBarry Smith #if 0 /* Do not currently support holes */
2279566063dSJacob Faibussowitsch   PetscCall(DMPlexCopyHoles(*dm, boundary));
2283a074057SBarry Smith #endif
2299566063dSJacob Faibussowitsch   PetscCall(FiniOutput_Triangle(&out));
2303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2313a074057SBarry Smith }
2323a074057SBarry Smith 
233d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexRefine_Triangle(DM dm, PetscReal *inmaxVolumes, DM *dmRefined)
234d71ae5a4SJacob Faibussowitsch {
2353a074057SBarry Smith   MPI_Comm             comm;
2363a074057SBarry Smith   PetscInt             dim       = 2;
2373a074057SBarry Smith   const char          *labelName = "marker";
2383a074057SBarry Smith   struct triangulateio in;
2393a074057SBarry Smith   struct triangulateio out;
2403a074057SBarry Smith   DMLabel              label;
2418457c761SMatthew G. Knepley   PetscInt             vStart, vEnd, v, gcStart, cStart, cEnd, c, depth, depthGlobal;
2423a074057SBarry Smith   PetscMPIInt          rank;
243800b850cSBarry Smith   double              *maxVolumes;
2443a074057SBarry Smith 
2453a074057SBarry Smith   PetscFunctionBegin;
2469566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2489566063dSJacob Faibussowitsch   PetscCall(InitInput_Triangle(&in));
2499566063dSJacob Faibussowitsch   PetscCall(InitOutput_Triangle(&out));
2509566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
251462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm));
2529566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2539566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, labelName, &label));
2543a074057SBarry Smith 
255784a7648SPierre Jolivet   PetscCall(PetscCIntCast(vEnd - vStart, &in.numberofpoints));
2563a074057SBarry Smith   if (in.numberofpoints > 0) {
2573a074057SBarry Smith     PetscSection coordSection;
2583a074057SBarry Smith     Vec          coordinates;
2593a074057SBarry Smith     PetscScalar *array;
2603a074057SBarry Smith 
2619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints * dim, &in.pointlist));
2629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints, &in.pointmarkerlist));
2639566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2649566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dm, &coordSection));
2659566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &array));
2663a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
2673a074057SBarry Smith       const PetscInt idx = v - vStart;
268469e3fe5SMatthew G. Knepley       PetscInt       off, d, val;
2693a074057SBarry Smith 
2709566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
271ad540459SPierre Jolivet       for (d = 0; d < dim; ++d) in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]);
272469e3fe5SMatthew G. Knepley       if (label) {
2739566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, v, &val));
274784a7648SPierre Jolivet         PetscCall(PetscCIntCast(val, &in.pointmarkerlist[idx]));
275469e3fe5SMatthew G. Knepley       }
2763a074057SBarry Smith     }
2779566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &array));
2783a074057SBarry Smith   }
2799566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2802827ebadSStefano Zampini   PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &gcStart, NULL));
2818457c761SMatthew G. Knepley   if (gcStart >= 0) cEnd = gcStart;
2823a074057SBarry Smith 
2833a074057SBarry Smith   in.numberofcorners = 3;
284784a7648SPierre Jolivet   PetscCall(PetscCIntCast(cEnd - cStart, &in.numberoftriangles));
2853a074057SBarry Smith 
286800b850cSBarry Smith #if !defined(PETSC_USE_REAL_DOUBLE)
2879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cEnd - cStart, &maxVolumes));
288800b850cSBarry Smith   for (c = 0; c < cEnd - cStart; ++c) maxVolumes[c] = (double)inmaxVolumes[c];
289800b850cSBarry Smith #else
290800b850cSBarry Smith   maxVolumes = inmaxVolumes;
291800b850cSBarry Smith #endif
292800b850cSBarry Smith 
2933a074057SBarry Smith   in.trianglearealist = (double *)maxVolumes;
2943a074057SBarry Smith   if (in.numberoftriangles > 0) {
2959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberoftriangles * in.numberofcorners, &in.trianglelist));
2963a074057SBarry Smith     for (c = cStart; c < cEnd; ++c) {
2973a074057SBarry Smith       const PetscInt idx     = c - cStart;
2983a074057SBarry Smith       PetscInt      *closure = NULL;
2993a074057SBarry Smith       PetscInt       closureSize;
3003a074057SBarry Smith 
3019566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
30263a3b9bcSJacob 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);
303784a7648SPierre Jolivet       for (v = 0; v < 3; ++v) PetscCall(PetscCIntCast(closure[(v + closureSize - 3) * 2] - vStart, &in.trianglelist[idx * in.numberofcorners + v]));
3049566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
3053a074057SBarry Smith     }
3063a074057SBarry Smith   }
3073a074057SBarry Smith   /* TODO: Segment markers are missing on input */
3083a074057SBarry Smith #if 0 /* Do not currently support holes */
3093a074057SBarry Smith   PetscReal *holeCoords;
3103a074057SBarry Smith   PetscInt   h, d;
3113a074057SBarry Smith 
3129566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHoles(boundary, &in.numberofholes, &holeCords));
3133a074057SBarry Smith   if (in.numberofholes > 0) {
3149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofholes*dim, &in.holelist));
3153a074057SBarry Smith     for (h = 0; h < in.numberofholes; ++h) {
3163a074057SBarry Smith       for (d = 0; d < dim; ++d) {
3173a074057SBarry Smith         in.holelist[h*dim+d] = holeCoords[h*dim+d];
3183a074057SBarry Smith       }
3193a074057SBarry Smith     }
3203a074057SBarry Smith   }
3213a074057SBarry Smith #endif
322dd400576SPatrick Sanan   if (rank == 0) {
3233a074057SBarry Smith     char args[32];
3243a074057SBarry Smith 
3253a074057SBarry Smith     /* Take away 'Q' for verbose output */
326c6a7a370SJeremy L Thompson     PetscCall(PetscStrncpy(args, "pqezQra", sizeof(args)));
3273a074057SBarry Smith     triangulate(args, &in, &out, NULL);
3283a074057SBarry Smith   }
3299566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointlist));
3309566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointmarkerlist));
3319566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentlist));
3329566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentmarkerlist));
3339566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.trianglelist));
3343a074057SBarry Smith 
3353a074057SBarry Smith   {
3363a074057SBarry Smith     DMLabel        rlabel      = NULL;
3373a074057SBarry Smith     const PetscInt numCorners  = 3;
3383a074057SBarry Smith     const PetscInt numCells    = out.numberoftriangles;
3393a074057SBarry Smith     const PetscInt numVertices = out.numberofpoints;
340a4a685f2SJacob Faibussowitsch     PetscInt      *cells;
341a4a685f2SJacob Faibussowitsch     PetscReal     *meshCoords;
3423a074057SBarry Smith     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
3433a074057SBarry Smith 
344a4a685f2SJacob Faibussowitsch     if (sizeof(PetscReal) == sizeof(out.pointlist[0])) {
345a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *)out.pointlist;
346a4a685f2SJacob Faibussowitsch     } else {
347a4a685f2SJacob Faibussowitsch       PetscInt i;
348a4a685f2SJacob Faibussowitsch 
3499566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dim * numVertices, &meshCoords));
350ad540459SPierre Jolivet       for (i = 0; i < dim * numVertices; i++) meshCoords[i] = (PetscReal)out.pointlist[i];
351a4a685f2SJacob Faibussowitsch     }
352a4a685f2SJacob Faibussowitsch     if (sizeof(PetscInt) == sizeof(out.trianglelist[0])) {
353a4a685f2SJacob Faibussowitsch       cells = (PetscInt *)out.trianglelist;
354a4a685f2SJacob Faibussowitsch     } else {
355a4a685f2SJacob Faibussowitsch       PetscInt i;
356a4a685f2SJacob Faibussowitsch 
3579566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(numCells * numCorners, &cells));
358ad540459SPierre Jolivet       for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out.trianglelist[i];
359a4a685f2SJacob Faibussowitsch     }
360a4a685f2SJacob Faibussowitsch 
3619566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined));
362a4a685f2SJacob Faibussowitsch     if (label) {
3639566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(*dmRefined, labelName));
3649566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dmRefined, labelName, &rlabel));
365a4a685f2SJacob Faibussowitsch     }
36648a46eb9SPierre Jolivet     if (sizeof(PetscReal) != sizeof(out.pointlist[0])) PetscCall(PetscFree(meshCoords));
36748a46eb9SPierre Jolivet     if (sizeof(PetscInt) != sizeof(out.trianglelist[0])) PetscCall(PetscFree(cells));
3683a074057SBarry Smith     /* Set labels */
3693a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
3703a074057SBarry Smith       if (out.pointmarkerlist[v]) {
3719566063dSJacob Faibussowitsch         if (rlabel) PetscCall(DMLabelSetValue(rlabel, v + numCells, out.pointmarkerlist[v]));
3723a074057SBarry Smith       }
3733a074057SBarry Smith     }
3743a074057SBarry Smith     if (interpolate) {
3753a074057SBarry Smith       PetscInt e;
3763a074057SBarry Smith 
3773a074057SBarry Smith       for (e = 0; e < out.numberofedges; e++) {
3783a074057SBarry Smith         if (out.edgemarkerlist[e]) {
3793a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells};
3803a074057SBarry Smith           const PetscInt *edges;
3813a074057SBarry Smith           PetscInt        numEdges;
3823a074057SBarry Smith 
3839566063dSJacob Faibussowitsch           PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges));
38463a3b9bcSJacob Faibussowitsch           PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
3859566063dSJacob Faibussowitsch           if (rlabel) PetscCall(DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]));
3869566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges));
3873a074057SBarry Smith         }
3883a074057SBarry Smith       }
3893a074057SBarry Smith     }
3909566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE));
3913a074057SBarry Smith   }
3923a074057SBarry Smith #if 0 /* Do not currently support holes */
3939566063dSJacob Faibussowitsch   PetscCall(DMPlexCopyHoles(*dm, boundary));
3943a074057SBarry Smith #endif
3959566063dSJacob Faibussowitsch   PetscCall(FiniOutput_Triangle(&out));
396800b850cSBarry Smith #if !defined(PETSC_USE_REAL_DOUBLE)
3979566063dSJacob Faibussowitsch   PetscCall(PetscFree(maxVolumes));
398800b850cSBarry Smith #endif
3993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4003a074057SBarry Smith }
401