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