xref: /petsc/src/dm/impls/plex/generators/triangle/trigenerate.c (revision 1c2dc1cbabe5212f80cf309c4ca5a26f0cadc660)
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 
83a074057SBarry Smith static PetscErrorCode InitInput_Triangle(struct triangulateio *inputCtx)
93a074057SBarry Smith {
103a074057SBarry Smith   PetscFunctionBegin;
113a074057SBarry Smith   inputCtx->numberofpoints             = 0;
123a074057SBarry Smith   inputCtx->numberofpointattributes    = 0;
133a074057SBarry Smith   inputCtx->pointlist                  = NULL;
143a074057SBarry Smith   inputCtx->pointattributelist         = NULL;
153a074057SBarry Smith   inputCtx->pointmarkerlist            = NULL;
163a074057SBarry Smith   inputCtx->numberofsegments           = 0;
173a074057SBarry Smith   inputCtx->segmentlist                = NULL;
183a074057SBarry Smith   inputCtx->segmentmarkerlist          = NULL;
193a074057SBarry Smith   inputCtx->numberoftriangleattributes = 0;
203a074057SBarry Smith   inputCtx->trianglelist               = NULL;
213a074057SBarry Smith   inputCtx->numberofholes              = 0;
223a074057SBarry Smith   inputCtx->holelist                   = NULL;
233a074057SBarry Smith   inputCtx->numberofregions            = 0;
243a074057SBarry Smith   inputCtx->regionlist                 = NULL;
253a074057SBarry Smith   PetscFunctionReturn(0);
263a074057SBarry Smith }
273a074057SBarry Smith 
283a074057SBarry Smith static PetscErrorCode InitOutput_Triangle(struct triangulateio *outputCtx)
293a074057SBarry Smith {
303a074057SBarry Smith   PetscFunctionBegin;
313a074057SBarry Smith   outputCtx->numberofpoints        = 0;
323a074057SBarry Smith   outputCtx->pointlist             = NULL;
333a074057SBarry Smith   outputCtx->pointattributelist    = NULL;
343a074057SBarry Smith   outputCtx->pointmarkerlist       = NULL;
353a074057SBarry Smith   outputCtx->numberoftriangles     = 0;
363a074057SBarry Smith   outputCtx->trianglelist          = NULL;
373a074057SBarry Smith   outputCtx->triangleattributelist = NULL;
383a074057SBarry Smith   outputCtx->neighborlist          = NULL;
393a074057SBarry Smith   outputCtx->segmentlist           = NULL;
403a074057SBarry Smith   outputCtx->segmentmarkerlist     = NULL;
413a074057SBarry Smith   outputCtx->numberofedges         = 0;
423a074057SBarry Smith   outputCtx->edgelist              = NULL;
433a074057SBarry Smith   outputCtx->edgemarkerlist        = NULL;
443a074057SBarry Smith   PetscFunctionReturn(0);
453a074057SBarry Smith }
463a074057SBarry Smith 
473a074057SBarry Smith static PetscErrorCode FiniOutput_Triangle(struct triangulateio *outputCtx)
483a074057SBarry Smith {
493a074057SBarry Smith   PetscFunctionBegin;
503a074057SBarry Smith   free(outputCtx->pointlist);
513a074057SBarry Smith   free(outputCtx->pointmarkerlist);
523a074057SBarry Smith   free(outputCtx->segmentlist);
533a074057SBarry Smith   free(outputCtx->segmentmarkerlist);
543a074057SBarry Smith   free(outputCtx->edgelist);
553a074057SBarry Smith   free(outputCtx->edgemarkerlist);
563a074057SBarry Smith   free(outputCtx->trianglelist);
573a074057SBarry Smith   free(outputCtx->neighborlist);
583a074057SBarry Smith   PetscFunctionReturn(0);
593a074057SBarry Smith }
603a074057SBarry Smith 
613a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerate_Triangle(DM boundary, PetscBool interpolate, DM *dm)
623a074057SBarry Smith {
633a074057SBarry Smith   MPI_Comm             comm;
643a074057SBarry Smith   DM_Plex             *mesh             = (DM_Plex *) boundary->data;
653a074057SBarry Smith   PetscInt             dim              = 2;
663a074057SBarry Smith   const PetscBool      createConvexHull = PETSC_FALSE;
673a074057SBarry Smith   const PetscBool      constrained      = PETSC_FALSE;
683a074057SBarry Smith   const char          *labelName        = "marker";
693a074057SBarry Smith   const char          *labelName2       = "Face Sets";
703a074057SBarry Smith   struct triangulateio in;
713a074057SBarry Smith   struct triangulateio out;
723a074057SBarry Smith   DMLabel              label, label2;
733a074057SBarry Smith   PetscInt             vStart, vEnd, v, eStart, eEnd, e;
743a074057SBarry Smith   PetscMPIInt          rank;
753a074057SBarry Smith 
763a074057SBarry Smith   PetscFunctionBegin;
779566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)boundary,&comm));
789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
799566063dSJacob Faibussowitsch   PetscCall(InitInput_Triangle(&in));
809566063dSJacob Faibussowitsch   PetscCall(InitOutput_Triangle(&out));
819566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd));
829566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(boundary, labelName,  &label));
839566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(boundary, labelName2, &label2));
843a074057SBarry Smith 
853a074057SBarry Smith   in.numberofpoints = vEnd - vStart;
863a074057SBarry Smith   if (in.numberofpoints > 0) {
873a074057SBarry Smith     PetscSection coordSection;
883a074057SBarry Smith     Vec          coordinates;
893a074057SBarry Smith     PetscScalar *array;
903a074057SBarry Smith 
919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints*dim, &in.pointlist));
929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints, &in.pointmarkerlist));
939566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(boundary, &coordinates));
949566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(boundary, &coordSection));
959566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &array));
963a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
973a074057SBarry Smith       const PetscInt idx = v - vStart;
98469e3fe5SMatthew G. Knepley       PetscInt       val, off, d;
993a074057SBarry Smith 
1009566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
1013a074057SBarry Smith       for (d = 0; d < dim; ++d) {
1023a074057SBarry Smith         in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
1033a074057SBarry Smith       }
104469e3fe5SMatthew G. Knepley       if (label) {
1059566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, v, &val));
106469e3fe5SMatthew G. Knepley         in.pointmarkerlist[idx] = val;
107469e3fe5SMatthew G. Knepley       }
1083a074057SBarry Smith     }
1099566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &array));
1103a074057SBarry Smith   }
1119566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(boundary, 0, &eStart, &eEnd));
1123a074057SBarry Smith   in.numberofsegments = eEnd - eStart;
1133a074057SBarry Smith   if (in.numberofsegments > 0) {
1149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofsegments*2, &in.segmentlist));
1159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofsegments, &in.segmentmarkerlist));
1163a074057SBarry Smith     for (e = eStart; e < eEnd; ++e) {
1173a074057SBarry Smith       const PetscInt  idx = e - eStart;
1183a074057SBarry Smith       const PetscInt *cone;
119469e3fe5SMatthew G. Knepley       PetscInt        val;
1203a074057SBarry Smith 
1219566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(boundary, e, &cone));
1223a074057SBarry Smith 
1233a074057SBarry Smith       in.segmentlist[idx*2+0] = cone[0] - vStart;
1243a074057SBarry Smith       in.segmentlist[idx*2+1] = cone[1] - vStart;
1253a074057SBarry Smith 
126469e3fe5SMatthew G. Knepley       if (label) {
1279566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, e, &val));
128469e3fe5SMatthew G. Knepley         in.segmentmarkerlist[idx] = val;
129469e3fe5SMatthew G. Knepley       }
1303a074057SBarry Smith     }
1313a074057SBarry Smith   }
1323a074057SBarry Smith #if 0 /* Do not currently support holes */
1333a074057SBarry Smith   PetscReal *holeCoords;
1343a074057SBarry Smith   PetscInt   h, d;
1353a074057SBarry Smith 
1369566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHoles(boundary, &in.numberofholes, &holeCords));
1373a074057SBarry Smith   if (in.numberofholes > 0) {
1389566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofholes*dim, &in.holelist));
1393a074057SBarry Smith     for (h = 0; h < in.numberofholes; ++h) {
1403a074057SBarry Smith       for (d = 0; d < dim; ++d) {
1413a074057SBarry Smith         in.holelist[h*dim+d] = holeCoords[h*dim+d];
1423a074057SBarry Smith       }
1433a074057SBarry Smith     }
1443a074057SBarry Smith   }
1453a074057SBarry Smith #endif
146dd400576SPatrick Sanan   if (rank == 0) {
1473a074057SBarry Smith     char args[32];
1483a074057SBarry Smith 
1493a074057SBarry Smith     /* Take away 'Q' for verbose output */
1509566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(args, "pqezQ"));
1519566063dSJacob Faibussowitsch     if (createConvexHull)   PetscCall(PetscStrcat(args, "c"));
1529566063dSJacob Faibussowitsch     if (constrained)        PetscCall(PetscStrcpy(args, "zepDQ"));
1533a074057SBarry Smith     if (mesh->triangleOpts) {triangulate(mesh->triangleOpts, &in, &out, NULL);}
1543a074057SBarry Smith     else                    {triangulate(args, &in, &out, NULL);}
1553a074057SBarry Smith   }
1569566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointlist));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointmarkerlist));
1589566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentlist));
1599566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentmarkerlist));
1609566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.holelist));
1613a074057SBarry Smith 
1623a074057SBarry Smith   {
1633a074057SBarry Smith     DMLabel          glabel      = NULL;
1643a074057SBarry Smith     DMLabel          glabel2     = NULL;
1653a074057SBarry Smith     const PetscInt   numCorners  = 3;
1663a074057SBarry Smith     const PetscInt   numCells    = out.numberoftriangles;
1673a074057SBarry Smith     const PetscInt   numVertices = out.numberofpoints;
168a4a685f2SJacob Faibussowitsch     PetscInt         *cells;
169a4a685f2SJacob Faibussowitsch     PetscReal        *meshCoords;
1703a074057SBarry Smith 
171a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) == sizeof (out.pointlist[0])) {
172a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *) out.pointlist;
173a4a685f2SJacob Faibussowitsch     } else {
174a4a685f2SJacob Faibussowitsch       PetscInt i;
175a4a685f2SJacob Faibussowitsch 
1769566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dim * numVertices,&meshCoords));
177a4a685f2SJacob Faibussowitsch       for (i = 0; i < dim * numVertices; i++) {
178a4a685f2SJacob Faibussowitsch         meshCoords[i] = (PetscReal) out.pointlist[i];
179a4a685f2SJacob Faibussowitsch       }
180a4a685f2SJacob Faibussowitsch     }
181a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) == sizeof (out.trianglelist[0])) {
182a4a685f2SJacob Faibussowitsch       cells = (PetscInt *) out.trianglelist;
183a4a685f2SJacob Faibussowitsch     } else {
184a4a685f2SJacob Faibussowitsch       PetscInt i;
185a4a685f2SJacob Faibussowitsch 
1869566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(numCells * numCorners, &cells));
187a4a685f2SJacob Faibussowitsch       for (i = 0; i < numCells * numCorners; i++) {
188a4a685f2SJacob Faibussowitsch         cells[i] = (PetscInt) out.trianglelist[i];
189a4a685f2SJacob Faibussowitsch       }
190a4a685f2SJacob Faibussowitsch     }
1919566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm));
192a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) != sizeof (out.pointlist[0])) {
1939566063dSJacob Faibussowitsch       PetscCall(PetscFree(meshCoords));
194a4a685f2SJacob Faibussowitsch     }
195a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) != sizeof (out.trianglelist[0])) {
1969566063dSJacob Faibussowitsch       PetscCall(PetscFree(cells));
197a4a685f2SJacob Faibussowitsch     }
198a4a685f2SJacob Faibussowitsch     if (label)  {
1999566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(*dm, labelName));
2009566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dm, labelName, &glabel));
201a4a685f2SJacob Faibussowitsch     }
202a4a685f2SJacob Faibussowitsch     if (label2) {
2039566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(*dm, labelName2));
2049566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dm, labelName2, &glabel2));
205a4a685f2SJacob Faibussowitsch     }
2063a074057SBarry Smith     /* Set labels */
2073a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
2083a074057SBarry Smith       if (out.pointmarkerlist[v]) {
2099566063dSJacob Faibussowitsch         if (glabel) PetscCall(DMLabelSetValue(glabel, v+numCells, out.pointmarkerlist[v]));
2103a074057SBarry Smith       }
2113a074057SBarry Smith     }
2123a074057SBarry Smith     if (interpolate) {
2133a074057SBarry Smith       for (e = 0; e < out.numberofedges; e++) {
2143a074057SBarry Smith         if (out.edgemarkerlist[e]) {
2153a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
2163a074057SBarry Smith           const PetscInt *edges;
2173a074057SBarry Smith           PetscInt        numEdges;
2183a074057SBarry Smith 
2199566063dSJacob Faibussowitsch           PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges));
2202c71b3e2SJacob Faibussowitsch           PetscCheckFalse(numEdges != 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
2219566063dSJacob Faibussowitsch           if (glabel)  PetscCall(DMLabelSetValue(glabel,  edges[0], out.edgemarkerlist[e]));
2229566063dSJacob Faibussowitsch           if (glabel2) PetscCall(DMLabelSetValue(glabel2, edges[0], out.edgemarkerlist[e]));
2239566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges));
2243a074057SBarry Smith         }
2253a074057SBarry Smith       }
2263a074057SBarry Smith     }
2279566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
2283a074057SBarry Smith   }
2293a074057SBarry Smith #if 0 /* Do not currently support holes */
2309566063dSJacob Faibussowitsch   PetscCall(DMPlexCopyHoles(*dm, boundary));
2313a074057SBarry Smith #endif
2329566063dSJacob Faibussowitsch   PetscCall(FiniOutput_Triangle(&out));
2333a074057SBarry Smith   PetscFunctionReturn(0);
2343a074057SBarry Smith }
2353a074057SBarry Smith 
236800b850cSBarry Smith PETSC_EXTERN PetscErrorCode DMPlexRefine_Triangle(DM dm, PetscReal *inmaxVolumes, DM *dmRefined)
2373a074057SBarry Smith {
2383a074057SBarry Smith   MPI_Comm             comm;
2393a074057SBarry Smith   PetscInt             dim       = 2;
2403a074057SBarry Smith   const char          *labelName = "marker";
2413a074057SBarry Smith   struct triangulateio in;
2423a074057SBarry Smith   struct triangulateio out;
2433a074057SBarry Smith   DMLabel              label;
2448457c761SMatthew G. Knepley   PetscInt             vStart, vEnd, v, gcStart, cStart, cEnd, c, depth, depthGlobal;
2453a074057SBarry Smith   PetscMPIInt          rank;
246800b850cSBarry Smith   double               *maxVolumes;
2473a074057SBarry Smith 
2483a074057SBarry Smith   PetscFunctionBegin;
2499566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
2509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2519566063dSJacob Faibussowitsch   PetscCall(InitInput_Triangle(&in));
2529566063dSJacob Faibussowitsch   PetscCall(InitOutput_Triangle(&out));
2539566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
254*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm));
2559566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2569566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, labelName, &label));
2573a074057SBarry Smith 
2583a074057SBarry Smith   in.numberofpoints = vEnd - vStart;
2593a074057SBarry Smith   if (in.numberofpoints > 0) {
2603a074057SBarry Smith     PetscSection coordSection;
2613a074057SBarry Smith     Vec          coordinates;
2623a074057SBarry Smith     PetscScalar *array;
2633a074057SBarry Smith 
2649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints*dim, &in.pointlist));
2659566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofpoints, &in.pointmarkerlist));
2669566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2679566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dm, &coordSection));
2689566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &array));
2693a074057SBarry Smith     for (v = vStart; v < vEnd; ++v) {
2703a074057SBarry Smith       const PetscInt idx = v - vStart;
271469e3fe5SMatthew G. Knepley       PetscInt       off, d, val;
2723a074057SBarry Smith 
2739566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
2743a074057SBarry Smith       for (d = 0; d < dim; ++d) {
2753a074057SBarry Smith         in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
2763a074057SBarry Smith       }
277469e3fe5SMatthew G. Knepley       if (label) {
2789566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, v, &val));
279469e3fe5SMatthew G. Knepley         in.pointmarkerlist[idx] = val;
280469e3fe5SMatthew G. Knepley       }
2813a074057SBarry Smith     }
2829566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &array));
2833a074057SBarry Smith   }
2849566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2859566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGhostCellStratum(dm, &gcStart, NULL));
2868457c761SMatthew G. Knepley   if (gcStart >= 0) cEnd = gcStart;
2873a074057SBarry Smith 
2883a074057SBarry Smith   in.numberofcorners   = 3;
2893a074057SBarry Smith   in.numberoftriangles = cEnd - cStart;
2903a074057SBarry Smith 
291800b850cSBarry Smith #if !defined(PETSC_USE_REAL_DOUBLE)
2929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cEnd - cStart,&maxVolumes));
293800b850cSBarry Smith   for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = (double)inmaxVolumes[c];
294800b850cSBarry Smith #else
295800b850cSBarry Smith   maxVolumes = inmaxVolumes;
296800b850cSBarry Smith #endif
297800b850cSBarry Smith 
2983a074057SBarry Smith   in.trianglearealist  = (double*) maxVolumes;
2993a074057SBarry Smith   if (in.numberoftriangles > 0) {
3009566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberoftriangles*in.numberofcorners, &in.trianglelist));
3013a074057SBarry Smith     for (c = cStart; c < cEnd; ++c) {
3023a074057SBarry Smith       const PetscInt idx      = c - cStart;
3033a074057SBarry Smith       PetscInt      *closure = NULL;
3043a074057SBarry Smith       PetscInt       closureSize;
3053a074057SBarry Smith 
3069566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
3072c71b3e2SJacob Faibussowitsch       PetscCheckFalse((closureSize != 4) && (closureSize != 7),comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a triangle, %D vertices in closure", closureSize);
3083a074057SBarry Smith       for (v = 0; v < 3; ++v) {
3093a074057SBarry Smith         in.trianglelist[idx*in.numberofcorners + v] = closure[(v+closureSize-3)*2] - vStart;
3103a074057SBarry Smith       }
3119566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
3123a074057SBarry Smith     }
3133a074057SBarry Smith   }
3143a074057SBarry Smith   /* TODO: Segment markers are missing on input */
3153a074057SBarry Smith #if 0 /* Do not currently support holes */
3163a074057SBarry Smith   PetscReal *holeCoords;
3173a074057SBarry Smith   PetscInt   h, d;
3183a074057SBarry Smith 
3199566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHoles(boundary, &in.numberofholes, &holeCords));
3203a074057SBarry Smith   if (in.numberofholes > 0) {
3219566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(in.numberofholes*dim, &in.holelist));
3223a074057SBarry Smith     for (h = 0; h < in.numberofholes; ++h) {
3233a074057SBarry Smith       for (d = 0; d < dim; ++d) {
3243a074057SBarry Smith         in.holelist[h*dim+d] = holeCoords[h*dim+d];
3253a074057SBarry Smith       }
3263a074057SBarry Smith     }
3273a074057SBarry Smith   }
3283a074057SBarry Smith #endif
329dd400576SPatrick Sanan   if (rank == 0) {
3303a074057SBarry Smith     char args[32];
3313a074057SBarry Smith 
3323a074057SBarry Smith     /* Take away 'Q' for verbose output */
3339566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(args, "pqezQra"));
3343a074057SBarry Smith     triangulate(args, &in, &out, NULL);
3353a074057SBarry Smith   }
3369566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointlist));
3379566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointmarkerlist));
3389566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentlist));
3399566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.segmentmarkerlist));
3409566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.trianglelist));
3413a074057SBarry Smith 
3423a074057SBarry Smith   {
3433a074057SBarry Smith     DMLabel          rlabel      = NULL;
3443a074057SBarry Smith     const PetscInt   numCorners  = 3;
3453a074057SBarry Smith     const PetscInt   numCells    = out.numberoftriangles;
3463a074057SBarry Smith     const PetscInt   numVertices = out.numberofpoints;
347a4a685f2SJacob Faibussowitsch     PetscInt         *cells;
348a4a685f2SJacob Faibussowitsch     PetscReal        *meshCoords;
3493a074057SBarry Smith     PetscBool        interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
3503a074057SBarry Smith 
351a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) == sizeof (out.pointlist[0])) {
352a4a685f2SJacob Faibussowitsch       meshCoords = (PetscReal *) out.pointlist;
353a4a685f2SJacob Faibussowitsch     } else {
354a4a685f2SJacob Faibussowitsch       PetscInt i;
355a4a685f2SJacob Faibussowitsch 
3569566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dim * numVertices,&meshCoords));
357a4a685f2SJacob Faibussowitsch       for (i = 0; i < dim * numVertices; i++) {
358a4a685f2SJacob Faibussowitsch         meshCoords[i] = (PetscReal) out.pointlist[i];
359a4a685f2SJacob Faibussowitsch       }
360a4a685f2SJacob Faibussowitsch     }
361a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) == sizeof (out.trianglelist[0])) {
362a4a685f2SJacob Faibussowitsch       cells = (PetscInt *) out.trianglelist;
363a4a685f2SJacob Faibussowitsch     } else {
364a4a685f2SJacob Faibussowitsch       PetscInt i;
365a4a685f2SJacob Faibussowitsch 
3669566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(numCells * numCorners, &cells));
367a4a685f2SJacob Faibussowitsch       for (i = 0; i < numCells * numCorners; i++) {
368a4a685f2SJacob Faibussowitsch         cells[i] = (PetscInt) out.trianglelist[i];
369a4a685f2SJacob Faibussowitsch       }
370a4a685f2SJacob Faibussowitsch     }
371a4a685f2SJacob Faibussowitsch 
3729566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined));
373a4a685f2SJacob Faibussowitsch     if (label) {
3749566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(*dmRefined, labelName));
3759566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dmRefined, labelName, &rlabel));
376a4a685f2SJacob Faibussowitsch     }
377a4a685f2SJacob Faibussowitsch     if (sizeof (PetscReal) != sizeof (out.pointlist[0])) {
3789566063dSJacob Faibussowitsch       PetscCall(PetscFree(meshCoords));
379a4a685f2SJacob Faibussowitsch     }
380a4a685f2SJacob Faibussowitsch     if (sizeof (PetscInt) != sizeof (out.trianglelist[0])) {
3819566063dSJacob Faibussowitsch       PetscCall(PetscFree(cells));
382a4a685f2SJacob Faibussowitsch     }
3833a074057SBarry Smith     /* Set labels */
3843a074057SBarry Smith     for (v = 0; v < numVertices; ++v) {
3853a074057SBarry Smith       if (out.pointmarkerlist[v]) {
3869566063dSJacob Faibussowitsch         if (rlabel) PetscCall(DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]));
3873a074057SBarry Smith       }
3883a074057SBarry Smith     }
3893a074057SBarry Smith     if (interpolate) {
3903a074057SBarry Smith       PetscInt e;
3913a074057SBarry Smith 
3923a074057SBarry Smith       for (e = 0; e < out.numberofedges; e++) {
3933a074057SBarry Smith         if (out.edgemarkerlist[e]) {
3943a074057SBarry Smith           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
3953a074057SBarry Smith           const PetscInt *edges;
3963a074057SBarry Smith           PetscInt        numEdges;
3973a074057SBarry Smith 
3989566063dSJacob Faibussowitsch           PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges));
3992c71b3e2SJacob Faibussowitsch           PetscCheckFalse(numEdges != 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
4009566063dSJacob Faibussowitsch           if (rlabel) PetscCall(DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]));
4019566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges));
4023a074057SBarry Smith         }
4033a074057SBarry Smith       }
4043a074057SBarry Smith     }
4059566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE));
4063a074057SBarry Smith   }
4073a074057SBarry Smith #if 0 /* Do not currently support holes */
4089566063dSJacob Faibussowitsch   PetscCall(DMPlexCopyHoles(*dm, boundary));
4093a074057SBarry Smith #endif
4109566063dSJacob Faibussowitsch   PetscCall(FiniOutput_Triangle(&out));
411800b850cSBarry Smith #if !defined(PETSC_USE_REAL_DOUBLE)
4129566063dSJacob Faibussowitsch   PetscCall(PetscFree(maxVolumes));
413800b850cSBarry Smith #endif
4143a074057SBarry Smith   PetscFunctionReturn(0);
4153a074057SBarry Smith }
416