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