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