xref: /petsc/src/dm/impls/plex/plexgenerate.c (revision 94ef8dde638caef1d0cd84a7dc8a2db65fcda8b6)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2d9deefdfSMatthew G. Knepley 
3d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexInvertCell_Internal(PetscInt dim, PetscInt numCorners, PetscInt cone[])
4d9deefdfSMatthew G. Knepley {
5d9deefdfSMatthew G. Knepley   int tmpc;
6d9deefdfSMatthew G. Knepley 
7d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
8d9deefdfSMatthew G. Knepley   if (dim != 3) PetscFunctionReturn(0);
9d9deefdfSMatthew G. Knepley   switch (numCorners) {
10d9deefdfSMatthew G. Knepley   case 4:
11d9deefdfSMatthew G. Knepley     tmpc    = cone[0];
12d9deefdfSMatthew G. Knepley     cone[0] = cone[1];
13d9deefdfSMatthew G. Knepley     cone[1] = tmpc;
14d9deefdfSMatthew G. Knepley     break;
15d9deefdfSMatthew G. Knepley   case 8:
16d9deefdfSMatthew G. Knepley     tmpc    = cone[1];
17d9deefdfSMatthew G. Knepley     cone[1] = cone[3];
18d9deefdfSMatthew G. Knepley     cone[3] = tmpc;
19d9deefdfSMatthew G. Knepley     break;
20d9deefdfSMatthew G. Knepley   default: break;
21d9deefdfSMatthew G. Knepley   }
22d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
23d9deefdfSMatthew G. Knepley }
24d9deefdfSMatthew G. Knepley 
25d9deefdfSMatthew G. Knepley /*@C
26d9deefdfSMatthew G. Knepley   DMPlexInvertCell - This flips tetrahedron and hexahedron orientation since Plex stores them internally with outward normals. Other cells are left untouched.
27d9deefdfSMatthew G. Knepley 
28d9deefdfSMatthew G. Knepley   Input Parameters:
29d9deefdfSMatthew G. Knepley + numCorners - The number of vertices in a cell
30d9deefdfSMatthew G. Knepley - cone - The incoming cone
31d9deefdfSMatthew G. Knepley 
32d9deefdfSMatthew G. Knepley   Output Parameter:
33d9deefdfSMatthew G. Knepley . cone - The inverted cone (in-place)
34d9deefdfSMatthew G. Knepley 
35d9deefdfSMatthew G. Knepley   Level: developer
36d9deefdfSMatthew G. Knepley 
37d9deefdfSMatthew G. Knepley .seealso: DMPlexGenerate()
38d9deefdfSMatthew G. Knepley @*/
39d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexInvertCell(PetscInt dim, PetscInt numCorners, int cone[])
40d9deefdfSMatthew G. Knepley {
41d9deefdfSMatthew G. Knepley   int tmpc;
42d9deefdfSMatthew G. Knepley 
43d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
44d9deefdfSMatthew G. Knepley   if (dim != 3) PetscFunctionReturn(0);
45d9deefdfSMatthew G. Knepley   switch (numCorners) {
46d9deefdfSMatthew G. Knepley   case 4:
47d9deefdfSMatthew G. Knepley     tmpc    = cone[0];
48d9deefdfSMatthew G. Knepley     cone[0] = cone[1];
49d9deefdfSMatthew G. Knepley     cone[1] = tmpc;
50d9deefdfSMatthew G. Knepley     break;
51d9deefdfSMatthew G. Knepley   case 8:
52d9deefdfSMatthew G. Knepley     tmpc    = cone[1];
53d9deefdfSMatthew G. Knepley     cone[1] = cone[3];
54d9deefdfSMatthew G. Knepley     cone[3] = tmpc;
55d9deefdfSMatthew G. Knepley     break;
56d9deefdfSMatthew G. Knepley   default: break;
57d9deefdfSMatthew G. Knepley   }
58d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
59d9deefdfSMatthew G. Knepley }
60d9deefdfSMatthew G. Knepley 
61d9deefdfSMatthew G. Knepley /* This is to fix the tetrahedron orientation from TetGen */
62d9deefdfSMatthew G. Knepley PETSC_UNUSED static PetscErrorCode DMPlexInvertCells_Internal(PetscInt dim, PetscInt numCells, PetscInt numCorners, int cells[])
63d9deefdfSMatthew G. Knepley {
64d9deefdfSMatthew G. Knepley   PetscInt       bound = numCells*numCorners, coff;
65d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
66d9deefdfSMatthew G. Knepley 
67d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
68d9deefdfSMatthew G. Knepley   for (coff = 0; coff < bound; coff += numCorners) {
69d9deefdfSMatthew G. Knepley     ierr = DMPlexInvertCell(dim, numCorners, &cells[coff]);CHKERRQ(ierr);
70d9deefdfSMatthew G. Knepley   }
71d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
72d9deefdfSMatthew G. Knepley }
73d9deefdfSMatthew G. Knepley 
74*94ef8ddeSSatish Balay /*@C
75d9deefdfSMatthew G. Knepley   DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator
76d9deefdfSMatthew G. Knepley 
77d9deefdfSMatthew G. Knepley   Not Collective
78d9deefdfSMatthew G. Knepley 
79d9deefdfSMatthew G. Knepley   Inputs Parameters:
80d9deefdfSMatthew G. Knepley + dm - The DMPlex object
81d9deefdfSMatthew G. Knepley - opts - The command line options
82d9deefdfSMatthew G. Knepley 
83d9deefdfSMatthew G. Knepley   Level: developer
84d9deefdfSMatthew G. Knepley 
85d9deefdfSMatthew G. Knepley .keywords: mesh, points
86d9deefdfSMatthew G. Knepley .seealso: DMPlexTetgenSetOptions(), DMPlexGenerate()
87d9deefdfSMatthew G. Knepley @*/
88d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts)
89d9deefdfSMatthew G. Knepley {
90d9deefdfSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
91d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
92d9deefdfSMatthew G. Knepley 
93d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
94d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
95d9deefdfSMatthew G. Knepley   PetscValidPointer(opts, 2);
96606ac3a5SMatthew G. Knepley   ierr = PetscFree(mesh->triangleOpts);CHKERRQ(ierr);
97606ac3a5SMatthew G. Knepley   ierr = PetscStrallocpy(opts, &mesh->triangleOpts);CHKERRQ(ierr);
98d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
99d9deefdfSMatthew G. Knepley }
100d9deefdfSMatthew G. Knepley 
101*94ef8ddeSSatish Balay /*@C
102d9deefdfSMatthew G. Knepley   DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator
103d9deefdfSMatthew G. Knepley 
104d9deefdfSMatthew G. Knepley   Not Collective
105d9deefdfSMatthew G. Knepley 
106d9deefdfSMatthew G. Knepley   Inputs Parameters:
107d9deefdfSMatthew G. Knepley + dm - The DMPlex object
108d9deefdfSMatthew G. Knepley - opts - The command line options
109d9deefdfSMatthew G. Knepley 
110d9deefdfSMatthew G. Knepley   Level: developer
111d9deefdfSMatthew G. Knepley 
112d9deefdfSMatthew G. Knepley .keywords: mesh, points
113d9deefdfSMatthew G. Knepley .seealso: DMPlexTriangleSetOptions(), DMPlexGenerate()
114d9deefdfSMatthew G. Knepley @*/
115d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts)
116d9deefdfSMatthew G. Knepley {
117d9deefdfSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
118d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
119d9deefdfSMatthew G. Knepley 
120d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
121d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
122d9deefdfSMatthew G. Knepley   PetscValidPointer(opts, 2);
123606ac3a5SMatthew G. Knepley   ierr = PetscFree(mesh->tetgenOpts);CHKERRQ(ierr);
124606ac3a5SMatthew G. Knepley   ierr = PetscStrallocpy(opts, &mesh->tetgenOpts);CHKERRQ(ierr);
125d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
126d9deefdfSMatthew G. Knepley }
127d9deefdfSMatthew G. Knepley 
128d9deefdfSMatthew G. Knepley #if defined(PETSC_HAVE_TRIANGLE)
129d9deefdfSMatthew G. Knepley #include <triangle.h>
130d9deefdfSMatthew G. Knepley 
131e9db6b00SMatthew G. Knepley static PetscErrorCode InitInput_Triangle(struct triangulateio *inputCtx)
132d9deefdfSMatthew G. Knepley {
133d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
134d9deefdfSMatthew G. Knepley   inputCtx->numberofpoints             = 0;
135d9deefdfSMatthew G. Knepley   inputCtx->numberofpointattributes    = 0;
136d9deefdfSMatthew G. Knepley   inputCtx->pointlist                  = NULL;
137d9deefdfSMatthew G. Knepley   inputCtx->pointattributelist         = NULL;
138d9deefdfSMatthew G. Knepley   inputCtx->pointmarkerlist            = NULL;
139d9deefdfSMatthew G. Knepley   inputCtx->numberofsegments           = 0;
140d9deefdfSMatthew G. Knepley   inputCtx->segmentlist                = NULL;
141d9deefdfSMatthew G. Knepley   inputCtx->segmentmarkerlist          = NULL;
142d9deefdfSMatthew G. Knepley   inputCtx->numberoftriangleattributes = 0;
143d9deefdfSMatthew G. Knepley   inputCtx->trianglelist               = NULL;
144d9deefdfSMatthew G. Knepley   inputCtx->numberofholes              = 0;
145d9deefdfSMatthew G. Knepley   inputCtx->holelist                   = NULL;
146d9deefdfSMatthew G. Knepley   inputCtx->numberofregions            = 0;
147d9deefdfSMatthew G. Knepley   inputCtx->regionlist                 = NULL;
148d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
149d9deefdfSMatthew G. Knepley }
150d9deefdfSMatthew G. Knepley 
151e9db6b00SMatthew G. Knepley static PetscErrorCode InitOutput_Triangle(struct triangulateio *outputCtx)
152d9deefdfSMatthew G. Knepley {
153d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
154d9deefdfSMatthew G. Knepley   outputCtx->numberofpoints        = 0;
155d9deefdfSMatthew G. Knepley   outputCtx->pointlist             = NULL;
156d9deefdfSMatthew G. Knepley   outputCtx->pointattributelist    = NULL;
157d9deefdfSMatthew G. Knepley   outputCtx->pointmarkerlist       = NULL;
158d9deefdfSMatthew G. Knepley   outputCtx->numberoftriangles     = 0;
159d9deefdfSMatthew G. Knepley   outputCtx->trianglelist          = NULL;
160d9deefdfSMatthew G. Knepley   outputCtx->triangleattributelist = NULL;
161d9deefdfSMatthew G. Knepley   outputCtx->neighborlist          = NULL;
162d9deefdfSMatthew G. Knepley   outputCtx->segmentlist           = NULL;
163d9deefdfSMatthew G. Knepley   outputCtx->segmentmarkerlist     = NULL;
164d9deefdfSMatthew G. Knepley   outputCtx->numberofedges         = 0;
165d9deefdfSMatthew G. Knepley   outputCtx->edgelist              = NULL;
166d9deefdfSMatthew G. Knepley   outputCtx->edgemarkerlist        = NULL;
167d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
168d9deefdfSMatthew G. Knepley }
169d9deefdfSMatthew G. Knepley 
170e9db6b00SMatthew G. Knepley static PetscErrorCode FiniOutput_Triangle(struct triangulateio *outputCtx)
171d9deefdfSMatthew G. Knepley {
172d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
173d9deefdfSMatthew G. Knepley   free(outputCtx->pointlist);
174d9deefdfSMatthew G. Knepley   free(outputCtx->pointmarkerlist);
175d9deefdfSMatthew G. Knepley   free(outputCtx->segmentlist);
176d9deefdfSMatthew G. Knepley   free(outputCtx->segmentmarkerlist);
177d9deefdfSMatthew G. Knepley   free(outputCtx->edgelist);
178d9deefdfSMatthew G. Knepley   free(outputCtx->edgemarkerlist);
179d9deefdfSMatthew G. Knepley   free(outputCtx->trianglelist);
180d9deefdfSMatthew G. Knepley   free(outputCtx->neighborlist);
181d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
182d9deefdfSMatthew G. Knepley }
183d9deefdfSMatthew G. Knepley 
184e9db6b00SMatthew G. Knepley static PetscErrorCode DMPlexGenerate_Triangle(DM boundary, PetscBool interpolate, DM *dm)
185d9deefdfSMatthew G. Knepley {
186d9deefdfSMatthew G. Knepley   MPI_Comm             comm;
187d9deefdfSMatthew G. Knepley   DM_Plex             *mesh             = (DM_Plex *) boundary->data;
188d9deefdfSMatthew G. Knepley   PetscInt             dim              = 2;
189d9deefdfSMatthew G. Knepley   const PetscBool      createConvexHull = PETSC_FALSE;
190d9deefdfSMatthew G. Knepley   const PetscBool      constrained      = PETSC_FALSE;
191d988aadeSMatthew G. Knepley   const char          *labelName        = "marker";
192d9deefdfSMatthew G. Knepley   struct triangulateio in;
193d9deefdfSMatthew G. Knepley   struct triangulateio out;
194d988aadeSMatthew G. Knepley   DMLabel              label;
195d9deefdfSMatthew G. Knepley   PetscInt             vStart, vEnd, v, eStart, eEnd, e;
196d9deefdfSMatthew G. Knepley   PetscMPIInt          rank;
197d9deefdfSMatthew G. Knepley   PetscErrorCode       ierr;
198d9deefdfSMatthew G. Knepley 
199d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
200d9deefdfSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr);
201d9deefdfSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
202d9deefdfSMatthew G. Knepley   ierr = InitInput_Triangle(&in);CHKERRQ(ierr);
203d9deefdfSMatthew G. Knepley   ierr = InitOutput_Triangle(&out);CHKERRQ(ierr);
204d9deefdfSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr);
205c58f1c22SToby Isaac   ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr);
206d9deefdfSMatthew G. Knepley 
207d9deefdfSMatthew G. Knepley   in.numberofpoints = vEnd - vStart;
208d9deefdfSMatthew G. Knepley   if (in.numberofpoints > 0) {
209d9deefdfSMatthew G. Knepley     PetscSection coordSection;
210d9deefdfSMatthew G. Knepley     Vec          coordinates;
211d9deefdfSMatthew G. Knepley     PetscScalar *array;
212d9deefdfSMatthew G. Knepley 
213d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in.numberofpoints*dim, &in.pointlist);CHKERRQ(ierr);
214d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in.numberofpoints, &in.pointmarkerlist);CHKERRQ(ierr);
215d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr);
216d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr);
217d9deefdfSMatthew G. Knepley     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
218d9deefdfSMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
219d9deefdfSMatthew G. Knepley       const PetscInt idx = v - vStart;
220d9deefdfSMatthew G. Knepley       PetscInt       off, d;
221d9deefdfSMatthew G. Knepley 
222d9deefdfSMatthew G. Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
223d9deefdfSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
224d9deefdfSMatthew G. Knepley         in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
225d9deefdfSMatthew G. Knepley       }
226d988aadeSMatthew G. Knepley       if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);}
227d9deefdfSMatthew G. Knepley     }
228d9deefdfSMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
229d9deefdfSMatthew G. Knepley   }
230d9deefdfSMatthew G. Knepley   ierr  = DMPlexGetHeightStratum(boundary, 0, &eStart, &eEnd);CHKERRQ(ierr);
231d9deefdfSMatthew G. Knepley   in.numberofsegments = eEnd - eStart;
232d9deefdfSMatthew G. Knepley   if (in.numberofsegments > 0) {
233d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in.numberofsegments*2, &in.segmentlist);CHKERRQ(ierr);
234d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in.numberofsegments, &in.segmentmarkerlist);CHKERRQ(ierr);
235d9deefdfSMatthew G. Knepley     for (e = eStart; e < eEnd; ++e) {
236d9deefdfSMatthew G. Knepley       const PetscInt  idx = e - eStart;
237d9deefdfSMatthew G. Knepley       const PetscInt *cone;
238d9deefdfSMatthew G. Knepley 
239d9deefdfSMatthew G. Knepley       ierr = DMPlexGetCone(boundary, e, &cone);CHKERRQ(ierr);
240d9deefdfSMatthew G. Knepley 
241d9deefdfSMatthew G. Knepley       in.segmentlist[idx*2+0] = cone[0] - vStart;
242d9deefdfSMatthew G. Knepley       in.segmentlist[idx*2+1] = cone[1] - vStart;
243d9deefdfSMatthew G. Knepley 
244d988aadeSMatthew G. Knepley       if (label) {ierr = DMLabelGetValue(label, e, &in.segmentmarkerlist[idx]);CHKERRQ(ierr);}
245d9deefdfSMatthew G. Knepley     }
246d9deefdfSMatthew G. Knepley   }
247d9deefdfSMatthew G. Knepley #if 0 /* Do not currently support holes */
248d9deefdfSMatthew G. Knepley   PetscReal *holeCoords;
249d9deefdfSMatthew G. Knepley   PetscInt   h, d;
250d9deefdfSMatthew G. Knepley 
251d9deefdfSMatthew G. Knepley   ierr = DMPlexGetHoles(boundary, &in.numberofholes, &holeCords);CHKERRQ(ierr);
252d9deefdfSMatthew G. Knepley   if (in.numberofholes > 0) {
253d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in.numberofholes*dim, &in.holelist);CHKERRQ(ierr);
254d9deefdfSMatthew G. Knepley     for (h = 0; h < in.numberofholes; ++h) {
255d9deefdfSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
256d9deefdfSMatthew G. Knepley         in.holelist[h*dim+d] = holeCoords[h*dim+d];
257d9deefdfSMatthew G. Knepley       }
258d9deefdfSMatthew G. Knepley     }
259d9deefdfSMatthew G. Knepley   }
260d9deefdfSMatthew G. Knepley #endif
261d9deefdfSMatthew G. Knepley   if (!rank) {
262d9deefdfSMatthew G. Knepley     char args[32];
263d9deefdfSMatthew G. Knepley 
264d9deefdfSMatthew G. Knepley     /* Take away 'Q' for verbose output */
265d9deefdfSMatthew G. Knepley     ierr = PetscStrcpy(args, "pqezQ");CHKERRQ(ierr);
266d9deefdfSMatthew G. Knepley     if (createConvexHull)   {ierr = PetscStrcat(args, "c");CHKERRQ(ierr);}
267d9deefdfSMatthew G. Knepley     if (constrained)        {ierr = PetscStrcpy(args, "zepDQ");CHKERRQ(ierr);}
268d9deefdfSMatthew G. Knepley     if (mesh->triangleOpts) {triangulate(mesh->triangleOpts, &in, &out, NULL);}
269d9deefdfSMatthew G. Knepley     else                    {triangulate(args, &in, &out, NULL);}
270d9deefdfSMatthew G. Knepley   }
271d9deefdfSMatthew G. Knepley   ierr = PetscFree(in.pointlist);CHKERRQ(ierr);
272d9deefdfSMatthew G. Knepley   ierr = PetscFree(in.pointmarkerlist);CHKERRQ(ierr);
273d9deefdfSMatthew G. Knepley   ierr = PetscFree(in.segmentlist);CHKERRQ(ierr);
274d9deefdfSMatthew G. Knepley   ierr = PetscFree(in.segmentmarkerlist);CHKERRQ(ierr);
275d9deefdfSMatthew G. Knepley   ierr = PetscFree(in.holelist);CHKERRQ(ierr);
276d9deefdfSMatthew G. Knepley 
277d9deefdfSMatthew G. Knepley   {
278d988aadeSMatthew G. Knepley     DMLabel        glabel      = NULL;
279d9deefdfSMatthew G. Knepley     const PetscInt numCorners  = 3;
280d9deefdfSMatthew G. Knepley     const PetscInt numCells    = out.numberoftriangles;
281d9deefdfSMatthew G. Knepley     const PetscInt numVertices = out.numberofpoints;
282d9deefdfSMatthew G. Knepley     const int     *cells      = out.trianglelist;
283d9deefdfSMatthew G. Knepley     const double  *meshCoords = out.pointlist;
284d9deefdfSMatthew G. Knepley 
285d9deefdfSMatthew G. Knepley     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr);
286c58f1c22SToby Isaac     if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);}
287d9deefdfSMatthew G. Knepley     /* Set labels */
288d9deefdfSMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
289d9deefdfSMatthew G. Knepley       if (out.pointmarkerlist[v]) {
290d988aadeSMatthew G. Knepley         if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);}
291d9deefdfSMatthew G. Knepley       }
292d9deefdfSMatthew G. Knepley     }
293d9deefdfSMatthew G. Knepley     if (interpolate) {
294d9deefdfSMatthew G. Knepley       for (e = 0; e < out.numberofedges; e++) {
295d9deefdfSMatthew G. Knepley         if (out.edgemarkerlist[e]) {
296d9deefdfSMatthew G. Knepley           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
297d9deefdfSMatthew G. Knepley           const PetscInt *edges;
298d9deefdfSMatthew G. Knepley           PetscInt        numEdges;
299d9deefdfSMatthew G. Knepley 
300d9deefdfSMatthew G. Knepley           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
301d9deefdfSMatthew G. Knepley           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
302d988aadeSMatthew G. Knepley           if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
303d9deefdfSMatthew G. Knepley           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
304d9deefdfSMatthew G. Knepley         }
305d9deefdfSMatthew G. Knepley       }
306d9deefdfSMatthew G. Knepley     }
307d9deefdfSMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
308d9deefdfSMatthew G. Knepley   }
309d9deefdfSMatthew G. Knepley #if 0 /* Do not currently support holes */
310d9deefdfSMatthew G. Knepley   ierr = DMPlexCopyHoles(*dm, boundary);CHKERRQ(ierr);
311d9deefdfSMatthew G. Knepley #endif
312d9deefdfSMatthew G. Knepley   ierr = FiniOutput_Triangle(&out);CHKERRQ(ierr);
313d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
314d9deefdfSMatthew G. Knepley }
315d9deefdfSMatthew G. Knepley 
316e9db6b00SMatthew G. Knepley static PetscErrorCode DMPlexRefine_Triangle(DM dm, double *maxVolumes, DM *dmRefined)
317d9deefdfSMatthew G. Knepley {
318d9deefdfSMatthew G. Knepley   MPI_Comm             comm;
319d9deefdfSMatthew G. Knepley   PetscInt             dim       = 2;
320d988aadeSMatthew G. Knepley   const char          *labelName = "marker";
321d9deefdfSMatthew G. Knepley   struct triangulateio in;
322d9deefdfSMatthew G. Knepley   struct triangulateio out;
323d988aadeSMatthew G. Knepley   DMLabel              label;
324d9deefdfSMatthew G. Knepley   PetscInt             vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal;
325d9deefdfSMatthew G. Knepley   PetscMPIInt          rank;
326d9deefdfSMatthew G. Knepley   PetscErrorCode       ierr;
327d9deefdfSMatthew G. Knepley 
328d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
329d9deefdfSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
330d9deefdfSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
331d9deefdfSMatthew G. Knepley   ierr = InitInput_Triangle(&in);CHKERRQ(ierr);
332d9deefdfSMatthew G. Knepley   ierr = InitOutput_Triangle(&out);CHKERRQ(ierr);
333d9deefdfSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
334b2566f29SBarry Smith   ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
335d9deefdfSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
336c58f1c22SToby Isaac   ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr);
337d9deefdfSMatthew G. Knepley 
338d9deefdfSMatthew G. Knepley   in.numberofpoints = vEnd - vStart;
339d9deefdfSMatthew G. Knepley   if (in.numberofpoints > 0) {
340d9deefdfSMatthew G. Knepley     PetscSection coordSection;
341d9deefdfSMatthew G. Knepley     Vec          coordinates;
342d9deefdfSMatthew G. Knepley     PetscScalar *array;
343d9deefdfSMatthew G. Knepley 
344d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in.numberofpoints*dim, &in.pointlist);CHKERRQ(ierr);
345d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in.numberofpoints, &in.pointmarkerlist);CHKERRQ(ierr);
346d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
347d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
348d9deefdfSMatthew G. Knepley     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
349d9deefdfSMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
350d9deefdfSMatthew G. Knepley       const PetscInt idx = v - vStart;
351d9deefdfSMatthew G. Knepley       PetscInt       off, d;
352d9deefdfSMatthew G. Knepley 
353d9deefdfSMatthew G. Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
354d9deefdfSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
355d9deefdfSMatthew G. Knepley         in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
356d9deefdfSMatthew G. Knepley       }
357d988aadeSMatthew G. Knepley       if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);}
358d9deefdfSMatthew G. Knepley     }
359d9deefdfSMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
360d9deefdfSMatthew G. Knepley   }
361d9deefdfSMatthew G. Knepley   ierr  = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
362d9deefdfSMatthew G. Knepley 
363d9deefdfSMatthew G. Knepley   in.numberofcorners   = 3;
364d9deefdfSMatthew G. Knepley   in.numberoftriangles = cEnd - cStart;
365d9deefdfSMatthew G. Knepley 
366d9deefdfSMatthew G. Knepley   in.trianglearealist  = (double*) maxVolumes;
367d9deefdfSMatthew G. Knepley   if (in.numberoftriangles > 0) {
368d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in.numberoftriangles*in.numberofcorners, &in.trianglelist);CHKERRQ(ierr);
369d9deefdfSMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
370d9deefdfSMatthew G. Knepley       const PetscInt idx      = c - cStart;
371d9deefdfSMatthew G. Knepley       PetscInt      *closure = NULL;
372d9deefdfSMatthew G. Knepley       PetscInt       closureSize;
373d9deefdfSMatthew G. Knepley 
374d9deefdfSMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
375d9deefdfSMatthew G. Knepley       if ((closureSize != 4) && (closureSize != 7)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a triangle, %D vertices in closure", closureSize);
376d9deefdfSMatthew G. Knepley       for (v = 0; v < 3; ++v) {
377d9deefdfSMatthew G. Knepley         in.trianglelist[idx*in.numberofcorners + v] = closure[(v+closureSize-3)*2] - vStart;
378d9deefdfSMatthew G. Knepley       }
379d9deefdfSMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
380d9deefdfSMatthew G. Knepley     }
381d9deefdfSMatthew G. Knepley   }
382d9deefdfSMatthew G. Knepley   /* TODO: Segment markers are missing on input */
383d9deefdfSMatthew G. Knepley #if 0 /* Do not currently support holes */
384d9deefdfSMatthew G. Knepley   PetscReal *holeCoords;
385d9deefdfSMatthew G. Knepley   PetscInt   h, d;
386d9deefdfSMatthew G. Knepley 
387d9deefdfSMatthew G. Knepley   ierr = DMPlexGetHoles(boundary, &in.numberofholes, &holeCords);CHKERRQ(ierr);
388d9deefdfSMatthew G. Knepley   if (in.numberofholes > 0) {
389d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in.numberofholes*dim, &in.holelist);CHKERRQ(ierr);
390d9deefdfSMatthew G. Knepley     for (h = 0; h < in.numberofholes; ++h) {
391d9deefdfSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
392d9deefdfSMatthew G. Knepley         in.holelist[h*dim+d] = holeCoords[h*dim+d];
393d9deefdfSMatthew G. Knepley       }
394d9deefdfSMatthew G. Knepley     }
395d9deefdfSMatthew G. Knepley   }
396d9deefdfSMatthew G. Knepley #endif
397d9deefdfSMatthew G. Knepley   if (!rank) {
398d9deefdfSMatthew G. Knepley     char args[32];
399d9deefdfSMatthew G. Knepley 
400d9deefdfSMatthew G. Knepley     /* Take away 'Q' for verbose output */
401d9deefdfSMatthew G. Knepley     ierr = PetscStrcpy(args, "pqezQra");CHKERRQ(ierr);
402d9deefdfSMatthew G. Knepley     triangulate(args, &in, &out, NULL);
403d9deefdfSMatthew G. Knepley   }
404d9deefdfSMatthew G. Knepley   ierr = PetscFree(in.pointlist);CHKERRQ(ierr);
405d9deefdfSMatthew G. Knepley   ierr = PetscFree(in.pointmarkerlist);CHKERRQ(ierr);
406d9deefdfSMatthew G. Knepley   ierr = PetscFree(in.segmentlist);CHKERRQ(ierr);
407d9deefdfSMatthew G. Knepley   ierr = PetscFree(in.segmentmarkerlist);CHKERRQ(ierr);
408d9deefdfSMatthew G. Knepley   ierr = PetscFree(in.trianglelist);CHKERRQ(ierr);
409d9deefdfSMatthew G. Knepley 
410d9deefdfSMatthew G. Knepley   {
411d988aadeSMatthew G. Knepley     DMLabel        rlabel      = NULL;
412d9deefdfSMatthew G. Knepley     const PetscInt numCorners  = 3;
413d9deefdfSMatthew G. Knepley     const PetscInt numCells    = out.numberoftriangles;
414d9deefdfSMatthew G. Knepley     const PetscInt numVertices = out.numberofpoints;
415d9deefdfSMatthew G. Knepley     const int     *cells      = out.trianglelist;
416d9deefdfSMatthew G. Knepley     const double  *meshCoords = out.pointlist;
417d9deefdfSMatthew G. Knepley     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
418d9deefdfSMatthew G. Knepley 
419d9deefdfSMatthew G. Knepley     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
420c58f1c22SToby Isaac     if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);}
421d9deefdfSMatthew G. Knepley     /* Set labels */
422d9deefdfSMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
423d9deefdfSMatthew G. Knepley       if (out.pointmarkerlist[v]) {
424d988aadeSMatthew G. Knepley         if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);}
425d9deefdfSMatthew G. Knepley       }
426d9deefdfSMatthew G. Knepley     }
427d9deefdfSMatthew G. Knepley     if (interpolate) {
428d9deefdfSMatthew G. Knepley       PetscInt e;
429d9deefdfSMatthew G. Knepley 
430d9deefdfSMatthew G. Knepley       for (e = 0; e < out.numberofedges; e++) {
431d9deefdfSMatthew G. Knepley         if (out.edgemarkerlist[e]) {
432d9deefdfSMatthew G. Knepley           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
433d9deefdfSMatthew G. Knepley           const PetscInt *edges;
434d9deefdfSMatthew G. Knepley           PetscInt        numEdges;
435d9deefdfSMatthew G. Knepley 
436d9deefdfSMatthew G. Knepley           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
437d9deefdfSMatthew G. Knepley           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
438d988aadeSMatthew G. Knepley           if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
439d9deefdfSMatthew G. Knepley           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
440d9deefdfSMatthew G. Knepley         }
441d9deefdfSMatthew G. Knepley       }
442d9deefdfSMatthew G. Knepley     }
443d9deefdfSMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
444d9deefdfSMatthew G. Knepley   }
445d9deefdfSMatthew G. Knepley #if 0 /* Do not currently support holes */
446d9deefdfSMatthew G. Knepley   ierr = DMPlexCopyHoles(*dm, boundary);CHKERRQ(ierr);
447d9deefdfSMatthew G. Knepley #endif
448d9deefdfSMatthew G. Knepley   ierr = FiniOutput_Triangle(&out);CHKERRQ(ierr);
449d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
450d9deefdfSMatthew G. Knepley }
451d9deefdfSMatthew G. Knepley #endif /* PETSC_HAVE_TRIANGLE */
452d9deefdfSMatthew G. Knepley 
453d9deefdfSMatthew G. Knepley #if defined(PETSC_HAVE_TETGEN)
454d9deefdfSMatthew G. Knepley #include <tetgen.h>
455e9db6b00SMatthew G. Knepley static PetscErrorCode DMPlexGenerate_Tetgen(DM boundary, PetscBool interpolate, DM *dm)
456d9deefdfSMatthew G. Knepley {
457d9deefdfSMatthew G. Knepley   MPI_Comm       comm;
4583d8f7108SMatthew G. Knepley   DM_Plex       *mesh      = (DM_Plex *) boundary->data;
459d9deefdfSMatthew G. Knepley   const PetscInt dim       = 3;
460d988aadeSMatthew G. Knepley   const char    *labelName = "marker";
461d9deefdfSMatthew G. Knepley   ::tetgenio     in;
462d9deefdfSMatthew G. Knepley   ::tetgenio     out;
463d988aadeSMatthew G. Knepley   DMLabel        label;
464d9deefdfSMatthew G. Knepley   PetscInt       vStart, vEnd, v, fStart, fEnd, f;
465d9deefdfSMatthew G. Knepley   PetscMPIInt    rank;
466d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
467d9deefdfSMatthew G. Knepley 
468d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
469d9deefdfSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr);
470d9deefdfSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
471d9deefdfSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr);
472c58f1c22SToby Isaac   ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr);
473d988aadeSMatthew G. Knepley 
474d9deefdfSMatthew G. Knepley   in.numberofpoints = vEnd - vStart;
475d9deefdfSMatthew G. Knepley   if (in.numberofpoints > 0) {
476d9deefdfSMatthew G. Knepley     PetscSection coordSection;
477d9deefdfSMatthew G. Knepley     Vec          coordinates;
478d9deefdfSMatthew G. Knepley     PetscScalar *array;
479d9deefdfSMatthew G. Knepley 
480d9deefdfSMatthew G. Knepley     in.pointlist       = new double[in.numberofpoints*dim];
481d9deefdfSMatthew G. Knepley     in.pointmarkerlist = new int[in.numberofpoints];
482d9deefdfSMatthew G. Knepley 
483d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr);
484d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr);
485d9deefdfSMatthew G. Knepley     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
486d9deefdfSMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
487d9deefdfSMatthew G. Knepley       const PetscInt idx = v - vStart;
488d9deefdfSMatthew G. Knepley       PetscInt       off, d;
489d9deefdfSMatthew G. Knepley 
490d9deefdfSMatthew G. Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
491d9deefdfSMatthew G. Knepley       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
492d988aadeSMatthew G. Knepley       if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);}
493d9deefdfSMatthew G. Knepley     }
494d9deefdfSMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
495d9deefdfSMatthew G. Knepley   }
496d9deefdfSMatthew G. Knepley   ierr  = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr);
497d9deefdfSMatthew G. Knepley 
498d9deefdfSMatthew G. Knepley   in.numberoffacets = fEnd - fStart;
499d9deefdfSMatthew G. Knepley   if (in.numberoffacets > 0) {
500d9deefdfSMatthew G. Knepley     in.facetlist       = new tetgenio::facet[in.numberoffacets];
501d9deefdfSMatthew G. Knepley     in.facetmarkerlist = new int[in.numberoffacets];
502d9deefdfSMatthew G. Knepley     for (f = fStart; f < fEnd; ++f) {
503d9deefdfSMatthew G. Knepley       const PetscInt idx     = f - fStart;
504d9deefdfSMatthew G. Knepley       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v;
505d9deefdfSMatthew G. Knepley 
506d9deefdfSMatthew G. Knepley       in.facetlist[idx].numberofpolygons = 1;
507d9deefdfSMatthew G. Knepley       in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
508d9deefdfSMatthew G. Knepley       in.facetlist[idx].numberofholes    = 0;
509d9deefdfSMatthew G. Knepley       in.facetlist[idx].holelist         = NULL;
510d9deefdfSMatthew G. Knepley 
511d9deefdfSMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
512d9deefdfSMatthew G. Knepley       for (p = 0; p < numPoints*2; p += 2) {
513d9deefdfSMatthew G. Knepley         const PetscInt point = points[p];
514d9deefdfSMatthew G. Knepley         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
515d9deefdfSMatthew G. Knepley       }
516d9deefdfSMatthew G. Knepley 
517d9deefdfSMatthew G. Knepley       tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
518d9deefdfSMatthew G. Knepley       poly->numberofvertices = numVertices;
519d9deefdfSMatthew G. Knepley       poly->vertexlist       = new int[poly->numberofvertices];
520d9deefdfSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) {
521d9deefdfSMatthew G. Knepley         const PetscInt vIdx = points[v] - vStart;
522d9deefdfSMatthew G. Knepley         poly->vertexlist[v] = vIdx;
523d9deefdfSMatthew G. Knepley       }
524d988aadeSMatthew G. Knepley       if (label) {ierr = DMLabelGetValue(label, f, &in.facetmarkerlist[idx]);CHKERRQ(ierr);}
525d9deefdfSMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
526d9deefdfSMatthew G. Knepley     }
527d9deefdfSMatthew G. Knepley   }
528d9deefdfSMatthew G. Knepley   if (!rank) {
529d9deefdfSMatthew G. Knepley     char args[32];
530d9deefdfSMatthew G. Knepley 
531d9deefdfSMatthew G. Knepley     /* Take away 'Q' for verbose output */
532d9deefdfSMatthew G. Knepley     ierr = PetscStrcpy(args, "pqezQ");CHKERRQ(ierr);
533d9deefdfSMatthew G. Knepley     if (mesh->tetgenOpts) {::tetrahedralize(mesh->tetgenOpts, &in, &out);}
534d9deefdfSMatthew G. Knepley     else                  {::tetrahedralize(args, &in, &out);}
535d9deefdfSMatthew G. Knepley   }
536d9deefdfSMatthew G. Knepley   {
537d988aadeSMatthew G. Knepley     DMLabel        glabel      = NULL;
538d9deefdfSMatthew G. Knepley     const PetscInt numCorners  = 4;
539d9deefdfSMatthew G. Knepley     const PetscInt numCells    = out.numberoftetrahedra;
540d9deefdfSMatthew G. Knepley     const PetscInt numVertices = out.numberofpoints;
541d9deefdfSMatthew G. Knepley     const double   *meshCoords = out.pointlist;
542d9deefdfSMatthew G. Knepley     int            *cells      = out.tetrahedronlist;
543d9deefdfSMatthew G. Knepley 
544d9deefdfSMatthew G. Knepley     ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr);
545d9deefdfSMatthew G. Knepley     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr);
546c58f1c22SToby Isaac     if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);}
547d9deefdfSMatthew G. Knepley     /* Set labels */
548d9deefdfSMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
549d9deefdfSMatthew G. Knepley       if (out.pointmarkerlist[v]) {
550d988aadeSMatthew G. Knepley         if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);}
551d9deefdfSMatthew G. Knepley       }
552d9deefdfSMatthew G. Knepley     }
553d9deefdfSMatthew G. Knepley     if (interpolate) {
554d9deefdfSMatthew G. Knepley       PetscInt e;
555d9deefdfSMatthew G. Knepley 
556d9deefdfSMatthew G. Knepley       for (e = 0; e < out.numberofedges; e++) {
557d9deefdfSMatthew G. Knepley         if (out.edgemarkerlist[e]) {
558d9deefdfSMatthew G. Knepley           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
559d9deefdfSMatthew G. Knepley           const PetscInt *edges;
560d9deefdfSMatthew G. Knepley           PetscInt        numEdges;
561d9deefdfSMatthew G. Knepley 
562d9deefdfSMatthew G. Knepley           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
563d9deefdfSMatthew G. Knepley           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
564d988aadeSMatthew G. Knepley           if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
565d9deefdfSMatthew G. Knepley           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
566d9deefdfSMatthew G. Knepley         }
567d9deefdfSMatthew G. Knepley       }
568d9deefdfSMatthew G. Knepley       for (f = 0; f < out.numberoftrifaces; f++) {
569d9deefdfSMatthew G. Knepley         if (out.trifacemarkerlist[f]) {
570d9deefdfSMatthew G. Knepley           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
571d9deefdfSMatthew G. Knepley           const PetscInt *faces;
572d9deefdfSMatthew G. Knepley           PetscInt        numFaces;
573d9deefdfSMatthew G. Knepley 
574d9deefdfSMatthew G. Knepley           ierr = DMPlexGetJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
575d9deefdfSMatthew G. Knepley           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
57635b814a0SMatthew G. Knepley           if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);}
577d9deefdfSMatthew G. Knepley           ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
578d9deefdfSMatthew G. Knepley         }
579d9deefdfSMatthew G. Knepley       }
580d9deefdfSMatthew G. Knepley     }
581d9deefdfSMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
582d9deefdfSMatthew G. Knepley   }
583d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
584d9deefdfSMatthew G. Knepley }
585d9deefdfSMatthew G. Knepley 
586e9db6b00SMatthew G. Knepley static PetscErrorCode DMPlexRefine_Tetgen(DM dm, double *maxVolumes, DM *dmRefined)
587d9deefdfSMatthew G. Knepley {
588d9deefdfSMatthew G. Knepley   MPI_Comm       comm;
589d9deefdfSMatthew G. Knepley   const PetscInt dim       = 3;
590d988aadeSMatthew G. Knepley   const char    *labelName = "marker";
591d9deefdfSMatthew G. Knepley   ::tetgenio     in;
592d9deefdfSMatthew G. Knepley   ::tetgenio     out;
593d988aadeSMatthew G. Knepley   DMLabel        label;
594d9deefdfSMatthew G. Knepley   PetscInt       vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal;
595d9deefdfSMatthew G. Knepley   PetscMPIInt    rank;
596d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
597d9deefdfSMatthew G. Knepley 
598d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
599d9deefdfSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
600d9deefdfSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
601d9deefdfSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
602b2566f29SBarry Smith   ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
603d9deefdfSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
604c58f1c22SToby Isaac   ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr);
605d9deefdfSMatthew G. Knepley 
606d9deefdfSMatthew G. Knepley   in.numberofpoints = vEnd - vStart;
607d9deefdfSMatthew G. Knepley   if (in.numberofpoints > 0) {
608d9deefdfSMatthew G. Knepley     PetscSection coordSection;
609d9deefdfSMatthew G. Knepley     Vec          coordinates;
610d9deefdfSMatthew G. Knepley     PetscScalar *array;
611d9deefdfSMatthew G. Knepley 
612d9deefdfSMatthew G. Knepley     in.pointlist       = new double[in.numberofpoints*dim];
613d9deefdfSMatthew G. Knepley     in.pointmarkerlist = new int[in.numberofpoints];
614d9deefdfSMatthew G. Knepley 
615d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
616d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
617d9deefdfSMatthew G. Knepley     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
618d9deefdfSMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
619d9deefdfSMatthew G. Knepley       const PetscInt idx = v - vStart;
620d9deefdfSMatthew G. Knepley       PetscInt       off, d;
621d9deefdfSMatthew G. Knepley 
622d9deefdfSMatthew G. Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
623d9deefdfSMatthew G. Knepley       for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
624d988aadeSMatthew G. Knepley       if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);}
625d9deefdfSMatthew G. Knepley     }
626d9deefdfSMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
627d9deefdfSMatthew G. Knepley   }
628d9deefdfSMatthew G. Knepley   ierr  = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
629d9deefdfSMatthew G. Knepley 
630d9deefdfSMatthew G. Knepley   in.numberofcorners       = 4;
631d9deefdfSMatthew G. Knepley   in.numberoftetrahedra    = cEnd - cStart;
632d9deefdfSMatthew G. Knepley   in.tetrahedronvolumelist = (double*) maxVolumes;
633d9deefdfSMatthew G. Knepley   if (in.numberoftetrahedra > 0) {
634d9deefdfSMatthew G. Knepley     in.tetrahedronlist = new int[in.numberoftetrahedra*in.numberofcorners];
635d9deefdfSMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
636d9deefdfSMatthew G. Knepley       const PetscInt idx      = c - cStart;
637d9deefdfSMatthew G. Knepley       PetscInt      *closure = NULL;
638d9deefdfSMatthew G. Knepley       PetscInt       closureSize;
639d9deefdfSMatthew G. Knepley 
640d9deefdfSMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
641d9deefdfSMatthew G. Knepley       if ((closureSize != 5) && (closureSize != 15)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %D vertices in closure", closureSize);
642d9deefdfSMatthew G. Knepley       for (v = 0; v < 4; ++v) {
643d9deefdfSMatthew G. Knepley         in.tetrahedronlist[idx*in.numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
644d9deefdfSMatthew G. Knepley       }
645d9deefdfSMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
646d9deefdfSMatthew G. Knepley     }
647d9deefdfSMatthew G. Knepley   }
648d9deefdfSMatthew G. Knepley   /* TODO: Put in boundary faces with markers */
649d9deefdfSMatthew G. Knepley   if (!rank) {
650d9deefdfSMatthew G. Knepley     char args[32];
651d9deefdfSMatthew G. Knepley 
652d9deefdfSMatthew G. Knepley     /* Take away 'Q' for verbose output */
653d9deefdfSMatthew G. Knepley     /*ierr = PetscStrcpy(args, "qezQra");CHKERRQ(ierr); */
654d9deefdfSMatthew G. Knepley     ierr = PetscStrcpy(args, "qezraVVVV");CHKERRQ(ierr);
655d9deefdfSMatthew G. Knepley     ::tetrahedralize(args, &in, &out);
656d9deefdfSMatthew G. Knepley   }
657d9deefdfSMatthew G. Knepley   in.tetrahedronvolumelist = NULL;
658d9deefdfSMatthew G. Knepley 
659d9deefdfSMatthew G. Knepley   {
660d988aadeSMatthew G. Knepley     DMLabel        rlabel      = NULL;
661d9deefdfSMatthew G. Knepley     const PetscInt numCorners  = 4;
662d9deefdfSMatthew G. Knepley     const PetscInt numCells    = out.numberoftetrahedra;
663d9deefdfSMatthew G. Knepley     const PetscInt numVertices = out.numberofpoints;
664d9deefdfSMatthew G. Knepley     const double   *meshCoords = out.pointlist;
665d9deefdfSMatthew G. Knepley     int            *cells      = out.tetrahedronlist;
666d9deefdfSMatthew G. Knepley 
667d9deefdfSMatthew G. Knepley     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
668d9deefdfSMatthew G. Knepley 
669d9deefdfSMatthew G. Knepley     ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr);
670d9deefdfSMatthew G. Knepley     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
671c58f1c22SToby Isaac     if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);}
672d9deefdfSMatthew G. Knepley     /* Set labels */
673d9deefdfSMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
674d9deefdfSMatthew G. Knepley       if (out.pointmarkerlist[v]) {
675d988aadeSMatthew G. Knepley         if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);}
676d9deefdfSMatthew G. Knepley       }
677d9deefdfSMatthew G. Knepley     }
678d9deefdfSMatthew G. Knepley     if (interpolate) {
679d9deefdfSMatthew G. Knepley       PetscInt e, f;
680d9deefdfSMatthew G. Knepley 
681d9deefdfSMatthew G. Knepley       for (e = 0; e < out.numberofedges; e++) {
682d9deefdfSMatthew G. Knepley         if (out.edgemarkerlist[e]) {
683d9deefdfSMatthew G. Knepley           const PetscInt  vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells};
684d9deefdfSMatthew G. Knepley           const PetscInt *edges;
685d9deefdfSMatthew G. Knepley           PetscInt        numEdges;
686d9deefdfSMatthew G. Knepley 
687d9deefdfSMatthew G. Knepley           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
688d9deefdfSMatthew G. Knepley           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
689d988aadeSMatthew G. Knepley           if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);}
690d9deefdfSMatthew G. Knepley           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
691d9deefdfSMatthew G. Knepley         }
692d9deefdfSMatthew G. Knepley       }
693d9deefdfSMatthew G. Knepley       for (f = 0; f < out.numberoftrifaces; f++) {
694d9deefdfSMatthew G. Knepley         if (out.trifacemarkerlist[f]) {
695d9deefdfSMatthew G. Knepley           const PetscInt  vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells};
696d9deefdfSMatthew G. Knepley           const PetscInt *faces;
697d9deefdfSMatthew G. Knepley           PetscInt        numFaces;
698d9deefdfSMatthew G. Knepley 
699d9deefdfSMatthew G. Knepley           ierr = DMPlexGetJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
700d9deefdfSMatthew G. Knepley           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
701d988aadeSMatthew G. Knepley           if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);}
702d9deefdfSMatthew G. Knepley           ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
703d9deefdfSMatthew G. Knepley         }
704d9deefdfSMatthew G. Knepley       }
705d9deefdfSMatthew G. Knepley     }
706d9deefdfSMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
707d9deefdfSMatthew G. Knepley   }
708d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
709d9deefdfSMatthew G. Knepley }
710d9deefdfSMatthew G. Knepley #endif /* PETSC_HAVE_TETGEN */
711d9deefdfSMatthew G. Knepley 
712d9deefdfSMatthew G. Knepley #if defined(PETSC_HAVE_CTETGEN)
713d9deefdfSMatthew G. Knepley #include <ctetgen.h>
714d9deefdfSMatthew G. Knepley 
715e9db6b00SMatthew G. Knepley static PetscErrorCode DMPlexGenerate_CTetgen(DM boundary, PetscBool interpolate, DM *dm)
716d9deefdfSMatthew G. Knepley {
717d9deefdfSMatthew G. Knepley   MPI_Comm       comm;
718d9deefdfSMatthew G. Knepley   const PetscInt dim       = 3;
719d988aadeSMatthew G. Knepley   const char    *labelName = "marker";
720d9deefdfSMatthew G. Knepley   PLC           *in, *out;
721d988aadeSMatthew G. Knepley   DMLabel        label;
722d9deefdfSMatthew G. Knepley   PetscInt       verbose = 0, vStart, vEnd, v, fStart, fEnd, f;
723d9deefdfSMatthew G. Knepley   PetscMPIInt    rank;
724d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
725d9deefdfSMatthew G. Knepley 
726d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
727d9deefdfSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr);
728c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,((PetscObject) boundary)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr);
729d9deefdfSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
730d9deefdfSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr);
731c58f1c22SToby Isaac   ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr);
732d9deefdfSMatthew G. Knepley   ierr = PLCCreate(&in);CHKERRQ(ierr);
733d9deefdfSMatthew G. Knepley   ierr = PLCCreate(&out);CHKERRQ(ierr);
734d9deefdfSMatthew G. Knepley 
735d9deefdfSMatthew G. Knepley   in->numberofpoints = vEnd - vStart;
736d9deefdfSMatthew G. Knepley   if (in->numberofpoints > 0) {
737d9deefdfSMatthew G. Knepley     PetscSection coordSection;
738d9deefdfSMatthew G. Knepley     Vec          coordinates;
739d9deefdfSMatthew G. Knepley     PetscScalar *array;
740d9deefdfSMatthew G. Knepley 
741d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr);
742d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in->numberofpoints,       &in->pointmarkerlist);CHKERRQ(ierr);
743d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr);
744d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr);
745d9deefdfSMatthew G. Knepley     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
746d9deefdfSMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
747d9deefdfSMatthew G. Knepley       const PetscInt idx = v - vStart;
748a9a51d2cSMatthew G. Knepley       PetscInt       off, d, m = -1;
749d9deefdfSMatthew G. Knepley 
750d9deefdfSMatthew G. Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
751d9deefdfSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
752d9deefdfSMatthew G. Knepley         in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
753d9deefdfSMatthew G. Knepley       }
754d988aadeSMatthew G. Knepley       if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);}
755d9deefdfSMatthew G. Knepley 
756d9deefdfSMatthew G. Knepley       in->pointmarkerlist[idx] = (int) m;
757d9deefdfSMatthew G. Knepley     }
758d9deefdfSMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
759d9deefdfSMatthew G. Knepley   }
760d9deefdfSMatthew G. Knepley   ierr  = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr);
761d9deefdfSMatthew G. Knepley 
762d9deefdfSMatthew G. Knepley   in->numberoffacets = fEnd - fStart;
763d9deefdfSMatthew G. Knepley   if (in->numberoffacets > 0) {
764d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in->numberoffacets, &in->facetlist);CHKERRQ(ierr);
765d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in->numberoffacets,   &in->facetmarkerlist);CHKERRQ(ierr);
766d9deefdfSMatthew G. Knepley     for (f = fStart; f < fEnd; ++f) {
767d9deefdfSMatthew G. Knepley       const PetscInt idx     = f - fStart;
768a9a51d2cSMatthew G. Knepley       PetscInt      *points = NULL, numPoints, p, numVertices = 0, v, m = -1;
769d9deefdfSMatthew G. Knepley       polygon       *poly;
770d9deefdfSMatthew G. Knepley 
771d9deefdfSMatthew G. Knepley       in->facetlist[idx].numberofpolygons = 1;
772d9deefdfSMatthew G. Knepley 
773d9deefdfSMatthew G. Knepley       ierr = PetscMalloc1(in->facetlist[idx].numberofpolygons, &in->facetlist[idx].polygonlist);CHKERRQ(ierr);
774d9deefdfSMatthew G. Knepley 
775d9deefdfSMatthew G. Knepley       in->facetlist[idx].numberofholes    = 0;
776d9deefdfSMatthew G. Knepley       in->facetlist[idx].holelist         = NULL;
777d9deefdfSMatthew G. Knepley 
778d9deefdfSMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
779d9deefdfSMatthew G. Knepley       for (p = 0; p < numPoints*2; p += 2) {
780d9deefdfSMatthew G. Knepley         const PetscInt point = points[p];
781d9deefdfSMatthew G. Knepley         if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
782d9deefdfSMatthew G. Knepley       }
783d9deefdfSMatthew G. Knepley 
784d9deefdfSMatthew G. Knepley       poly                   = in->facetlist[idx].polygonlist;
785d9deefdfSMatthew G. Knepley       poly->numberofvertices = numVertices;
786d9deefdfSMatthew G. Knepley       ierr                   = PetscMalloc1(poly->numberofvertices, &poly->vertexlist);CHKERRQ(ierr);
787d9deefdfSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) {
788d9deefdfSMatthew G. Knepley         const PetscInt vIdx = points[v] - vStart;
789d9deefdfSMatthew G. Knepley         poly->vertexlist[v] = vIdx;
790d9deefdfSMatthew G. Knepley       }
791d988aadeSMatthew G. Knepley       if (label) {ierr = DMLabelGetValue(label, f, &m);CHKERRQ(ierr);}
792d9deefdfSMatthew G. Knepley       in->facetmarkerlist[idx] = (int) m;
793d9deefdfSMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
794d9deefdfSMatthew G. Knepley     }
795d9deefdfSMatthew G. Knepley   }
796d9deefdfSMatthew G. Knepley   if (!rank) {
797d9deefdfSMatthew G. Knepley     TetGenOpts t;
798d9deefdfSMatthew G. Knepley 
799d9deefdfSMatthew G. Knepley     ierr        = TetGenOptsInitialize(&t);CHKERRQ(ierr);
800d9deefdfSMatthew G. Knepley     t.in        = boundary; /* Should go away */
801d9deefdfSMatthew G. Knepley     t.plc       = 1;
802d9deefdfSMatthew G. Knepley     t.quality   = 1;
803d9deefdfSMatthew G. Knepley     t.edgesout  = 1;
804d9deefdfSMatthew G. Knepley     t.zeroindex = 1;
805d9deefdfSMatthew G. Knepley     t.quiet     = 1;
806d9deefdfSMatthew G. Knepley     t.verbose   = verbose;
807d9deefdfSMatthew G. Knepley     ierr        = TetGenCheckOpts(&t);CHKERRQ(ierr);
808d9deefdfSMatthew G. Knepley     ierr        = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr);
809d9deefdfSMatthew G. Knepley   }
810d9deefdfSMatthew G. Knepley   {
811d988aadeSMatthew G. Knepley     DMLabel        glabel      = NULL;
812d9deefdfSMatthew G. Knepley     const PetscInt numCorners  = 4;
813d9deefdfSMatthew G. Knepley     const PetscInt numCells    = out->numberoftetrahedra;
814d9deefdfSMatthew G. Knepley     const PetscInt numVertices = out->numberofpoints;
815e9ccc142SToby Isaac     double         *meshCoords;
816d9deefdfSMatthew G. Knepley     int            *cells      = out->tetrahedronlist;
817d9deefdfSMatthew G. Knepley 
818e9ccc142SToby Isaac     if (sizeof (PetscReal) == sizeof (double)) {
819e9ccc142SToby Isaac       meshCoords = (double *) out->pointlist;
820e9ccc142SToby Isaac     }
821e9ccc142SToby Isaac     else {
822e9ccc142SToby Isaac       PetscInt i;
823e9ccc142SToby Isaac 
824e9ccc142SToby Isaac       ierr = PetscMalloc1(3 * numVertices,&meshCoords);CHKERRQ(ierr);
825e9ccc142SToby Isaac       for (i = 0; i < 3 * numVertices; i++) {
826e9ccc142SToby Isaac         meshCoords[i] = (PetscReal) out->pointlist[i];
827e9ccc142SToby Isaac       }
828e9ccc142SToby Isaac     }
829e9ccc142SToby Isaac 
830d9deefdfSMatthew G. Knepley     ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr);
831d9deefdfSMatthew G. Knepley     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr);
832e9ccc142SToby Isaac     if (sizeof (PetscReal) != sizeof (double)) {
833e9ccc142SToby Isaac       ierr = PetscFree(meshCoords);CHKERRQ(ierr);
834e9ccc142SToby Isaac     }
835c58f1c22SToby Isaac     if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);}
836d9deefdfSMatthew G. Knepley     /* Set labels */
837d9deefdfSMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
838d9deefdfSMatthew G. Knepley       if (out->pointmarkerlist[v]) {
839d988aadeSMatthew G. Knepley         if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);}
840d9deefdfSMatthew G. Knepley       }
841d9deefdfSMatthew G. Knepley     }
842d9deefdfSMatthew G. Knepley     if (interpolate) {
843d9deefdfSMatthew G. Knepley       PetscInt e;
844d9deefdfSMatthew G. Knepley 
845d9deefdfSMatthew G. Knepley       for (e = 0; e < out->numberofedges; e++) {
846d9deefdfSMatthew G. Knepley         if (out->edgemarkerlist[e]) {
847d9deefdfSMatthew G. Knepley           const PetscInt  vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells};
848d9deefdfSMatthew G. Knepley           const PetscInt *edges;
849d9deefdfSMatthew G. Knepley           PetscInt        numEdges;
850d9deefdfSMatthew G. Knepley 
851d9deefdfSMatthew G. Knepley           ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
852d9deefdfSMatthew G. Knepley           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
853d988aadeSMatthew G. Knepley           if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);}
854d9deefdfSMatthew G. Knepley           ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
855d9deefdfSMatthew G. Knepley         }
856d9deefdfSMatthew G. Knepley       }
857d9deefdfSMatthew G. Knepley       for (f = 0; f < out->numberoftrifaces; f++) {
858d9deefdfSMatthew G. Knepley         if (out->trifacemarkerlist[f]) {
859d9deefdfSMatthew G. Knepley           const PetscInt  vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells};
860d9deefdfSMatthew G. Knepley           const PetscInt *faces;
861d9deefdfSMatthew G. Knepley           PetscInt        numFaces;
862d9deefdfSMatthew G. Knepley 
863d9deefdfSMatthew G. Knepley           ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
864d9deefdfSMatthew G. Knepley           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
865d988aadeSMatthew G. Knepley           if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);}
866d9deefdfSMatthew G. Knepley           ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
867d9deefdfSMatthew G. Knepley         }
868d9deefdfSMatthew G. Knepley       }
869d9deefdfSMatthew G. Knepley     }
870d9deefdfSMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
871d9deefdfSMatthew G. Knepley   }
872d9deefdfSMatthew G. Knepley 
873d9deefdfSMatthew G. Knepley   ierr = PLCDestroy(&in);CHKERRQ(ierr);
874d9deefdfSMatthew G. Knepley   ierr = PLCDestroy(&out);CHKERRQ(ierr);
875d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
876d9deefdfSMatthew G. Knepley }
877d9deefdfSMatthew G. Knepley 
878e9db6b00SMatthew G. Knepley static PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined)
879d9deefdfSMatthew G. Knepley {
880d9deefdfSMatthew G. Knepley   MPI_Comm       comm;
881d9deefdfSMatthew G. Knepley   const PetscInt dim       = 3;
882d988aadeSMatthew G. Knepley   const char    *labelName = "marker";
883d9deefdfSMatthew G. Knepley   PLC           *in, *out;
884d988aadeSMatthew G. Knepley   DMLabel        label;
885d9deefdfSMatthew G. Knepley   PetscInt       verbose = 0, vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal;
886d9deefdfSMatthew G. Knepley   PetscMPIInt    rank;
887d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
888d9deefdfSMatthew G. Knepley 
889d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
890d9deefdfSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
891c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,((PetscObject) dm)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr);
892d9deefdfSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
893d9deefdfSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
894b2566f29SBarry Smith   ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
895d9deefdfSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
896c58f1c22SToby Isaac   ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr);
897d9deefdfSMatthew G. Knepley   ierr = PLCCreate(&in);CHKERRQ(ierr);
898d9deefdfSMatthew G. Knepley   ierr = PLCCreate(&out);CHKERRQ(ierr);
899d9deefdfSMatthew G. Knepley 
900d9deefdfSMatthew G. Knepley   in->numberofpoints = vEnd - vStart;
901d9deefdfSMatthew G. Knepley   if (in->numberofpoints > 0) {
902d9deefdfSMatthew G. Knepley     PetscSection coordSection;
903d9deefdfSMatthew G. Knepley     Vec          coordinates;
904d9deefdfSMatthew G. Knepley     PetscScalar *array;
905d9deefdfSMatthew G. Knepley 
906d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr);
907d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in->numberofpoints,       &in->pointmarkerlist);CHKERRQ(ierr);
908d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
909d9deefdfSMatthew G. Knepley     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
910d9deefdfSMatthew G. Knepley     ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr);
911d9deefdfSMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
912d9deefdfSMatthew G. Knepley       const PetscInt idx = v - vStart;
913a9a51d2cSMatthew G. Knepley       PetscInt       off, d, m = -1;
914d9deefdfSMatthew G. Knepley 
915d9deefdfSMatthew G. Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
916d9deefdfSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
917d9deefdfSMatthew G. Knepley         in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]);
918d9deefdfSMatthew G. Knepley       }
919d988aadeSMatthew G. Knepley       if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);}
920d9deefdfSMatthew G. Knepley 
921d9deefdfSMatthew G. Knepley       in->pointmarkerlist[idx] = (int) m;
922d9deefdfSMatthew G. Knepley     }
923d9deefdfSMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr);
924d9deefdfSMatthew G. Knepley   }
925d9deefdfSMatthew G. Knepley   ierr  = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
926d9deefdfSMatthew G. Knepley 
927d9deefdfSMatthew G. Knepley   in->numberofcorners       = 4;
928d9deefdfSMatthew G. Knepley   in->numberoftetrahedra    = cEnd - cStart;
929d9deefdfSMatthew G. Knepley   in->tetrahedronvolumelist = maxVolumes;
930d9deefdfSMatthew G. Knepley   if (in->numberoftetrahedra > 0) {
931d9deefdfSMatthew G. Knepley     ierr = PetscMalloc1(in->numberoftetrahedra*in->numberofcorners, &in->tetrahedronlist);CHKERRQ(ierr);
932d9deefdfSMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
933d9deefdfSMatthew G. Knepley       const PetscInt idx      = c - cStart;
934d9deefdfSMatthew G. Knepley       PetscInt      *closure = NULL;
935d9deefdfSMatthew G. Knepley       PetscInt       closureSize;
936d9deefdfSMatthew G. Knepley 
937d9deefdfSMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
938d9deefdfSMatthew G. Knepley       if ((closureSize != 5) && (closureSize != 15)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %D vertices in closure", closureSize);
939d9deefdfSMatthew G. Knepley       for (v = 0; v < 4; ++v) {
940d9deefdfSMatthew G. Knepley         in->tetrahedronlist[idx*in->numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart;
941d9deefdfSMatthew G. Knepley       }
942d9deefdfSMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
943d9deefdfSMatthew G. Knepley     }
944d9deefdfSMatthew G. Knepley   }
945d9deefdfSMatthew G. Knepley   if (!rank) {
946d9deefdfSMatthew G. Knepley     TetGenOpts t;
947d9deefdfSMatthew G. Knepley 
948d9deefdfSMatthew G. Knepley     ierr = TetGenOptsInitialize(&t);CHKERRQ(ierr);
949d9deefdfSMatthew G. Knepley 
950d9deefdfSMatthew G. Knepley     t.in        = dm; /* Should go away */
951d9deefdfSMatthew G. Knepley     t.refine    = 1;
952d9deefdfSMatthew G. Knepley     t.varvolume = 1;
953d9deefdfSMatthew G. Knepley     t.quality   = 1;
954d9deefdfSMatthew G. Knepley     t.edgesout  = 1;
955d9deefdfSMatthew G. Knepley     t.zeroindex = 1;
956d9deefdfSMatthew G. Knepley     t.quiet     = 1;
957d9deefdfSMatthew G. Knepley     t.verbose   = verbose; /* Change this */
958d9deefdfSMatthew G. Knepley 
959d9deefdfSMatthew G. Knepley     ierr = TetGenCheckOpts(&t);CHKERRQ(ierr);
960d9deefdfSMatthew G. Knepley     ierr = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr);
961d9deefdfSMatthew G. Knepley   }
962d9deefdfSMatthew G. Knepley   {
963d988aadeSMatthew G. Knepley     DMLabel        rlabel      = NULL;
964d9deefdfSMatthew G. Knepley     const PetscInt numCorners  = 4;
965d9deefdfSMatthew G. Knepley     const PetscInt numCells    = out->numberoftetrahedra;
966d9deefdfSMatthew G. Knepley     const PetscInt numVertices = out->numberofpoints;
967e9ccc142SToby Isaac     double         *meshCoords;
968d9deefdfSMatthew G. Knepley     int            *cells      = out->tetrahedronlist;
969d9deefdfSMatthew G. Knepley     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
970d9deefdfSMatthew G. Knepley 
971e9ccc142SToby Isaac     if (sizeof (PetscReal) == sizeof (double)) {
972e9ccc142SToby Isaac       meshCoords = (double *) out->pointlist;
973e9ccc142SToby Isaac     }
974e9ccc142SToby Isaac     else {
975e9ccc142SToby Isaac       PetscInt i;
976e9ccc142SToby Isaac 
977e9ccc142SToby Isaac       ierr = PetscMalloc1(3 * numVertices,&meshCoords);CHKERRQ(ierr);
978e9ccc142SToby Isaac       for (i = 0; i < 3 * numVertices; i++) {
979e9ccc142SToby Isaac         meshCoords[i] = (PetscReal) out->pointlist[i];
980e9ccc142SToby Isaac       }
981e9ccc142SToby Isaac     }
982e9ccc142SToby Isaac 
983d9deefdfSMatthew G. Knepley     ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr);
984d9deefdfSMatthew G. Knepley     ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr);
985e9ccc142SToby Isaac     if (sizeof (PetscReal) != sizeof (double)) {
986e9ccc142SToby Isaac       ierr = PetscFree(meshCoords);CHKERRQ(ierr);
987e9ccc142SToby Isaac     }
988c58f1c22SToby Isaac     if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);}
989d9deefdfSMatthew G. Knepley     /* Set labels */
990d9deefdfSMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
991d9deefdfSMatthew G. Knepley       if (out->pointmarkerlist[v]) {
992d988aadeSMatthew G. Knepley         if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);}
993d9deefdfSMatthew G. Knepley       }
994d9deefdfSMatthew G. Knepley     }
995d9deefdfSMatthew G. Knepley     if (interpolate) {
996d9deefdfSMatthew G. Knepley       PetscInt e, f;
997d9deefdfSMatthew G. Knepley 
998d9deefdfSMatthew G. Knepley       for (e = 0; e < out->numberofedges; e++) {
999d9deefdfSMatthew G. Knepley         if (out->edgemarkerlist[e]) {
1000d9deefdfSMatthew G. Knepley           const PetscInt  vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells};
1001d9deefdfSMatthew G. Knepley           const PetscInt *edges;
1002d9deefdfSMatthew G. Knepley           PetscInt        numEdges;
1003d9deefdfSMatthew G. Knepley 
1004d9deefdfSMatthew G. Knepley           ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1005d9deefdfSMatthew G. Knepley           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges);
1006d988aadeSMatthew G. Knepley           if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);}
1007d9deefdfSMatthew G. Knepley           ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr);
1008d9deefdfSMatthew G. Knepley         }
1009d9deefdfSMatthew G. Knepley       }
1010d9deefdfSMatthew G. Knepley       for (f = 0; f < out->numberoftrifaces; f++) {
1011d9deefdfSMatthew G. Knepley         if (out->trifacemarkerlist[f]) {
1012d9deefdfSMatthew G. Knepley           const PetscInt  vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells};
1013d9deefdfSMatthew G. Knepley           const PetscInt *faces;
1014d9deefdfSMatthew G. Knepley           PetscInt        numFaces;
1015d9deefdfSMatthew G. Knepley 
1016d9deefdfSMatthew G. Knepley           ierr = DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1017d9deefdfSMatthew G. Knepley           if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces);
1018d988aadeSMatthew G. Knepley           if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);}
1019d9deefdfSMatthew G. Knepley           ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr);
1020d9deefdfSMatthew G. Knepley         }
1021d9deefdfSMatthew G. Knepley       }
1022d9deefdfSMatthew G. Knepley     }
1023d9deefdfSMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr);
1024d9deefdfSMatthew G. Knepley   }
1025d9deefdfSMatthew G. Knepley   ierr = PLCDestroy(&in);CHKERRQ(ierr);
1026d9deefdfSMatthew G. Knepley   ierr = PLCDestroy(&out);CHKERRQ(ierr);
1027d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
1028d9deefdfSMatthew G. Knepley }
1029d9deefdfSMatthew G. Knepley #endif /* PETSC_HAVE_CTETGEN */
1030d9deefdfSMatthew G. Knepley 
1031d9deefdfSMatthew G. Knepley /*@C
1032d9deefdfSMatthew G. Knepley   DMPlexGenerate - Generates a mesh.
1033d9deefdfSMatthew G. Knepley 
1034d9deefdfSMatthew G. Knepley   Not Collective
1035d9deefdfSMatthew G. Knepley 
1036d9deefdfSMatthew G. Knepley   Input Parameters:
1037d9deefdfSMatthew G. Knepley + boundary - The DMPlex boundary object
1038d9deefdfSMatthew G. Knepley . name - The mesh generation package name
1039d9deefdfSMatthew G. Knepley - interpolate - Flag to create intermediate mesh elements
1040d9deefdfSMatthew G. Knepley 
1041d9deefdfSMatthew G. Knepley   Output Parameter:
1042d9deefdfSMatthew G. Knepley . mesh - The DMPlex object
1043d9deefdfSMatthew G. Knepley 
1044d9deefdfSMatthew G. Knepley   Level: intermediate
1045d9deefdfSMatthew G. Knepley 
1046d9deefdfSMatthew G. Knepley .keywords: mesh, elements
1047d9deefdfSMatthew G. Knepley .seealso: DMPlexCreate(), DMRefine()
1048d9deefdfSMatthew G. Knepley @*/
1049d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
1050d9deefdfSMatthew G. Knepley {
1051d9deefdfSMatthew G. Knepley   PetscInt       dim;
1052d9deefdfSMatthew G. Knepley   char           genname[1024];
1053d9deefdfSMatthew G. Knepley   PetscBool      isTriangle = PETSC_FALSE, isTetgen = PETSC_FALSE, isCTetgen = PETSC_FALSE, flg;
1054d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
1055d9deefdfSMatthew G. Knepley 
1056d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
1057d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(boundary, DM_CLASSID, 1);
1058d9deefdfSMatthew G. Knepley   PetscValidLogicalCollectiveBool(boundary, interpolate, 2);
1059d9deefdfSMatthew G. Knepley   ierr = DMGetDimension(boundary, &dim);CHKERRQ(ierr);
1060c5929fdfSBarry Smith   ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr);
1061d9deefdfSMatthew G. Knepley   if (flg) name = genname;
1062d9deefdfSMatthew G. Knepley   if (name) {
1063d9deefdfSMatthew G. Knepley     ierr = PetscStrcmp(name, "triangle", &isTriangle);CHKERRQ(ierr);
1064d9deefdfSMatthew G. Knepley     ierr = PetscStrcmp(name, "tetgen",   &isTetgen);CHKERRQ(ierr);
1065d9deefdfSMatthew G. Knepley     ierr = PetscStrcmp(name, "ctetgen",  &isCTetgen);CHKERRQ(ierr);
1066d9deefdfSMatthew G. Knepley   }
1067d9deefdfSMatthew G. Knepley   switch (dim) {
1068d9deefdfSMatthew G. Knepley   case 1:
1069d9deefdfSMatthew G. Knepley     if (!name || isTriangle) {
1070d9deefdfSMatthew G. Knepley #if defined(PETSC_HAVE_TRIANGLE)
1071d9deefdfSMatthew G. Knepley       ierr = DMPlexGenerate_Triangle(boundary, interpolate, mesh);CHKERRQ(ierr);
1072d9deefdfSMatthew G. Knepley #else
1073d9deefdfSMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Mesh generation needs external package support.\nPlease reconfigure with --download-triangle.");
1074d9deefdfSMatthew G. Knepley #endif
1075d9deefdfSMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Unknown 2D mesh generation package %s", name);
1076d9deefdfSMatthew G. Knepley     break;
1077d9deefdfSMatthew G. Knepley   case 2:
1078d9deefdfSMatthew G. Knepley     if (!name || isCTetgen) {
1079d9deefdfSMatthew G. Knepley #if defined(PETSC_HAVE_CTETGEN)
1080d9deefdfSMatthew G. Knepley       ierr = DMPlexGenerate_CTetgen(boundary, interpolate, mesh);CHKERRQ(ierr);
1081d9deefdfSMatthew G. Knepley #else
1082d9deefdfSMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "CTetgen needs external package support.\nPlease reconfigure with --download-ctetgen.");
1083d9deefdfSMatthew G. Knepley #endif
1084d9deefdfSMatthew G. Knepley     } else if (isTetgen) {
1085d9deefdfSMatthew G. Knepley #if defined(PETSC_HAVE_TETGEN)
1086d9deefdfSMatthew G. Knepley       ierr = DMPlexGenerate_Tetgen(boundary, interpolate, mesh);CHKERRQ(ierr);
1087d9deefdfSMatthew G. Knepley #else
10882ca110d9SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Tetgen needs external package support.\nPlease reconfigure with --download-tetgen.");
1089d9deefdfSMatthew G. Knepley #endif
1090d9deefdfSMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Unknown 3D mesh generation package %s", name);
1091d9deefdfSMatthew G. Knepley     break;
1092d9deefdfSMatthew G. Knepley   default:
1093d9deefdfSMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Mesh generation for a dimension %d boundary is not supported.", dim);
1094d9deefdfSMatthew G. Knepley   }
1095d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
1096d9deefdfSMatthew G. Knepley }
1097d9deefdfSMatthew G. Knepley 
1098755f5a09SToby Isaac #if defined(PETSC_HAVE_TRIANGLE) || defined(PETSC_HAVE_CTETGEN) || defined(PETSC_HAVE_TETGEN)
1099713918a9SToby Isaac static PetscErrorCode DMRefine_Plex_Label(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal maxVolumes[])
1100713918a9SToby Isaac {
1101713918a9SToby Isaac   PetscInt       dim, c;
1102713918a9SToby Isaac   PetscReal      refRatio;
1103713918a9SToby Isaac   PetscErrorCode ierr;
1104713918a9SToby Isaac 
1105713918a9SToby Isaac   PetscFunctionBegin;
1106713918a9SToby Isaac   ierr     = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1107713918a9SToby Isaac   refRatio = (PetscReal) ((PetscInt) 1 << dim);
1108713918a9SToby Isaac   for (c = cStart; c < cEnd; c++) {
1109713918a9SToby Isaac     PetscReal vol;
1110713918a9SToby Isaac     PetscInt  i, closureSize = 0;
1111713918a9SToby Isaac     PetscInt  *closure = NULL;
1112713918a9SToby Isaac     PetscBool anyRefine  = PETSC_FALSE;
1113713918a9SToby Isaac     PetscBool anyCoarsen = PETSC_FALSE;
1114cd3c525cSToby Isaac     PetscBool anyKeep    = PETSC_FALSE;
1115713918a9SToby Isaac 
1116713918a9SToby Isaac     ierr = DMPlexComputeCellGeometryFVM(dm,c,&vol,NULL,NULL);CHKERRQ(ierr);
1117cd3c525cSToby Isaac     maxVolumes[c - cStart] = vol;
1118a19b9e93SToby Isaac     ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1119713918a9SToby Isaac     for (i = 0; i < closureSize; i++) {
1120713918a9SToby Isaac       PetscInt point = closure[2 * i], refFlag;
1121713918a9SToby Isaac 
1122713918a9SToby Isaac       ierr = DMLabelGetValue(adaptLabel,point,&refFlag);CHKERRQ(ierr);
1123713918a9SToby Isaac       switch (refFlag) {
1124713918a9SToby Isaac       case DM_ADAPT_REFINE:
1125713918a9SToby Isaac         anyRefine = PETSC_TRUE;
1126713918a9SToby Isaac         break;
1127713918a9SToby Isaac       case DM_ADAPT_COARSEN:
1128713918a9SToby Isaac         anyCoarsen = PETSC_TRUE;
1129713918a9SToby Isaac         break;
1130713918a9SToby Isaac       case DM_ADAPT_KEEP:
1131cd3c525cSToby Isaac         anyKeep = PETSC_TRUE;
1132cd3c525cSToby Isaac         break;
1133cd3c525cSToby Isaac       case DM_ADAPT_DETERMINE:
1134713918a9SToby Isaac         break;
1135713918a9SToby Isaac       default:
1136713918a9SToby Isaac         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n",refFlag);
1137713918a9SToby Isaac         break;
1138713918a9SToby Isaac       }
1139713918a9SToby Isaac       if (anyRefine) break;
1140713918a9SToby Isaac     }
1141a19b9e93SToby Isaac     ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1142713918a9SToby Isaac     if (anyRefine) {
1143713918a9SToby Isaac       maxVolumes[c - cStart] = vol / refRatio;
1144cd3c525cSToby Isaac     } else if (anyKeep) {
1145cd3c525cSToby Isaac       maxVolumes[c - cStart] = vol;
1146713918a9SToby Isaac     } else if (anyCoarsen) {
1147713918a9SToby Isaac       maxVolumes[c - cStart] = vol * refRatio;
1148713918a9SToby Isaac     }
1149713918a9SToby Isaac   }
1150713918a9SToby Isaac   PetscFunctionReturn(0);
1151713918a9SToby Isaac }
1152755f5a09SToby Isaac #endif
1153713918a9SToby Isaac 
1154713918a9SToby Isaac static PetscErrorCode DMRefine_Plex_Internal(DM dm, MPI_Comm comm, DMLabel adaptLabel, DM *dmRefined)
1155d9deefdfSMatthew G. Knepley {
1156b28003e6SMatthew G. Knepley   PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *);
1157d9deefdfSMatthew G. Knepley   PetscReal        refinementLimit;
1158d9deefdfSMatthew G. Knepley   PetscInt         dim, cStart, cEnd;
1159d9deefdfSMatthew G. Knepley   char             genname[1024], *name = NULL;
116036447a5eSToby Isaac   PetscBool        isUniform, isTriangle = PETSC_FALSE, isTetgen = PETSC_FALSE, isCTetgen = PETSC_FALSE, flg, localized;
1161d9deefdfSMatthew G. Knepley   PetscErrorCode   ierr;
1162d9deefdfSMatthew G. Knepley 
1163d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
116436447a5eSToby Isaac   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1165d9deefdfSMatthew G. Knepley   ierr = DMPlexGetRefinementUniform(dm, &isUniform);CHKERRQ(ierr);
1166d9deefdfSMatthew G. Knepley   if (isUniform) {
1167d9deefdfSMatthew G. Knepley     CellRefiner cellRefiner;
1168d9deefdfSMatthew G. Knepley 
1169d9deefdfSMatthew G. Knepley     ierr = DMPlexGetCellRefiner_Internal(dm, &cellRefiner);CHKERRQ(ierr);
1170d9deefdfSMatthew G. Knepley     ierr = DMPlexRefineUniform_Internal(dm, cellRefiner, dmRefined);CHKERRQ(ierr);
1171a6ba4734SToby Isaac     ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr);
117236447a5eSToby Isaac     if (localized) {
117336447a5eSToby Isaac       ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);
117436447a5eSToby Isaac     }
1175d9deefdfSMatthew G. Knepley     PetscFunctionReturn(0);
1176d9deefdfSMatthew G. Knepley   }
1177d9deefdfSMatthew G. Knepley   ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr);
1178b28003e6SMatthew G. Knepley   ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr);
1179b28003e6SMatthew G. Knepley   if (refinementLimit == 0.0 && !refinementFunc) PetscFunctionReturn(0);
1180d9deefdfSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1181d9deefdfSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1182c5929fdfSBarry Smith   ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr);
1183d9deefdfSMatthew G. Knepley   if (flg) name = genname;
1184d9deefdfSMatthew G. Knepley   if (name) {
1185d9deefdfSMatthew G. Knepley     ierr = PetscStrcmp(name, "triangle", &isTriangle);CHKERRQ(ierr);
1186d9deefdfSMatthew G. Knepley     ierr = PetscStrcmp(name, "tetgen",   &isTetgen);CHKERRQ(ierr);
1187d9deefdfSMatthew G. Knepley     ierr = PetscStrcmp(name, "ctetgen",  &isCTetgen);CHKERRQ(ierr);
1188d9deefdfSMatthew G. Knepley   }
1189d9deefdfSMatthew G. Knepley   switch (dim) {
1190d9deefdfSMatthew G. Knepley   case 2:
1191d9deefdfSMatthew G. Knepley     if (!name || isTriangle) {
1192d9deefdfSMatthew G. Knepley #if defined(PETSC_HAVE_TRIANGLE)
1193d9deefdfSMatthew G. Knepley       double  *maxVolumes;
1194d9deefdfSMatthew G. Knepley       PetscInt c;
1195d9deefdfSMatthew G. Knepley 
1196854ce69bSBarry Smith       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
1197713918a9SToby Isaac       if (adaptLabel) {
1198713918a9SToby Isaac         ierr = DMRefine_Plex_Label(dm,adaptLabel,cStart,cEnd,maxVolumes);CHKERRQ(ierr);
1199713918a9SToby Isaac       } else if (refinementFunc) {
1200b28003e6SMatthew G. Knepley         for (c = cStart; c < cEnd; ++c) {
1201b28003e6SMatthew G. Knepley           PetscReal vol, centroid[3];
1202e9ccc142SToby Isaac           PetscReal maxVol;
1203b28003e6SMatthew G. Knepley 
1204b28003e6SMatthew G. Knepley           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
1205e9ccc142SToby Isaac           ierr = (*refinementFunc)(centroid, &maxVol);CHKERRQ(ierr);
1206e9ccc142SToby Isaac           maxVolumes[c - cStart] = (double) maxVol;
1207b28003e6SMatthew G. Knepley         }
1208b28003e6SMatthew G. Knepley       } else {
1209d9deefdfSMatthew G. Knepley         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
1210b28003e6SMatthew G. Knepley       }
1211d9deefdfSMatthew G. Knepley       ierr = DMPlexRefine_Triangle(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
1212d9deefdfSMatthew G. Knepley       ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
1213d9deefdfSMatthew G. Knepley #else
1214d9deefdfSMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement needs external package support.\nPlease reconfigure with --download-triangle.");
1215d9deefdfSMatthew G. Knepley #endif
1216d9deefdfSMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 2D mesh generation package %s", name);
1217d9deefdfSMatthew G. Knepley     break;
1218d9deefdfSMatthew G. Knepley   case 3:
1219d9deefdfSMatthew G. Knepley     if (!name || isCTetgen) {
1220d9deefdfSMatthew G. Knepley #if defined(PETSC_HAVE_CTETGEN)
1221d9deefdfSMatthew G. Knepley       PetscReal *maxVolumes;
1222d9deefdfSMatthew G. Knepley       PetscInt   c;
1223d9deefdfSMatthew G. Knepley 
1224854ce69bSBarry Smith       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
1225713918a9SToby Isaac       if (adaptLabel) {
1226713918a9SToby Isaac         ierr = DMRefine_Plex_Label(dm,adaptLabel,cStart,cEnd,maxVolumes);CHKERRQ(ierr);
1227713918a9SToby Isaac       } else if (refinementFunc) {
1228b28003e6SMatthew G. Knepley         for (c = cStart; c < cEnd; ++c) {
1229b28003e6SMatthew G. Knepley           PetscReal vol, centroid[3];
1230b28003e6SMatthew G. Knepley 
1231b28003e6SMatthew G. Knepley           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
1232b28003e6SMatthew G. Knepley           ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr);
1233b28003e6SMatthew G. Knepley         }
1234b28003e6SMatthew G. Knepley       } else {
1235d9deefdfSMatthew G. Knepley         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
1236b28003e6SMatthew G. Knepley       }
1237d9deefdfSMatthew G. Knepley       ierr = DMPlexRefine_CTetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
1238d9deefdfSMatthew G. Knepley #else
1239d9deefdfSMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "CTetgen needs external package support.\nPlease reconfigure with --download-ctetgen.");
1240d9deefdfSMatthew G. Knepley #endif
1241d9deefdfSMatthew G. Knepley     } else if (isTetgen) {
1242d9deefdfSMatthew G. Knepley #if defined(PETSC_HAVE_TETGEN)
1243d9deefdfSMatthew G. Knepley       double  *maxVolumes;
1244d9deefdfSMatthew G. Knepley       PetscInt c;
1245d9deefdfSMatthew G. Knepley 
1246854ce69bSBarry Smith       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
1247713918a9SToby Isaac       if (adaptLabel) {
1248713918a9SToby Isaac         ierr = DMRefine_Plex_Label(dm,adaptLabel,cStart,cEnd,maxVolumes);CHKERRQ(ierr);
1249713918a9SToby Isaac       } else if (refinementFunc) {
1250b28003e6SMatthew G. Knepley         for (c = cStart; c < cEnd; ++c) {
1251b28003e6SMatthew G. Knepley           PetscReal vol, centroid[3];
1252b28003e6SMatthew G. Knepley 
1253b28003e6SMatthew G. Knepley           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
1254b28003e6SMatthew G. Knepley           ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr);
1255b28003e6SMatthew G. Knepley         }
1256b28003e6SMatthew G. Knepley       } else {
1257d9deefdfSMatthew G. Knepley         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
1258b28003e6SMatthew G. Knepley       }
1259d9deefdfSMatthew G. Knepley       ierr = DMPlexRefine_Tetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
1260d9deefdfSMatthew G. Knepley #else
1261d9deefdfSMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tetgen needs external package support.\nPlease reconfigure with --with-c-language=cxx --download-tetgen.");
1262d9deefdfSMatthew G. Knepley #endif
1263d9deefdfSMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 3D mesh generation package %s", name);
1264d9deefdfSMatthew G. Knepley     break;
1265d9deefdfSMatthew G. Knepley   default:
1266d9deefdfSMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %d is not supported.", dim);
1267d9deefdfSMatthew G. Knepley   }
1268a6ba4734SToby Isaac   ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr);
126936447a5eSToby Isaac   if (localized) {
127036447a5eSToby Isaac     ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);
127136447a5eSToby Isaac   }
1272d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
1273d9deefdfSMatthew G. Knepley }
1274d9deefdfSMatthew G. Knepley 
1275713918a9SToby Isaac PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined)
1276713918a9SToby Isaac {
1277713918a9SToby Isaac   PetscErrorCode ierr;
1278713918a9SToby Isaac 
1279713918a9SToby Isaac   PetscFunctionBegin;
1280713918a9SToby Isaac   ierr = DMRefine_Plex_Internal(dm,comm,NULL,dmRefined);CHKERRQ(ierr);
1281713918a9SToby Isaac   PetscFunctionReturn(0);
1282713918a9SToby Isaac }
1283713918a9SToby Isaac 
1284713918a9SToby Isaac PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted)
1285713918a9SToby Isaac {
1286713918a9SToby Isaac   PetscInt       defFlag, minFlag, maxFlag, numFlags, i;
1287713918a9SToby Isaac   const PetscInt *flags;
1288713918a9SToby Isaac   IS             flagIS;
1289713918a9SToby Isaac   PetscErrorCode ierr;
1290713918a9SToby Isaac 
1291713918a9SToby Isaac   PetscFunctionBegin;
1292713918a9SToby Isaac   ierr = DMLabelGetDefaultValue(adaptLabel,&defFlag);CHKERRQ(ierr);
1293713918a9SToby Isaac   minFlag = defFlag;
1294713918a9SToby Isaac   maxFlag = defFlag;
1295713918a9SToby Isaac   ierr = DMLabelGetValueIS(adaptLabel,&flagIS);CHKERRQ(ierr);
1296713918a9SToby Isaac   ierr = ISGetLocalSize(flagIS,&numFlags);CHKERRQ(ierr);
1297713918a9SToby Isaac   ierr = ISGetIndices(flagIS,&flags);CHKERRQ(ierr);
1298713918a9SToby Isaac   for (i = 0; i < numFlags; i++) {
1299713918a9SToby Isaac     PetscInt flag = flags[i];
1300713918a9SToby Isaac 
1301713918a9SToby Isaac     minFlag = PetscMin(minFlag,flag);
1302713918a9SToby Isaac     maxFlag = PetscMax(maxFlag,flag);
1303713918a9SToby Isaac   }
1304713918a9SToby Isaac   ierr = ISRestoreIndices(flagIS,&flags);CHKERRQ(ierr);
1305713918a9SToby Isaac   ierr = ISDestroy(&flagIS);CHKERRQ(ierr);
1306713918a9SToby Isaac   {
13070e419529SToby Isaac     PetscInt minMaxFlag[2], minMaxFlagGlobal[2];
1308713918a9SToby Isaac 
13090e419529SToby Isaac     minMaxFlag[0] = minFlag;
13100e419529SToby Isaac     minMaxFlag[1] = -maxFlag;
1311713918a9SToby Isaac     ierr = MPI_Allreduce(minMaxFlag,minMaxFlagGlobal,2,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
1312713918a9SToby Isaac     minFlag = minMaxFlagGlobal[0];
1313713918a9SToby Isaac     maxFlag = -minMaxFlagGlobal[1];
1314713918a9SToby Isaac   }
1315713918a9SToby Isaac   if (minFlag == maxFlag) {
1316713918a9SToby Isaac     switch (minFlag) {
1317cd3c525cSToby Isaac     case DM_ADAPT_DETERMINE:
1318713918a9SToby Isaac       *dmAdapted = NULL;
1319713918a9SToby Isaac       break;
1320713918a9SToby Isaac     case DM_ADAPT_REFINE:
1321713918a9SToby Isaac       ierr = DMPlexSetRefinementUniform(dm,PETSC_TRUE);CHKERRQ(ierr);
1322713918a9SToby Isaac       ierr = DMRefine(dm,MPI_COMM_NULL,dmAdapted);CHKERRQ(ierr);
1323713918a9SToby Isaac       break;
1324713918a9SToby Isaac     case DM_ADAPT_COARSEN:
1325713918a9SToby Isaac       ierr = DMCoarsen(dm,MPI_COMM_NULL,dmAdapted);CHKERRQ(ierr);
1326713918a9SToby Isaac       break;
1327713918a9SToby Isaac     default:
1328713918a9SToby Isaac       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n",minFlag);
1329713918a9SToby Isaac       break;
1330713918a9SToby Isaac     }
1331713918a9SToby Isaac   } else {
1332713918a9SToby Isaac     ierr = DMPlexSetRefinementUniform(dm,PETSC_FALSE);CHKERRQ(ierr);
1333713918a9SToby Isaac     ierr = DMRefine_Plex_Internal(dm,MPI_COMM_NULL,adaptLabel,dmAdapted);CHKERRQ(ierr);
1334713918a9SToby Isaac   }
1335713918a9SToby Isaac   PetscFunctionReturn(0);
1336713918a9SToby Isaac }
1337713918a9SToby Isaac 
1338d9deefdfSMatthew G. Knepley PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[])
1339d9deefdfSMatthew G. Knepley {
1340d9deefdfSMatthew G. Knepley   DM             cdm = dm;
1341d9deefdfSMatthew G. Knepley   PetscInt       r;
134236447a5eSToby Isaac   PetscBool      isUniform, localized;
1343d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
1344d9deefdfSMatthew G. Knepley 
1345d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
1346d9deefdfSMatthew G. Knepley   ierr = DMPlexGetRefinementUniform(dm, &isUniform);CHKERRQ(ierr);
134736447a5eSToby Isaac   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1348d9deefdfSMatthew G. Knepley   if (!isUniform) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Non-uniform refinement is incompatible with the hierarchy");
1349d9deefdfSMatthew G. Knepley   for (r = 0; r < nlevels; ++r) {
1350d9deefdfSMatthew G. Knepley     CellRefiner cellRefiner;
1351d9deefdfSMatthew G. Knepley 
1352d9deefdfSMatthew G. Knepley     ierr = DMPlexGetCellRefiner_Internal(cdm, &cellRefiner);CHKERRQ(ierr);
1353d9deefdfSMatthew G. Knepley     ierr = DMPlexRefineUniform_Internal(cdm, cellRefiner, &dmRefined[r]);CHKERRQ(ierr);
1354a6ba4734SToby Isaac     ierr = DMCopyBoundary(cdm, dmRefined[r]);CHKERRQ(ierr);
135536447a5eSToby Isaac     if (localized) {
135636447a5eSToby Isaac       ierr = DMLocalizeCoordinates(dmRefined[r]);CHKERRQ(ierr);
135736447a5eSToby Isaac     }
1358a8fb8f29SToby Isaac     ierr = DMSetCoarseDM(dmRefined[r], cdm);CHKERRQ(ierr);
13590aef6b92SMatthew G. Knepley     ierr = DMPlexSetRegularRefinement(dmRefined[r], PETSC_TRUE);CHKERRQ(ierr);
1360d9deefdfSMatthew G. Knepley     cdm  = dmRefined[r];
1361d9deefdfSMatthew G. Knepley   }
1362d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
1363d9deefdfSMatthew G. Knepley }
1364