xref: /petsc/src/dm/impls/plex/tests/ex1.c (revision efbe7e8a80d07327753dbe0b33efee01e046af3f)
1 static char help[] = "Tests various DMPlex routines to construct, refine and distribute a mesh.\n\n";
2 
3 #include <petscdmplex.h>
4 #include <petscsf.h>
5 
6 typedef enum {BOX, CYLINDER, SPHERE, BALL} DomainShape;
7 enum {STAGE_LOAD, STAGE_DISTRIBUTE, STAGE_REFINE, STAGE_OVERLAP};
8 
9 typedef struct {
10   DM            dm;                /* REQUIRED in order to use SNES evaluation functions */
11   PetscInt      debug;             /* The debugging level */
12   PetscLogEvent createMeshEvent;
13   PetscLogStage stages[4];
14   /* Domain and mesh definition */
15   PetscInt      dim;                             /* The topological mesh dimension */
16   PetscBool     interpolate;                     /* Generate intermediate mesh elements */
17   PetscReal     refinementLimit;                 /* The largest allowable cell volume */
18   PetscBool     cellSimplex;                     /* Use simplices or hexes */
19   PetscBool     cellWedge;                       /* Use wedges */
20   DomainShape   domainShape;                     /* Shape of the region to be meshed */
21   PetscInt      *domainBoxSizes;                 /* Sizes of the box mesh */
22   PetscReal     *domainBoxL,*domainBoxU;         /* Lower left, upper right corner of the box mesh */
23   DMBoundaryType periodicity[3];                 /* The domain periodicity */
24   char          filename[PETSC_MAX_PATH_LEN];    /* Import mesh from file */
25   char          bdfilename[PETSC_MAX_PATH_LEN];  /* Import mesh boundary from file */
26   char          extfilename[PETSC_MAX_PATH_LEN]; /* Import 2D mesh to be extruded from file */
27   PetscBool     testPartition;                   /* Use a fixed partitioning for testing */
28   PetscInt      overlap;                         /* The cell overlap to use during partitioning */
29   PetscReal     extrude_thickness;               /* Thickness of extrusion */
30   PetscInt      extrude_layers;                  /* Layers to be extruded */
31   PetscBool     extrude_hfirst;                  /* New numbering: height first? */
32   PetscBool     testp4est[2];
33   PetscBool     redistribute;
34   PetscBool     final_ref;                       /* Run refinement at the end */
35   PetscBool     final_diagnostics;               /* Run diagnostics on the final mesh */
36 } AppCtx;
37 
38 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
39 {
40   const char       *dShapes[4] = {"box", "cylinder", "sphere", "ball"};
41   PetscInt         shape, bd, n;
42   static PetscInt  domainBoxSizes[3] = {1,1,1};
43   static PetscReal domainBoxL[3] = {0.,0.,0.};
44   static PetscReal domainBoxU[3] = {1.,1.,1.};
45   PetscBool        flg;
46   PetscErrorCode   ierr;
47 
48   PetscFunctionBegin;
49   options->debug             = 0;
50   options->dim               = 2;
51   options->interpolate       = PETSC_FALSE;
52   options->refinementLimit   = 0.0;
53   options->cellSimplex       = PETSC_TRUE;
54   options->cellWedge         = PETSC_FALSE;
55   options->domainShape       = BOX;
56   options->domainBoxSizes    = NULL;
57   options->domainBoxL        = NULL;
58   options->domainBoxU        = NULL;
59   options->periodicity[0]    = DM_BOUNDARY_NONE;
60   options->periodicity[1]    = DM_BOUNDARY_NONE;
61   options->periodicity[2]    = DM_BOUNDARY_NONE;
62   options->filename[0]       = '\0';
63   options->bdfilename[0]     = '\0';
64   options->extfilename[0]    = '\0';
65   options->testPartition     = PETSC_FALSE;
66   options->overlap           = 0;
67   options->extrude_layers    = 0;
68   options->extrude_thickness = 0.1;
69   options->extrude_hfirst    = PETSC_TRUE;
70   options->testp4est[0]      = PETSC_FALSE;
71   options->testp4est[1]      = PETSC_FALSE;
72   options->redistribute      = PETSC_FALSE;
73   options->final_ref         = PETSC_FALSE;
74   options->final_diagnostics = PETSC_TRUE;
75 
76   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
77   ierr = PetscOptionsBoundedInt("-debug", "The debugging level", "ex1.c", options->debug, &options->debug, NULL,0);CHKERRQ(ierr);
78   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex1.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
79   ierr = PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex1.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr);
80   ierr = PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex1.c", options->refinementLimit, &options->refinementLimit, NULL);CHKERRQ(ierr);
81   ierr = PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex1.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr);
82   ierr = PetscOptionsBool("-cell_wedge", "Use wedges if true", "ex1.c", options->cellWedge, &options->cellWedge, NULL);CHKERRQ(ierr);
83   shape = options->domainShape;
84   ierr = PetscOptionsEList("-domain_shape","The shape of the domain","ex1.c", dShapes, 4, dShapes[options->domainShape], &shape, NULL);CHKERRQ(ierr);
85   options->domainShape = (DomainShape) shape;
86   ierr = PetscOptionsIntArray("-domain_box_sizes","The sizes of the box domain","ex1.c", domainBoxSizes, (n=3,&n), &flg);CHKERRQ(ierr);
87   if (flg) { options->domainShape = BOX; options->domainBoxSizes = domainBoxSizes;}
88   ierr = PetscOptionsRealArray("-domain_box_ll","Coordinates of the lower left corner of the box domain","ex1.c", domainBoxL, (n=3,&n), &flg);CHKERRQ(ierr);
89   if (flg) { options->domainBoxL = domainBoxL;}
90   ierr = PetscOptionsRealArray("-domain_box_ur","Coordinates of the upper right corner of the box domain","ex1.c", domainBoxU, (n=3,&n), &flg);CHKERRQ(ierr);
91   if (flg) { options->domainBoxU = domainBoxU;}
92   bd = options->periodicity[0];
93   ierr = PetscOptionsEList("-x_periodicity", "The x-boundary periodicity", "ex1.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[0]], &bd, NULL);CHKERRQ(ierr);
94   options->periodicity[0] = (DMBoundaryType) bd;
95   bd = options->periodicity[1];
96   ierr = PetscOptionsEList("-y_periodicity", "The y-boundary periodicity", "ex1.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[1]], &bd, NULL);CHKERRQ(ierr);
97   options->periodicity[1] = (DMBoundaryType) bd;
98   bd = options->periodicity[2];
99   ierr = PetscOptionsEList("-z_periodicity", "The z-boundary periodicity", "ex1.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[2]], &bd, NULL);CHKERRQ(ierr);
100   options->periodicity[2] = (DMBoundaryType) bd;
101   ierr = PetscOptionsString("-filename", "The mesh file", "ex1.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
102   ierr = PetscOptionsString("-bd_filename", "The mesh boundary file", "ex1.c", options->bdfilename, options->bdfilename, sizeof(options->bdfilename), NULL);CHKERRQ(ierr);
103   ierr = PetscOptionsString("-ext_filename", "The 2D mesh file to be extruded", "ex1.c", options->extfilename, options->extfilename, sizeof(options->extfilename), &flg);CHKERRQ(ierr);
104   if (flg) options->extrude_layers = 2;
105   ierr = PetscOptionsBoundedInt("-ext_layers", "The number of layers to extrude", "ex1.c", options->extrude_layers, &options->extrude_layers, NULL,0);CHKERRQ(ierr);
106   ierr = PetscOptionsReal("-ext_thickness", "The thickness of the layer to be extruded", "ex1.c", options->extrude_thickness, &options->extrude_thickness, NULL);CHKERRQ(ierr);
107   ierr = PetscOptionsBool("-ext_hfirst", "Order the cells in the height first", "ex1.c", options->extrude_hfirst, &options->extrude_hfirst, NULL);CHKERRQ(ierr);
108   ierr = PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex1.c", options->testPartition, &options->testPartition, NULL);CHKERRQ(ierr);
109   ierr = PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex1.c", options->overlap, &options->overlap, NULL,0);CHKERRQ(ierr);
110   ierr = PetscOptionsBool("-test_p4est_seq", "Test p4est with sequential base DM", "ex1.c", options->testp4est[0], &options->testp4est[0], NULL);CHKERRQ(ierr);
111   ierr = PetscOptionsBool("-test_p4est_par", "Test p4est with parallel base DM", "ex1.c", options->testp4est[1], &options->testp4est[1], NULL);CHKERRQ(ierr);
112   ierr = PetscOptionsBool("-test_redistribute", "Test redistribution", "ex1.c", options->redistribute, &options->redistribute, NULL);CHKERRQ(ierr);
113   ierr = PetscOptionsBool("-final_ref", "Run uniform refinement on the final mesh", "ex1.c", options->final_ref, &options->final_ref, NULL);CHKERRQ(ierr);
114   ierr = PetscOptionsBool("-final_diagnostics", "Run diagnostics on the final mesh", "ex1.c", options->final_diagnostics, &options->final_diagnostics, NULL);CHKERRQ(ierr);
115   ierr = PetscOptionsEnd();
116 
117   ierr = PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);CHKERRQ(ierr);
118   ierr = PetscLogStageRegister("MeshLoad",       &options->stages[STAGE_LOAD]);CHKERRQ(ierr);
119   ierr = PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr);
120   ierr = PetscLogStageRegister("MeshRefine",     &options->stages[STAGE_REFINE]);CHKERRQ(ierr);
121   ierr = PetscLogStageRegister("MeshOverlap",    &options->stages[STAGE_OVERLAP]);CHKERRQ(ierr);
122   PetscFunctionReturn(0);
123 }
124 
125 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
126 {
127   PetscInt       dim                  = user->dim;
128   PetscBool      interpolate          = user->interpolate;
129   PetscReal      refinementLimit      = user->refinementLimit;
130   PetscBool      cellSimplex          = user->cellSimplex;
131   PetscBool      cellWedge            = user->cellWedge;
132   const char    *filename             = user->filename;
133   const char    *bdfilename           = user->bdfilename;
134   const char    *extfilename          = user->extfilename;
135   PetscBool      testp4est_seq        = user->testp4est[0];
136   PetscBool      testp4est_par        = user->testp4est[1];
137   PetscInt       triSizes_n2[2]       = {4, 4};
138   PetscInt       triPoints_n2[8]      = {3, 5, 6, 7, 0, 1, 2, 4};
139   PetscInt       triSizes_n8[8]       = {1, 1, 1, 1, 1, 1, 1, 1};
140   PetscInt       triPoints_n8[8]      = {0, 1, 2, 3, 4, 5, 6, 7};
141   PetscInt       quadSizes[2]         = {2, 2};
142   PetscInt       quadPoints[4]        = {2, 3, 0, 1};
143   PetscInt       gmshSizes_n3[3]      = {14, 14, 14};
144   PetscInt       gmshPoints_n3[42]    = {1, 2,  4,  5,  9, 10, 11, 15, 16, 20, 21, 27, 28, 29,
145                                          3, 8, 12, 30, 31, 32, 33, 34, 35, 36, 37, 38, 39, 40,
146                                          0, 6,  7, 13, 14, 17, 18, 19, 22, 23, 24, 25, 26, 41};
147   PetscInt       fluentSizes_n3[3]    = {50, 50, 50};
148   PetscInt       fluentPoints_n3[150] = { 5,  6,  7,  8, 12, 14, 16,  34,  36,  37,  38,  39,  40,  41,  42,  43,  44,  45,  46,  48,  50,  51,  80,  81,  89,
149                                          91, 93, 94, 95, 96, 97, 98,  99, 100, 101, 104, 121, 122, 124, 125, 126, 127, 128, 129, 131, 133, 143, 144, 145, 147,
150                                           1,  3,  4,  9, 10, 17, 18,  19,  24,  25,  26,  27,  28,  29,  30,  31,  32,  33,  35,  47,  61,  71,  72,  73,  74,
151                                          75, 76, 77, 78, 79, 86, 87,  88,  90,  92, 113, 115, 116, 117, 118, 119, 120, 123, 138, 140, 141, 142, 146, 148, 149,
152                                           0,  2, 11, 13, 15, 20, 21,  22,  23,  49,  52,  53,  54,  55,  56,  57,  58,  59,  60,  62,  63,  64,  65,  66,  67,
153                                          68, 69, 70, 82, 83, 84, 85, 102, 103, 105, 106, 107, 108, 109, 110, 111, 112, 114, 130, 132, 134, 135, 136, 137, 139};
154   size_t         len, bdlen, extlen;
155   PetscMPIInt    rank, size;
156   PetscBool      periodic;
157   PetscErrorCode ierr;
158 
159   PetscFunctionBegin;
160   ierr = PetscLogEventBegin(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr);
161   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
162   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
163   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
164   ierr = PetscStrlen(bdfilename, &bdlen);CHKERRQ(ierr);
165   ierr = PetscStrlen(extfilename, &extlen);CHKERRQ(ierr);
166   ierr = PetscLogStagePush(user->stages[STAGE_LOAD]);CHKERRQ(ierr);
167   if (len) {
168     ierr = DMPlexCreateFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
169   } else if (bdlen) {
170     DM boundary;
171 
172     ierr = DMPlexCreateFromFile(comm, bdfilename, interpolate, &boundary);CHKERRQ(ierr);
173     ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr);
174     ierr = DMDestroy(&boundary);CHKERRQ(ierr);
175   } else if (extlen) {
176     DM edm;
177 
178     ierr = DMPlexCreateFromFile(comm, extfilename, interpolate, &edm);CHKERRQ(ierr);
179     ierr = DMPlexExtrude(edm, user->extrude_layers, user->extrude_thickness, user->extrude_hfirst, NULL, interpolate, dm);CHKERRQ(ierr);
180     ierr = DMDestroy(&edm);CHKERRQ(ierr);
181   } else {
182     switch (user->domainShape) {
183     case BOX:
184       if (cellWedge) {
185         if (dim != 3) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Dimension must be 3 for a wedge mesh, not %D", dim);
186         ierr = DMPlexCreateWedgeBoxMesh(comm, user->domainBoxSizes, user->domainBoxL, user->domainBoxU, user->periodicity, PETSC_FALSE, interpolate, dm);CHKERRQ(ierr);
187       } else {
188         ierr = DMPlexCreateBoxMesh(comm, dim, cellSimplex, user->domainBoxSizes, user->domainBoxL, user->domainBoxU, user->periodicity, interpolate, dm);CHKERRQ(ierr);
189       }
190       break;
191     case CYLINDER:
192       if (cellSimplex) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot mesh a cylinder with simplices");
193       if (dim != 3)    SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Dimension must be 3 for a cylinder mesh, not %D", dim);
194       if (cellWedge) {
195         ierr = DMPlexCreateWedgeCylinderMesh(comm, 6, interpolate, dm);CHKERRQ(ierr);
196       } else {
197         ierr = DMPlexCreateHexCylinderMesh(comm, 1, user->periodicity[2], dm);CHKERRQ(ierr);
198       }
199       break;
200     case SPHERE:
201       ierr = DMPlexCreateSphereMesh(comm, dim, cellSimplex, 1.0, dm);CHKERRQ(ierr);
202       break;
203     case BALL:
204       ierr = DMPlexCreateBallMesh(comm, dim, 1.0, dm);CHKERRQ(ierr);
205       break;
206     default: SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Unknown domain shape %D", user->domainShape);
207     }
208   }
209   if (!extlen && user->extrude_layers > 0) {
210     DM edm;
211 
212     ierr = DMPlexExtrude(*dm, user->extrude_layers, user->extrude_thickness, user->extrude_hfirst, NULL, interpolate, &edm);CHKERRQ(ierr);
213     ierr = DMDestroy(dm);CHKERRQ(ierr);
214     *dm  = edm;
215   }
216 
217   /* For topologically periodic meshes, we first localize coordinates,
218      and then remove any information related with the
219      automatic computation of localized vertices.
220      This way, refinement operations and conversions to p4est
221      will preserve the shape of the domain in physical space */
222   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
223   ierr = DMGetPeriodicity(*dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr);
224   if (periodic) {
225     ierr = DMSetPeriodicity(*dm,PETSC_TRUE,NULL,NULL,NULL);CHKERRQ(ierr);
226   }
227 
228   ierr = DMViewFromOptions(*dm,NULL,"-init_dm_view");CHKERRQ(ierr);
229   ierr = DMGetDimension(*dm,&dim);CHKERRQ(ierr);
230 
231   if (testp4est_seq) {
232 #if defined(PETSC_HAVE_P4EST)
233     DM dmConv = NULL;
234 
235     ierr = DMPlexCheckSymmetry(*dm);CHKERRQ(ierr);
236     ierr = DMPlexCheckSkeleton(*dm, 0);CHKERRQ(ierr);
237     ierr = DMPlexCheckFaces(*dm, 0);CHKERRQ(ierr);
238     ierr = DMPlexCheckGeometry(*dm);CHKERRQ(ierr);
239     ierr = DMPlexCheckPointSF(*dm);CHKERRQ(ierr);
240     ierr = DMPlexCheckInterfaceCones(*dm);CHKERRQ(ierr);
241     ierr = DMPlexSetRefinementUniform(*dm, PETSC_TRUE);CHKERRQ(ierr);
242     ierr = DMPlexSetCellRefinerType(*dm, DM_REFINER_TO_BOX);CHKERRQ(ierr);
243     ierr = DMRefine(*dm, PETSC_COMM_WORLD, &dmConv);CHKERRQ(ierr);
244     if (dmConv) {
245       ierr = DMDestroy(dm);CHKERRQ(ierr);
246       *dm  = dmConv;
247     }
248     ierr = DMViewFromOptions(*dm,NULL,"-initref_dm_view");CHKERRQ(ierr);
249     ierr = DMPlexCheckSymmetry(*dm);CHKERRQ(ierr);
250     ierr = DMPlexCheckSkeleton(*dm, 0);CHKERRQ(ierr);
251     ierr = DMPlexCheckFaces(*dm, 0);CHKERRQ(ierr);
252     ierr = DMPlexCheckGeometry(*dm);CHKERRQ(ierr);
253     ierr = DMPlexCheckPointSF(*dm);CHKERRQ(ierr);
254     ierr = DMPlexCheckInterfaceCones(*dm);CHKERRQ(ierr);
255     user->cellSimplex = PETSC_FALSE;
256 
257     ierr = DMConvert(*dm,dim == 2 ? DMP4EST : DMP8EST,&dmConv);CHKERRQ(ierr);
258     if (dmConv) {
259       ierr = PetscObjectSetOptionsPrefix((PetscObject) dmConv, "conv_seq_1_");CHKERRQ(ierr);
260       ierr = DMSetFromOptions(dmConv);CHKERRQ(ierr);
261       ierr = DMDestroy(dm);CHKERRQ(ierr);
262       *dm  = dmConv;
263     }
264     ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "conv_seq_1_");CHKERRQ(ierr);
265     ierr = DMSetUp(*dm);CHKERRQ(ierr);
266     ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
267     ierr = DMConvert(*dm,DMPLEX,&dmConv);CHKERRQ(ierr);
268     if (dmConv) {
269       ierr = PetscObjectSetOptionsPrefix((PetscObject) dmConv, "conv_seq_2_");CHKERRQ(ierr);
270       ierr = DMSetFromOptions(dmConv);CHKERRQ(ierr);
271       ierr = DMDestroy(dm);CHKERRQ(ierr);
272       *dm  = dmConv;
273     }
274     ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "conv_seq_2_");CHKERRQ(ierr);
275     ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
276     ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL);CHKERRQ(ierr);
277 #else
278     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with --download-p4est");
279 #endif
280   }
281 
282   ierr = PetscLogStagePop();CHKERRQ(ierr);
283   if (!testp4est_seq) {
284     DM refinedMesh     = NULL;
285     DM distributedMesh = NULL;
286 
287     if (user->testPartition) {
288       const PetscInt  *sizes = NULL;
289       const PetscInt  *points = NULL;
290       PetscPartitioner part;
291 
292       if (!rank) {
293         if (dim == 2 && cellSimplex && size == 2) {
294            sizes = triSizes_n2; points = triPoints_n2;
295         } else if (dim == 2 && cellSimplex && size == 8) {
296           sizes = triSizes_n8; points = triPoints_n8;
297         } else if (dim == 2 && !cellSimplex && size == 2) {
298           sizes = quadSizes; points = quadPoints;
299         } else if (dim == 2 && size == 3) {
300           PetscInt Nc;
301 
302           ierr = DMPlexGetHeightStratum(*dm, 0, NULL, &Nc);CHKERRQ(ierr);
303           if (Nc == 42) { /* Gmsh 3 & 4 */
304             sizes = gmshSizes_n3; points = gmshPoints_n3;
305           } else if (Nc == 150) { /* Fluent 1 */
306             sizes = fluentSizes_n3; points = fluentPoints_n3;
307           } else if (Nc == 42) { /* Med 1 */
308           } else if (Nc == 161) { /* Med 3 */
309           }
310         }
311       }
312       ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
313       ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr);
314       ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr);
315     } else {
316       PetscPartitioner part;
317 
318       ierr = DMPlexGetPartitioner(*dm,&part);CHKERRQ(ierr);
319       ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
320     }
321     /* Distribute mesh over processes */
322     ierr = PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr);
323     ierr = DMViewFromOptions(*dm, NULL, "-dm_pre_dist_view");CHKERRQ(ierr);
324     ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
325     if (distributedMesh) {
326       ierr = DMDestroy(dm);CHKERRQ(ierr);
327       *dm  = distributedMesh;
328     }
329     ierr = PetscLogStagePop();CHKERRQ(ierr);
330     ierr = DMViewFromOptions(*dm, NULL, "-distributed_dm_view");CHKERRQ(ierr);
331     /* Refine mesh using a volume constraint */
332     ierr = PetscLogStagePush(user->stages[STAGE_REFINE]);CHKERRQ(ierr);
333     ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr);
334     ierr = DMPlexSetRefinementLimit(*dm, refinementLimit);CHKERRQ(ierr);
335     ierr = DMRefine(*dm, comm, &refinedMesh);CHKERRQ(ierr);
336     if (refinedMesh) {
337       ierr = DMDestroy(dm);CHKERRQ(ierr);
338       *dm  = refinedMesh;
339     }
340     ierr = PetscLogStagePop();CHKERRQ(ierr);
341   }
342   ierr = PetscLogStagePush(user->stages[STAGE_REFINE]);CHKERRQ(ierr);
343   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
344   ierr = PetscLogStagePop();CHKERRQ(ierr);
345 
346   if (testp4est_par) {
347 #if defined(PETSC_HAVE_P4EST)
348     DM dmConv = NULL;
349 
350     ierr = DMViewFromOptions(*dm, NULL, "-dm_tobox_view");CHKERRQ(ierr);
351     ierr = DMPlexSetRefinementUniform(*dm, PETSC_TRUE);CHKERRQ(ierr);
352     ierr = DMPlexSetCellRefinerType(*dm, DM_REFINER_TO_BOX);CHKERRQ(ierr);
353     ierr = DMRefine(*dm, PETSC_COMM_WORLD, &dmConv);CHKERRQ(ierr);
354     if (dmConv) {
355       ierr = DMDestroy(dm);CHKERRQ(ierr);
356       *dm  = dmConv;
357     }
358     user->cellSimplex = PETSC_FALSE;
359     ierr = DMViewFromOptions(*dm, NULL, "-dm_tobox_view");CHKERRQ(ierr);
360     ierr = DMPlexCheckSymmetry(*dm);CHKERRQ(ierr);
361     ierr = DMPlexCheckSkeleton(*dm, 0);CHKERRQ(ierr);
362     ierr = DMPlexCheckFaces(*dm, 0);CHKERRQ(ierr);
363     ierr = DMPlexCheckGeometry(*dm);CHKERRQ(ierr);
364     ierr = DMPlexCheckPointSF(*dm);CHKERRQ(ierr);
365     ierr = DMPlexCheckInterfaceCones(*dm);CHKERRQ(ierr);
366 
367     ierr = DMConvert(*dm,dim == 2 ? DMP4EST : DMP8EST,&dmConv);CHKERRQ(ierr);
368     if (dmConv) {
369       ierr = PetscObjectSetOptionsPrefix((PetscObject) dmConv, "conv_par_1_");CHKERRQ(ierr);
370       ierr = DMSetFromOptions(dmConv);CHKERRQ(ierr);
371       ierr = DMDestroy(dm);CHKERRQ(ierr);
372       *dm  = dmConv;
373     }
374     ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "conv_par_1_");CHKERRQ(ierr);
375     ierr = DMSetUp(*dm);CHKERRQ(ierr);
376     ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
377     ierr = DMConvert(*dm, DMPLEX, &dmConv);CHKERRQ(ierr);
378     if (dmConv) {
379       ierr = PetscObjectSetOptionsPrefix((PetscObject) dmConv, "conv_par_2_");CHKERRQ(ierr);
380       ierr = DMSetFromOptions(dmConv);CHKERRQ(ierr);
381       ierr = DMDestroy(dm);CHKERRQ(ierr);
382       *dm  = dmConv;
383     }
384     ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "conv_par_2_");CHKERRQ(ierr);
385     ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
386     ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL);CHKERRQ(ierr);
387 #else
388     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with --download-p4est");
389 #endif
390   }
391 
392   /* test redistribution of an already distributed mesh */
393   if (user->redistribute) {
394     DM       distributedMesh;
395     PetscSF  sf;
396     PetscInt nranks;
397 
398     ierr = DMViewFromOptions(*dm, NULL, "-dm_pre_redist_view");CHKERRQ(ierr);
399     ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
400     if (distributedMesh) {
401       ierr = DMGetPointSF(distributedMesh, &sf);CHKERRQ(ierr);
402       ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
403       ierr = DMGetNeighbors(distributedMesh, &nranks, NULL);CHKERRQ(ierr);
404       ierr = MPI_Allreduce(MPI_IN_PLACE, &nranks, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)*dm));CHKERRQ(ierr);
405       ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)*dm)), "Minimum number of neighbors: %D\n", nranks);CHKERRQ(ierr);
406       ierr = DMDestroy(dm);CHKERRQ(ierr);
407       *dm  = distributedMesh;
408     }
409     ierr = DMViewFromOptions(*dm, NULL, "-dm_post_redist_view");CHKERRQ(ierr);
410   }
411 
412   if (user->overlap) {
413     DM overlapMesh = NULL;
414 
415     /* Add the overlap to refined mesh */
416     ierr = PetscLogStagePush(user->stages[STAGE_OVERLAP]);CHKERRQ(ierr);
417     ierr = DMViewFromOptions(*dm, NULL, "-dm_pre_overlap_view");CHKERRQ(ierr);
418     ierr = DMPlexDistributeOverlap(*dm, user->overlap, NULL, &overlapMesh);CHKERRQ(ierr);
419     if (overlapMesh) {
420       PetscInt overlap;
421       ierr = DMPlexGetOverlap(overlapMesh, &overlap);CHKERRQ(ierr);
422       ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "Overlap: %D\n", overlap);CHKERRQ(ierr);
423       ierr = DMDestroy(dm);CHKERRQ(ierr);
424       *dm = overlapMesh;
425     }
426     ierr = DMViewFromOptions(*dm, NULL, "-dm_post_overlap_view");CHKERRQ(ierr);
427     ierr = PetscLogStagePop();CHKERRQ(ierr);
428   }
429   if (user->final_ref) {
430     DM refinedMesh = NULL;
431 
432     ierr = DMPlexSetRefinementUniform(*dm, PETSC_TRUE);CHKERRQ(ierr);
433     ierr = DMRefine(*dm, comm, &refinedMesh);CHKERRQ(ierr);
434     if (refinedMesh) {
435       ierr = DMDestroy(dm);CHKERRQ(ierr);
436       *dm  = refinedMesh;
437     }
438   }
439 
440   ierr = PetscObjectSetName((PetscObject) *dm, "Simplicial Mesh");CHKERRQ(ierr);
441   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
442   if (user->final_diagnostics) {
443     DMPlexInterpolatedFlag interpolated;
444     PetscInt  dim, depth;
445 
446     ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
447     ierr = DMPlexGetDepth(*dm, &depth);CHKERRQ(ierr);
448     ierr = DMPlexIsInterpolatedCollective(*dm, &interpolated);CHKERRQ(ierr);
449 
450     ierr = DMPlexCheckSymmetry(*dm);CHKERRQ(ierr);
451     if (interpolated == DMPLEX_INTERPOLATED_FULL) {
452       ierr = DMPlexCheckFaces(*dm, 0);CHKERRQ(ierr);
453     }
454     ierr = DMPlexCheckSkeleton(*dm, 0);CHKERRQ(ierr);
455     ierr = DMPlexCheckGeometry(*dm);CHKERRQ(ierr);
456   }
457   ierr = PetscLogEventEnd(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr);
458   user->dm = *dm;
459   PetscFunctionReturn(0);
460 }
461 
462 int main(int argc, char **argv)
463 {
464   AppCtx         user;                 /* user-defined work context */
465   PetscErrorCode ierr;
466 
467   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
468   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
469   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &user.dm);CHKERRQ(ierr);
470   ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
471   ierr = PetscFinalize();
472   return ierr;
473 }
474 
475 /*TEST
476 
477   # CTetGen 0-1
478   test:
479     suffix: 0
480     requires: ctetgen
481     args: -dim 3 -ctetgen_verbose 4 -dm_view ascii::ascii_info_detail -info :~sys
482   test:
483     suffix: 1
484     requires: ctetgen
485     args: -dim 3 -ctetgen_verbose 4 -refinement_limit 0.0625 -dm_view ascii::ascii_info_detail -info :~sys
486 
487 
488   # 2D LaTex and ASCII output 2-9
489   test:
490     suffix: 2
491     requires: triangle
492     args: -dim 2 -dm_view ascii::ascii_latex
493   test:
494     suffix: 3
495     requires: triangle
496     args: -dim 2 -dm_refine 1 -interpolate 1 -dm_view ascii::ascii_info_detail
497   test:
498     suffix: 4
499     requires: triangle
500     nsize: 2
501     args: -dim 2 -dm_refine 1 -interpolate 1 -test_partition -dm_view ascii::ascii_info_detail
502   test:
503     suffix: 5
504     requires: triangle
505     nsize: 2
506     args: -dim 2 -dm_refine 1 -interpolate 1 -test_partition -dm_view ascii::ascii_latex
507   test:
508     suffix: 6
509     args: -dim 2 -cell_simplex 0 -interpolate -dm_view ascii::ascii_info_detail
510   test:
511     suffix: 7
512     args: -dim 2 -cell_simplex 0 -interpolate -dm_refine 1 -dm_view ascii::ascii_info_detail
513   test:
514     suffix: 8
515     nsize: 2
516     args: -dim 2 -cell_simplex 0 -interpolate -dm_refine 1 -interpolate 1 -test_partition -dm_view ascii::ascii_latex
517 
518   # 1D ASCII output
519   test:
520     suffix: 1d_0
521     args: -dim 1 -domain_shape box -dm_view ascii::ascii_info_detail
522   test:
523     suffix: 1d_1
524     args: -dim 1 -domain_shape box -dm_refine 2 -dm_view ascii::ascii_info_detail
525   test:
526     suffix: 1d_2
527     args: -dim 1 -domain_box_sizes 5 -x_periodicity periodic -dm_view ascii::ascii_info_detail -dm_plex_check_all
528 
529   # Parallel refinement tests with overlap
530   test:
531     suffix: refine_overlap_1d
532     nsize: 2
533     args: -dim 1 -domain_box_sizes 4 -dm_refine 1 -overlap {{0 1 2}separate output} -petscpartitioner_type simple -dm_view ascii::ascii_info
534   test:
535     suffix: refine_overlap_2d
536     requires: triangle
537     nsize: {{2 8}separate output}
538     args: -dim 2 -cell_simplex 1 -dm_refine 1 -interpolate 1 -test_partition -overlap {{0 1 2}separate output} -dm_view ascii::ascii_info
539 
540   # Parallel simple partitioner tests
541   test:
542     suffix: part_simple_0
543     requires: triangle
544     nsize: 2
545     args: -dim 2 -cell_simplex 1 -dm_refine 0 -interpolate 0 -petscpartitioner_type simple -partition_view -dm_view ascii::ascii_info_detail
546   test:
547     suffix: part_simple_1
548     requires: triangle
549     nsize: 8
550     args: -dim 2 -cell_simplex 1 -dm_refine 1 -interpolate 1 -petscpartitioner_type simple -partition_view -dm_view ascii::ascii_info_detail
551 
552   # Parallel partitioner tests
553   test:
554     suffix: part_parmetis_0
555     requires: parmetis
556     nsize: 2
557     args: -dim 2 -cell_simplex 0 -dm_refine 1 -interpolate 1 -petscpartitioner_type parmetis -dm_view -petscpartitioner_view -test_redistribute -dm_plex_csr_via_mat {{0 1}} -dm_pre_redist_view ::load_balance -dm_post_redist_view ::load_balance -petscpartitioner_view_graph
558   test:
559     suffix: part_ptscotch_0
560     requires: ptscotch
561     nsize: 2
562     args: -dim 2 -cell_simplex 0 -dm_refine 0 -interpolate 1 -petscpartitioner_type ptscotch -petscpartitioner_view -petscpartitioner_ptscotch_strategy quality -test_redistribute -dm_plex_csr_via_mat {{0 1}} -dm_pre_redist_view ::load_balance -dm_post_redist_view ::load_balance -petscpartitioner_view_graph
563   test:
564     suffix: part_ptscotch_1
565     requires: ptscotch
566     nsize: 8
567     args: -dim 2 -cell_simplex 0 -dm_refine 1 -interpolate 1 -petscpartitioner_type ptscotch -petscpartitioner_view -petscpartitioner_ptscotch_imbalance 0.1
568 
569   # CGNS reader tests 10-11 (need to find smaller test meshes)
570   test:
571     suffix: cgns_0
572     requires: cgns
573     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/tut21.cgns -interpolate 1 -dm_view
574 
575   # Gmsh mesh reader tests
576   test:
577     suffix: gmsh_0
578     requires: !single
579     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -interpolate 1 -dm_view
580   test:
581     suffix: gmsh_1
582     requires: !single
583     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -interpolate 1 -dm_view
584   test:
585     suffix: gmsh_2
586     requires: !single
587     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin.msh -interpolate 1 -dm_view
588   test:
589     suffix: gmsh_3
590     nsize: 3
591     requires: !single
592     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -test_partition -interpolate 1 -dm_view
593   test:
594     suffix: gmsh_4
595     nsize: 3
596     requires: !single
597     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin.msh -test_partition -interpolate 1 -dm_view
598   test:
599     suffix: gmsh_5
600     requires: !single
601     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_quad.msh -interpolate 1 -dm_view
602   # TODO: it seems the mesh is not a valid gmsh (inverted cell)
603   test:
604     suffix: gmsh_6
605     requires: !single
606     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin_physnames.msh -interpolate 1 -dm_view -final_diagnostics 0
607   test:
608     suffix: gmsh_7
609     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere_bin.msh -dm_view ::ascii_info_detail -interpolate -dm_plex_check_all
610   test:
611     suffix: gmsh_8
612     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere.msh -dm_view ::ascii_info_detail -interpolate -dm_plex_check_all
613   testset:
614     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic_bin.msh -dm_view ::ascii_info_detail -interpolate -dm_plex_check_all
615     test:
616       suffix: gmsh_9
617     test:
618       suffix: gmsh_9_periodic_0
619       args: -dm_plex_gmsh_periodic 0
620   testset:
621     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_view ::ascii_info_detail -interpolate -dm_plex_check_all
622     test:
623       suffix: gmsh_10
624     test:
625       suffix: gmsh_10_periodic_0
626       args: -dm_plex_gmsh_periodic 0
627   testset:
628     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_view ::ascii_info_detail -interpolate -dm_plex_check_all -dm_refine 1
629     test:
630       suffix: gmsh_11
631     test:
632       suffix: gmsh_11_periodic_0
633       args: -dm_plex_gmsh_periodic 0
634   # TODO: it seems the mesh is not a valid gmsh (inverted cell)
635   test:
636     suffix: gmsh_12
637     nsize: 4
638     requires: !single mpiio
639     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin_physnames.msh -viewer_binary_mpiio -petscpartitioner_type simple -interpolate 1 -dm_view -final_diagnostics 0
640   test:
641     suffix: gmsh_13_hybs2t
642     nsize: 4
643     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_triquad.msh -petscpartitioner_type simple -interpolate 1 -dm_view -dm_refine 1 -dm_plex_cell_refiner tobox -dm_plex_check_all
644   test:
645     suffix: gmsh_14_ext
646     requires: !single
647     args: -ext_layers 2 -ext_thickness 1.5 -ext_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin.msh -dm_view -interpolate -dm_plex_check_all
648   test:
649     suffix: gmsh_14_ext_s2t
650     requires: !single
651     args: -ext_layers 2 -ext_thickness 1.5 -ext_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin.msh -dm_view -interpolate -dm_plex_check_all -dm_refine 1 -dm_plex_cell_refiner tobox
652   test:
653     suffix: gmsh_15_hyb3d
654     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_tetwedge.msh -dm_view -interpolate -dm_plex_check_all
655   test:
656     suffix: gmsh_15_hyb3d_vtk
657     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_tetwedge.msh -dm_view vtk: -interpolate -dm_plex_gmsh_hybrid -dm_plex_check_all
658   test:
659     suffix: gmsh_15_hyb3d_s2t
660     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_tetwedge.msh -dm_view -interpolate -dm_plex_check_all -dm_refine 1 -dm_plex_cell_refiner tobox
661   test:
662     suffix: gmsh_16_spheresurface
663     nsize : 4
664     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3 -dm_plex_check_all -dm_view -interpolate -petscpartitioner_type simple
665   test:
666     suffix: gmsh_16_spheresurface_s2t
667     nsize : 4
668     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3 -dm_refine 1 -dm_plex_cell_refiner tobox -dm_plex_check_all -dm_view -interpolate -petscpartitioner_type simple
669   test:
670     suffix: gmsh_16_spheresurface_extruded
671     nsize : 4
672     args: -ext_layers 3 -ext_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3 -dm_plex_check_all -dm_view -interpolate -petscpartitioner_type simple
673   test:
674     suffix: gmsh_16_spheresurface_extruded_s2t
675     nsize : 4
676     args: -ext_layers 3 -ext_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3 -dm_refine 1 -dm_plex_cell_refiner tobox -dm_plex_check_all -dm_view -interpolate -petscpartitioner_type simple
677   test:
678     suffix: gmsh_17_hyb3d_interp_ascii
679     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_hexwedge.msh -dm_view -interpolate -dm_plex_check_all
680   test:
681     suffix: exodus_17_hyb3d_interp_ascii
682     requires: exodusii
683     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_hexwedge.exo -dm_view -interpolate -dm_plex_check_all
684 
685   # Legacy Gmsh v22/v40 ascii/binary reader tests
686   testset:
687     output_file: output/ex1_gmsh_3d_legacy.out
688     args: -dm_view ::ascii_info_detail -interpolate -dm_plex_check_all
689     test:
690       suffix: gmsh_3d_ascii_v22
691       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii.msh2
692     test:
693       suffix: gmsh_3d_ascii_v40
694       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii.msh4
695     test:
696       suffix: gmsh_3d_binary_v22
697       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary.msh2
698     test:
699       suffix: gmsh_3d_binary_v40
700       requires: long64
701       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary.msh4
702 
703   # Gmsh v41 ascii/binary reader tests
704   testset: # 32bit mesh, sequential
705     args: -dm_view ::ascii_info_detail -interpolate -dm_plex_check_all
706     output_file: output/ex1_gmsh_3d_32.out
707     test:
708       suffix: gmsh_3d_ascii_v41_32
709       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii-32.msh
710     test:
711       suffix: gmsh_3d_binary_v41_32
712       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-32.msh
713     test:
714       suffix: gmsh_3d_binary_v41_32_mpiio
715       requires: define(PETSC_HAVE_MPIIO)
716       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-32.msh -viewer_binary_mpiio
717   testset:  # 32bit mesh, parallel
718     args:  -petscpartitioner_type simple -dm_view ::ascii_info_detail -interpolate -dm_plex_check_all
719     nsize: 2
720     output_file: output/ex1_gmsh_3d_32_np2.out
721     test:
722       suffix: gmsh_3d_ascii_v41_32_np2
723       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii-32.msh
724     test:
725       suffix: gmsh_3d_binary_v41_32_np2
726       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-32.msh
727     test:
728       suffix: gmsh_3d_binary_v41_32_np2_mpiio
729       requires: define(PETSC_HAVE_MPIIO)
730       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-32.msh -viewer_binary_mpiio
731   testset: # 64bit mesh, sequential
732     args: -dm_view ::ascii_info_detail -interpolate -dm_plex_check_all
733     output_file: output/ex1_gmsh_3d_64.out
734     test:
735       suffix: gmsh_3d_ascii_v41_64
736       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii-64.msh
737     test:
738       suffix: gmsh_3d_binary_v41_64
739       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-64.msh
740     test:
741       suffix: gmsh_3d_binary_v41_64_mpiio
742       requires: define(PETSC_HAVE_MPIIO)
743       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-64.msh -viewer_binary_mpiio
744   testset:  # 64bit mesh, parallel
745     args:  -petscpartitioner_type simple -dm_view ::ascii_info_detail -interpolate -dm_plex_check_all
746     nsize: 2
747     output_file: output/ex1_gmsh_3d_64_np2.out
748     test:
749       suffix: gmsh_3d_ascii_v41_64_np2
750       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii-64.msh
751     test:
752       suffix: gmsh_3d_binary_v41_64_np2
753       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-64.msh
754     test:
755       suffix: gmsh_3d_binary_v41_64_np2_mpiio
756       requires: define(PETSC_HAVE_MPIIO)
757       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-64.msh -viewer_binary_mpiio
758 
759   # Fluent mesh reader tests
760   # TODO: Geometry checks fail
761   test:
762     suffix: fluent_0
763     requires: !complex
764     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.cas -interpolate 1 -dm_view -final_diagnostics 0
765   test:
766     suffix: fluent_1
767     nsize: 3
768     requires: !complex
769     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.cas -interpolate 1 -test_partition -dm_view -final_diagnostics 0
770   test:
771     suffix: fluent_2
772     requires: !complex
773     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/cube_5tets_ascii.cas -interpolate 1 -dm_view -final_diagnostics 0
774   test:
775     suffix: fluent_3
776     requires: !complex
777     TODO: Fails on non-linux: fseek(), fileno() ? https://gitlab.com/petsc/petsc/merge_requests/2206#note_238166382
778     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/cube_5tets.cas -interpolate 1 -dm_view -final_diagnostics 0
779 
780   # Med mesh reader tests, including parallel file reads
781   test:
782     suffix: med_0
783     requires: med
784     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.med -interpolate 1 -dm_view
785   test:
786     suffix: med_1
787     requires: med
788     nsize: 3
789     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.med -interpolate 1 -petscpartitioner_type simple -dm_view
790   test:
791     suffix: med_2
792     requires: med
793     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/cylinder.med -interpolate 1 -dm_view
794   test:
795     suffix: med_3
796     requires: med
797     TODO: MED
798     nsize: 3
799     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/cylinder.med -interpolate 1 -petscpartitioner_type simple -dm_view
800 
801   # Test shape quality
802   test:
803     suffix: test_shape
804     requires: ctetgen
805     args: -dim 3 -interpolate -dm_refine_hierarchy 3 -dm_plex_check_all -dm_plex_check_cell_shape
806 
807   # Test simplex to tensor conversion
808   test:
809     suffix: s2t2
810     requires: triangle
811     args: -dim 2 -dm_refine 1 -interpolate -dm_plex_cell_refiner tobox -refinement_limit 0.0625 -dm_view ascii::ascii_info_detail
812 
813   test:
814     suffix: s2t3
815     requires: ctetgen
816     args: -dim 3 -dm_refine 1 -interpolate -dm_plex_cell_refiner tobox -refinement_limit 0.0625 -dm_view ascii::ascii_info_detail
817 
818   # Test domain shapes
819   test:
820     suffix: cylinder
821     args: -dim 3 -cell_simplex 0 -interpolate -domain_shape cylinder -dm_plex_check_all -dm_view
822 
823   test:
824     suffix: cylinder_per
825     args: -dim 3 -cell_simplex 0 -interpolate -domain_shape cylinder -z_periodicity periodic -dm_plex_check_all -dm_view
826 
827   test:
828     suffix: cylinder_wedge
829     args: -dim 3 -cell_simplex 0 -interpolate -cell_wedge -domain_shape cylinder -dm_view vtk: -dm_plex_check_all
830 
831   test:
832     suffix: cylinder_wedge_int
833     output_file: output/ex1_cylinder_wedge.out
834     args: -dim 3 -cell_simplex 0 -interpolate -cell_wedge -domain_shape cylinder -dm_view vtk: -dm_plex_check_all
835 
836   test:
837     suffix: box_2d
838     args: -dim 2 -cell_simplex 0 -interpolate -domain_shape box -dm_refine 2 -dm_plex_check_all -dm_view
839 
840   test:
841     suffix: box_2d_per
842     args: -dim 2 -cell_simplex 0 -interpolate -domain_shape box -dm_refine 2 -dm_plex_check_all -dm_view
843 
844   test:
845     suffix: box_2d_per_unint
846     args: -dim 2 -cell_simplex 0 -interpolate 0 -domain_shape box -domain_box_sizes 3,3 -dm_plex_check_all -dm_view ::ascii_info_detail
847 
848   test:
849     suffix: box_3d
850     args: -dim 3 -cell_simplex 0 -interpolate -domain_shape box -dm_refine 3 -dm_plex_check_all -dm_view
851 
852   test:
853     requires: triangle
854     suffix: box_wedge
855     args: -dim 3 -cell_simplex 0 -interpolate -cell_wedge -domain_shape box -dm_view vtk: -dm_plex_check_all
856 
857   testset:
858     requires: triangle
859     args: -dim 3 -cell_simplex 0 -interpolate -cell_wedge -domain_shape box -domain_box_sizes 2,3,1 -dm_view -dm_plex_check_all -dm_refine 1 -dm_plex_cell_refiner tobox
860     test:
861       suffix: box_wedge_s2t
862     test:
863       nsize: 3
864       args: -petscpartitioner_type simple
865       suffix: box_wedge_s2t_parallel
866 
867   # Test GLVis output
868   test:
869     suffix: glvis_2d_tet
870     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0 -dm_view glvis:
871 
872   test:
873     suffix: glvis_2d_tet_per
874     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary 0
875 
876   test:
877     suffix: glvis_2d_tet_per_mfem
878     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -viewer_glvis_dm_plex_enable_boundary -viewer_glvis_dm_plex_enable_mfem -dm_view glvis: -interpolate
879 
880   test:
881     suffix: glvis_2d_quad
882     args: -dim 2 -cell_simplex 0 -interpolate -domain_shape box -domain_box_sizes 3,3 -dm_view glvis:
883 
884   test:
885     suffix: glvis_2d_quad_per
886     args: -dim 2 -cell_simplex 0 -interpolate -domain_shape box -domain_box_sizes 3,3 -x_periodicity periodic -y_periodicity periodic -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary
887 
888   test:
889     suffix: glvis_2d_quad_per_mfem
890     args: -dim 2 -cell_simplex 0 -interpolate -domain_shape box -domain_box_sizes 3,3 -x_periodicity periodic -y_periodicity periodic -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary -viewer_glvis_dm_plex_enable_mfem
891 
892   test:
893     suffix: glvis_3d_tet
894     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere_bin.msh -dm_plex_gmsh_periodic 0 -dm_view glvis:
895 
896   test:
897     suffix: glvis_3d_tet_per
898     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere_bin.msh -dm_view glvis: -interpolate -viewer_glvis_dm_plex_enable_boundary
899 
900   test:
901     suffix: glvis_3d_tet_per_mfem
902     TODO: broken
903     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere_bin.msh -viewer_glvis_dm_plex_enable_mfem -dm_view glvis: -interpolate
904 
905   test:
906     suffix: glvis_3d_hex
907     args: -dim 3 -cell_simplex 0 -interpolate -domain_shape box -domain_box_sizes 3,3,3 -dm_view glvis:
908 
909   test:
910     suffix: glvis_3d_hex_per
911     args: -dim 3 -cell_simplex 0 -interpolate -domain_shape box -domain_box_sizes 3,3,3 -x_periodicity periodic -y_periodicity periodic -z_periodicity periodic -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary 0
912 
913   test:
914     suffix: glvis_3d_hex_per_mfem
915     args: -dim 3 -cell_simplex 0 -domain_shape box -domain_box_sizes 3,3,3 -x_periodicity periodic -y_periodicity periodic -z_periodicity periodic -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary -viewer_glvis_dm_plex_enable_mfem -interpolate
916 
917   # Test P4EST
918   testset:
919     requires: p4est
920     args: -interpolate -dm_view -test_p4est_seq -conv_seq_2_dm_plex_check_all -conv_seq_1_dm_forest_minimum_refinement 1
921     test:
922       suffix: p4est_periodic
923       args: -dim 2 -domain_shape box -cell_simplex 0 -x_periodicity periodic -y_periodicity periodic -domain_box_sizes 3,5 -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 2 -conv_seq_1_dm_p4est_refine_pattern hash
924     test:
925       suffix: p4est_periodic_3d
926       args: -dim 3 -domain_shape box -cell_simplex 0 -x_periodicity periodic -y_periodicity periodic -z_periodicity -domain_box_sizes 3,5,4 -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 2 -conv_seq_1_dm_p4est_refine_pattern hash
927     test:
928       suffix: p4est_gmsh_periodic
929       args: -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 1 -conv_seq_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh
930     test:
931       suffix: p4est_gmsh_surface
932       args: -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 1 -conv_seq_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3
933     test:
934       suffix: p4est_gmsh_surface_parallel
935       nsize: 2
936       args: -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 1 -conv_seq_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3 -petscpartitioner_type simple -dm_view ::load_balance
937     test:
938       suffix: p4est_hyb_2d
939       args: -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 1 -conv_seq_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_triquad.msh
940     test:
941       suffix: p4est_hyb_3d
942       args: -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 1 -conv_seq_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_tetwedge.msh
943     test:
944       requires: ctetgen
945       suffix: p4est_s2t_bugfaces_3d
946       args: -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 0 -dim 3 -domain_box_sizes 1,1 -cell_simplex
947     test:
948       suffix: p4est_bug_overlapsf
949       nsize: 3
950       args: -dim 3 -cell_simplex 0 -domain_box_sizes 2,2,1 -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 1 -conv_seq_1_dm_p4est_refine_pattern hash  -petscpartitioner_type simple
951     test:
952       suffix: p4est_redistribute
953       nsize: 3
954       args: -dim 3 -cell_simplex 0 -domain_box_sizes 2,2,1 -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 1 -conv_seq_1_dm_p4est_refine_pattern hash  -petscpartitioner_type simple -test_redistribute -dm_plex_csr_via_mat {{0 1}} -dm_view ::load_balance
955     test:
956       suffix: p4est_gmsh_s2t_3d
957       args: -conv_seq_1_dm_forest_initial_refinement 1 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh
958     test:
959       suffix: p4est_gmsh_s2t_3d_hash
960       args: -conv_seq_1_dm_forest_initial_refinement 1 -conv_seq_1_dm_forest_maximum_refinement 2 -conv_seq_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh
961     test:
962       requires: long_runtime
963       suffix: p4est_gmsh_periodic_3d
964       args: -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 1 -conv_seq_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere.msh
965 
966   testset:
967     requires: p4est
968     nsize: 6
969     args: -interpolate -test_p4est_par -conv_par_2_dm_plex_check_all -conv_par_1_dm_forest_minimum_refinement 1 -conv_par_1_dm_forest_partition_overlap 0
970     test:
971       TODO: interface cones do not conform
972       suffix: p4est_par_periodic
973       args: -dim 2 -domain_shape box -cell_simplex 0 -x_periodicity periodic -y_periodicity periodic -domain_box_sizes 3,5 -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 2 -conv_par_1_dm_p4est_refine_pattern hash
974     test:
975       TODO: interface cones do not conform
976       suffix: p4est_par_periodic_3d
977       args: -dim 3 -domain_shape box -cell_simplex 0 -x_periodicity periodic -y_periodicity periodic -z_periodicity periodic -domain_box_sizes 3,5,4 -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 2 -conv_par_1_dm_p4est_refine_pattern hash
978     test:
979       TODO: interface cones do not conform
980       suffix: p4est_par_gmsh_periodic
981       args: -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 1 -conv_par_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh
982     test:
983       suffix: p4est_par_gmsh_surface
984       args: -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 1 -conv_par_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3
985     test:
986       suffix: p4est_par_gmsh_s2t_3d
987       args: -conv_par_1_dm_forest_initial_refinement 1 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh
988     test:
989       TODO: interface cones do not conform
990       suffix: p4est_par_gmsh_s2t_3d_hash
991       args: -conv_par_1_dm_forest_initial_refinement 1 -conv_par_1_dm_forest_maximum_refinement 2 -conv_par_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh
992     test:
993       requires: long_runtime
994       suffix: p4est_par_gmsh_periodic_3d
995       args: -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 1 -conv_par_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere.msh
996 
997   testset:
998     requires: p4est
999     nsize: 6
1000     args: -interpolate -test_p4est_par -conv_par_2_dm_plex_check_all -conv_par_1_dm_forest_minimum_refinement 1 -conv_par_1_dm_forest_partition_overlap 1 -petscpartitioner_type simple
1001     test:
1002       suffix: p4est_par_ovl_periodic
1003       args: -dim 2 -domain_shape box -cell_simplex 0 -x_periodicity periodic -y_periodicity periodic -domain_box_sizes 3,5 -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 2 -conv_par_1_dm_p4est_refine_pattern hash
1004     #TODO Mesh cell 201 is inverted, vol = 0. (FVM Volume. Is it correct? -> Diagnostics disabled)
1005     test:
1006       suffix: p4est_par_ovl_periodic_3d
1007       args: -dim 3 -domain_shape box -cell_simplex 0 -x_periodicity periodic -y_periodicity periodic -z_periodicity -domain_box_sizes 3,5,4 -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 2 -conv_par_1_dm_p4est_refine_pattern hash -final_diagnostics 0
1008     test:
1009       suffix: p4est_par_ovl_gmsh_periodic
1010       args: -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 1 -conv_par_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh
1011     test:
1012       suffix: p4est_par_ovl_gmsh_surface
1013       args: -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 1 -conv_par_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3
1014     test:
1015       suffix: p4est_par_ovl_gmsh_s2t_3d
1016       args: -conv_par_1_dm_forest_initial_refinement 1 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh
1017     test:
1018       suffix: p4est_par_ovl_gmsh_s2t_3d_hash
1019       args: -conv_par_1_dm_forest_initial_refinement 1 -conv_par_1_dm_forest_maximum_refinement 2 -conv_par_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh
1020     test:
1021       requires: long_runtime
1022       suffix: p4est_par_ovl_gmsh_periodic_3d
1023       args: -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 1 -conv_par_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere.msh
1024     test:
1025       suffix: p4est_par_ovl_hyb_2d
1026       args: -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 1 -conv_par_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_triquad.msh
1027     test:
1028       suffix: p4est_par_ovl_hyb_3d
1029       args: -conv_par_1_dm_forest_initial_refinement 0 -conv_par_1_dm_forest_maximum_refinement 1 -conv_par_1_dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_tetwedge.msh
1030 
1031   test:
1032     TODO: broken
1033     requires: p4est
1034     nsize: 2
1035     suffix: p4est_bug_labels_noovl
1036     args: -interpolate -test_p4est_seq -dm_plex_check_all -dm_forest_minimum_refinement 0 -dm_forest_partition_overlap 1 -dim 2 -domain_shape box -cell_simplex 0 -domain_box_sizes 3,3 -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 2 -dm_p4est_refine_pattern hash -petscpartitioner_type simple -dm_forest_print_label_error
1037 
1038   test:
1039     requires: p4est
1040     nsize: 2
1041     suffix: p4est_bug_distribute_overlap
1042     args: -interpolate -test_p4est_seq -conv_seq_2_dm_plex_check_all -conv_seq_1_dm_forest_minimum_refinement 0 -conv_seq_1_dm_forest_partition_overlap 0 -dim 2 -domain_shape box -cell_simplex 0 -domain_box_sizes 3,3 -conv_seq_1_dm_forest_initial_refinement 0 -conv_seq_1_dm_forest_maximum_refinement 2 -conv_seq_1_dm_p4est_refine_pattern hash -petscpartitioner_type simple -overlap 1 -dm_view ::load_balance
1043     args: -dm_post_overlap_view
1044 
1045   test:
1046     suffix: glvis_2d_hyb
1047     args: -dim 2 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_triquad.msh -interpolate -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary -petscpartitioner_type simple
1048 
1049   test:
1050     suffix: glvis_3d_hyb
1051     args: -dim 3 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_tetwedge.msh -interpolate -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary -petscpartitioner_type simple
1052 
1053   test:
1054     suffix: glvis_3d_hyb_s2t
1055     args: -dim 3 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_3d_cube.msh -interpolate -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary -petscpartitioner_type simple -dm_refine 1 -dm_plex_cell_refiner tobox -dm_plex_check_all
1056 
1057   test:
1058     suffix: ref_alfeld2d_0
1059     requires: triangle
1060     args: -dim 2 -domain_shape box -cell_simplex 1 -domain_box_sizes 5 -dm_view -interpolate -dm_plex_check_all -dm_refine 1 -dm_plex_cell_refiner alfeld2d -final_diagnostics
1061   test:
1062     suffix: ref_alfeld3d_0
1063     requires: ctetgen
1064     args: -dim 3 -domain_shape box -cell_simplex 1 -domain_box_sizes 5 -dm_view -interpolate -dm_plex_check_all -dm_refine 1 -dm_plex_cell_refiner alfeld3d -final_diagnostics
1065 
1066   # Boundary layer refiners
1067   test:
1068     suffix: ref_bl_1
1069     args: -dim 1 -domain_shape box -cell_simplex 0 -domain_box_sizes 5 -dm_view -interpolate -dm_plex_check_all 0 -dm_refine 1 -dm_plex_cell_refiner boundarylayer -ext_layers 2 -final_diagnostics -dm_plex_refine_boundarylayer_splits 3 -ext_hfirst {{0 1}}
1070   test:
1071     suffix: ref_bl_2_tri
1072     requires: triangle
1073     args: -dim 2 -domain_shape box -cell_simplex 1 -domain_box_sizes 5 -dm_view -interpolate -dm_plex_check_all 0 -dm_refine 1 -dm_plex_cell_refiner boundarylayer -ext_layers 3 -final_diagnostics -dm_plex_refine_boundarylayer_splits 4 -ext_hfirst {{0 1}}
1074   test:
1075     suffix: ref_bl_3_quad
1076     args: -dim 2 -domain_shape box -cell_simplex 0 -domain_box_sizes 5 -dm_view -interpolate -dm_plex_check_all 0 -dm_refine 1 -dm_plex_cell_refiner boundarylayer -ext_layers 3 -final_diagnostics -dm_plex_refine_boundarylayer_splits 4 -ext_hfirst {{0 1}}
1077   test:
1078     suffix: ref_bl_spheresurface_extruded
1079     nsize : 4
1080     args: -ext_layers 3 -ext_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3 -dm_plex_check_all -dm_view -interpolate -petscpartitioner_type simple -ext_hfirst {{0 1}separate output} -final_diagnostics -dm_refine 1 -dm_plex_cell_refiner boundarylayer -dm_plex_refine_boundarylayer_splits 2
1081   test:
1082     suffix: ref_bl_3d_hyb
1083     nsize : 4
1084     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_3d_cube.msh -dm_plex_check_all -dm_view -interpolate -petscpartitioner_type simple -final_diagnostics -dm_refine 1 -dm_plex_cell_refiner boundarylayer -dm_plex_refine_boundarylayer_splits 4 -dm_plex_refine_boundarylayer_progression 3.1
1085 
1086   test:
1087     suffix: sphere_0
1088     args: -dim 2 -domain_shape sphere -dm_plex_check_all -dm_view ::ascii_info_detail
1089 
1090   test:
1091     suffix: sphere_1
1092     args: -dim 2 -domain_shape sphere -dm_plex_check_all -dm_refine 2 -dm_view
1093 
1094   test:
1095     suffix: sphere_2
1096     args: -dim 2 -cell_simplex 0 -domain_shape sphere -dm_plex_check_all -dm_view ::ascii_info_detail
1097 
1098   test:
1099     suffix: sphere_3
1100     args: -dim 2 -cell_simplex 0 -domain_shape sphere -dm_plex_check_all -dm_refine 2 -dm_view
1101 
1102   test:
1103     suffix: ball_0
1104     requires: ctetgen
1105     args: -dim 3 -domain_shape ball -dm_plex_check_all -dm_view
1106 
1107   test:
1108     suffix: ball_1
1109     requires: ctetgen
1110     args: -dim 3 -domain_shape ball -bd_dm_refine 2 -dm_plex_check_all -dm_view
1111 
1112 TEST*/
1113