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