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