Actual source code: ex1.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: static char help[] = "Tests various DMPlex routines to construct, refine and distribute a mesh.\n\n";

  3:  #include <petscdmplex.h>

  5: typedef enum {BOX, CYLINDER} DomainShape;
  6: enum {STAGE_LOAD, STAGE_DISTRIBUTE, STAGE_REFINE, STAGE_OVERLAP};

  8: typedef struct {
  9:   DM            dm;                /* REQUIRED in order to use SNES evaluation functions */
 10:   PetscInt      debug;             /* The debugging level */
 11:   PetscLogEvent createMeshEvent;
 12:   PetscLogStage stages[4];
 13:   /* Domain and mesh definition */
 14:   PetscInt      dim;                             /* The topological mesh dimension */
 15:   PetscBool     interpolate;                     /* Generate intermediate mesh elements */
 16:   PetscReal     refinementLimit;                 /* The largest allowable cell volume */
 17:   PetscBool     cellSimplex;                     /* Use simplices or hexes */
 18:   PetscBool     cellWedge;                       /* Use wedges */
 19:   PetscBool     simplex2tensor;                  /* Refine simplicials in hexes */
 20:   DomainShape   domainShape;                     /* Shape of the region to be meshed */
 21:   PetscInt      *domainBoxSizes;                 /* Sizes of the box mesh */
 22:   PetscReal     *domainBoxL,*domainBoxU;         /* Lower left, upper right corner of the box mesh */
 23:   DMBoundaryType periodicity[3];                 /* The domain periodicity */
 24:   char          filename[PETSC_MAX_PATH_LEN];    /* Import mesh from file */
 25:   char          bdfilename[PETSC_MAX_PATH_LEN];  /* Import mesh boundary from file */
 26:   char          extfilename[PETSC_MAX_PATH_LEN]; /* Import 2D mesh to be extruded from file */
 27:   PetscBool     testPartition;                   /* Use a fixed partitioning for testing */
 28:   PetscInt      overlap;                         /* The cell overlap to use during partitioning */
 29:   PetscBool     testShape;                       /* Test the cell shape quality */
 30:   PetscReal     extrude_thickness;               /* Thickness of extrusion */
 31:   PetscInt      extrude_layers;                  /* Layers to be extruded */
 32:   PetscBool     testp4est[2];
 33: } AppCtx;

 35: PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
 36: {
 37:   const char       *dShapes[2] = {"box", "cylinder"};
 38:   PetscInt         shape, bd, n;
 39:   static PetscInt  domainBoxSizes[3] = {1,1,1};
 40:   static PetscReal domainBoxL[3] = {0.,0.,0.};
 41:   static PetscReal domainBoxU[3] = {1.,1.,1.};
 42:   PetscBool        flg;
 43:   PetscErrorCode   ierr;

 46:   options->debug             = 0;
 47:   options->dim               = 2;
 48:   options->interpolate       = PETSC_FALSE;
 49:   options->refinementLimit   = 0.0;
 50:   options->cellSimplex       = PETSC_TRUE;
 51:   options->cellWedge         = PETSC_FALSE;
 52:   options->domainShape       = BOX;
 53:   options->domainBoxSizes    = NULL;
 54:   options->domainBoxL        = NULL;
 55:   options->domainBoxU        = NULL;
 56:   options->periodicity[0]    = DM_BOUNDARY_NONE;
 57:   options->periodicity[1]    = DM_BOUNDARY_NONE;
 58:   options->periodicity[2]    = DM_BOUNDARY_NONE;
 59:   options->filename[0]       = '\0';
 60:   options->bdfilename[0]     = '\0';
 61:   options->extfilename[0]    = '\0';
 62:   options->testPartition     = PETSC_FALSE;
 63:   options->overlap           = PETSC_FALSE;
 64:   options->testShape         = PETSC_FALSE;
 65:   options->simplex2tensor    = PETSC_FALSE;
 66:   options->extrude_layers    = 2;
 67:   options->extrude_thickness = 0.1;
 68:   options->testp4est[0]      = PETSC_FALSE;
 69:   options->testp4est[1]      = PETSC_FALSE;

 71:   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
 72:   PetscOptionsInt("-debug", "The debugging level", "ex1.c", options->debug, &options->debug, NULL);
 73:   PetscOptionsInt("-dim", "The topological mesh dimension", "ex1.c", options->dim, &options->dim, NULL);
 74:   PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex1.c", options->interpolate, &options->interpolate, NULL);
 75:   PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex1.c", options->refinementLimit, &options->refinementLimit, NULL);
 76:   PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex1.c", options->cellSimplex, &options->cellSimplex, NULL);
 77:   PetscOptionsBool("-cell_wedge", "Use wedges if true", "ex1.c", options->cellWedge, &options->cellWedge, NULL);
 78:   PetscOptionsBool("-simplex2tensor", "Refine simplicial cells in tensor product cells", "ex1.c", options->simplex2tensor, &options->simplex2tensor, NULL);
 79:   if (options->simplex2tensor) options->interpolate = PETSC_TRUE;
 80:   shape = options->domainShape;
 81:   PetscOptionsEList("-domain_shape","The shape of the domain","ex1.c", dShapes, 2, dShapes[options->domainShape], &shape, NULL);
 82:   options->domainShape = (DomainShape) shape;
 83:   PetscOptionsIntArray("-domain_box_sizes","The sizes of the box domain","ex1.c", domainBoxSizes, (n=3,&n), &flg);
 84:   if (flg) { options->domainShape = BOX; options->domainBoxSizes = domainBoxSizes;}
 85:   PetscOptionsRealArray("-domain_box_ll","Coordinates of the lower left corner of the box domain","ex1.c", domainBoxL, (n=3,&n), &flg);
 86:   if (flg) { options->domainBoxL = domainBoxL;}
 87:   PetscOptionsRealArray("-domain_box_ur","Coordinates of the upper right corner of the box domain","ex1.c", domainBoxU, (n=3,&n), &flg);
 88:   if (flg) { options->domainBoxU = domainBoxU;}
 89:   bd = options->periodicity[0];
 90:   PetscOptionsEList("-x_periodicity", "The x-boundary periodicity", "ex1.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[0]], &bd, NULL);
 91:   options->periodicity[0] = (DMBoundaryType) bd;
 92:   bd = options->periodicity[1];
 93:   PetscOptionsEList("-y_periodicity", "The y-boundary periodicity", "ex1.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[1]], &bd, NULL);
 94:   options->periodicity[1] = (DMBoundaryType) bd;
 95:   bd = options->periodicity[2];
 96:   PetscOptionsEList("-z_periodicity", "The z-boundary periodicity", "ex1.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[2]], &bd, NULL);
 97:   options->periodicity[2] = (DMBoundaryType) bd;
 98:   PetscOptionsString("-filename", "The mesh file", "ex1.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);
 99:   PetscOptionsString("-bd_filename", "The mesh boundary file", "ex1.c", options->bdfilename, options->bdfilename, PETSC_MAX_PATH_LEN, NULL);
100:   PetscOptionsString("-ext_filename", "The 2D mesh file to be extruded", "ex1.c", options->extfilename, options->extfilename, PETSC_MAX_PATH_LEN, NULL);
101:   PetscOptionsInt("-ext_layers", "The number of layers to extrude", "ex1.c", options->extrude_layers, &options->extrude_layers, NULL);
102:   PetscOptionsReal("-ext_thickness", "The thickness of the layer to be extruded", "ex1.c", options->extrude_thickness, &options->extrude_thickness, NULL);
103:   PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex1.c", options->testPartition, &options->testPartition, NULL);
104:   PetscOptionsInt("-overlap", "The cell overlap for partitioning", "ex1.c", options->overlap, &options->overlap, NULL);
105:   PetscOptionsBool("-test_shape", "Report cell shape qualities (Jacobian condition numbers)", "ex1.c", options->testShape, &options->testShape, NULL);
106:   PetscOptionsBool("-test_p4est_seq", "Test p4est with sequential base DM", "ex1.c", options->testp4est[0], &options->testp4est[0], NULL);
107:   PetscOptionsBool("-test_p4est_par", "Test p4est with parallel base DM", "ex1.c", options->testp4est[1], &options->testp4est[1], NULL);
108:   PetscOptionsEnd();

110:   PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);
111:   PetscLogStageRegister("MeshLoad",       &options->stages[STAGE_LOAD]);
112:   PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]);
113:   PetscLogStageRegister("MeshRefine",     &options->stages[STAGE_REFINE]);
114:   PetscLogStageRegister("MeshOverlap",    &options->stages[STAGE_OVERLAP]);
115:   return(0);
116: }

118: PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
119: {
120:   PetscInt       dim                  = user->dim;
121:   PetscBool      interpolate          = user->interpolate;
122:   PetscReal      refinementLimit      = user->refinementLimit;
123:   PetscBool      cellSimplex          = user->cellSimplex;
124:   PetscBool      cellWedge            = user->cellWedge;
125:   PetscBool      simplex2tensor       = user->simplex2tensor;
126:   const char    *filename             = user->filename;
127:   const char    *bdfilename           = user->bdfilename;
128:   const char    *extfilename          = user->extfilename;
129:   PetscBool      testp4est_seq        = user->testp4est[0];
130:   PetscBool      testp4est_par        = user->testp4est[1];
131:   PetscInt       triSizes_n2[2]       = {4, 4};
132:   PetscInt       triPoints_n2[8]      = {3, 5, 6, 7, 0, 1, 2, 4};
133:   PetscInt       triSizes_n8[8]       = {1, 1, 1, 1, 1, 1, 1, 1};
134:   PetscInt       triPoints_n8[8]      = {0, 1, 2, 3, 4, 5, 6, 7};
135:   PetscInt       quadSizes[2]         = {2, 2};
136:   PetscInt       quadPoints[4]        = {2, 3, 0, 1};
137:   PetscInt       gmshSizes_n3[3]      = {14, 14, 14};
138:   PetscInt       gmshPoints_n3[42]    = {1, 2,  4,  5,  9, 10, 11, 15, 16, 20, 21, 27, 28, 29,
139:                                          3, 8, 12, 30, 31, 32, 33, 34, 35, 36, 37, 38, 39, 40,
140:                                          0, 6,  7, 13, 14, 17, 18, 19, 22, 23, 24, 25, 26, 41};
141:   PetscInt       fluentSizes_n3[3]    = {50, 50, 50};
142:   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,
143:                                          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,
144:                                           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,
145:                                          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,
146:                                           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,
147:                                          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};
148:   size_t         len, bdlen, extlen;
149:   PetscMPIInt    rank, size;

153:   PetscLogEventBegin(user->createMeshEvent,0,0,0,0);
154:   MPI_Comm_rank(comm, &rank);
155:   MPI_Comm_size(comm, &size);
156:   PetscStrlen(filename, &len);
157:   PetscStrlen(bdfilename, &bdlen);
158:   PetscStrlen(extfilename, &extlen);
159:   PetscLogStagePush(user->stages[STAGE_LOAD]);
160:   if (len) {
161:     DMPlexCreateFromFile(comm, filename, interpolate, dm);
162:   } else if (bdlen) {
163:     DM boundary;

165:     DMPlexCreateFromFile(comm, bdfilename, interpolate, &boundary);
166:     DMPlexGenerate(boundary, NULL, interpolate, dm);
167:     DMDestroy(&boundary);
168:   } else if (extlen) {
169:     DM edm;

171:     DMPlexCreateFromFile(comm, extfilename, interpolate, &edm);
172:     DMPlexExtrude(edm, user->extrude_layers, user->extrude_thickness, PETSC_TRUE, interpolate, dm);
173:     DMDestroy(&edm);
174:   } else {
175:     switch (user->domainShape) {
176:     case BOX:
177:       if (cellWedge) {
178:         if (dim != 3) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Dimension must be 3 for a wedge mesh, not %D", dim);
179:         DMPlexCreateWedgeBoxMesh(comm, user->domainBoxSizes, user->domainBoxL, user->domainBoxU, user->periodicity, PETSC_FALSE, interpolate, dm);
180:       } else {
181:         DMPlexCreateBoxMesh(comm, dim, cellSimplex, user->domainBoxSizes, user->domainBoxL, user->domainBoxU, user->periodicity, interpolate, dm);
182:       }
183:       break;
184:     case CYLINDER:
185:       if (cellSimplex) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot mesh a cylinder with simplices");
186:       if (dim != 3)    SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Dimension must be 3 for a cylinder mesh, not %D", dim);
187:       if (cellWedge) {
188:         DMPlexCreateWedgeCylinderMesh(comm, 6, interpolate, dm);
189:       } else {
190:         DMPlexCreateHexCylinderMesh(comm, 3, user->periodicity[2], dm);
191:       }
192:       break;
193:     default: SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Unknown domain shape %D", user->domainShape);
194:     }
195:   }
196:   DMLocalizeCoordinates(*dm); /* needed for periodic */
197:   DMViewFromOptions(*dm,NULL,"-init_dm_view");
198:   DMGetDimension(*dm,&dim);

200:   if (testp4est_seq) {
201: #if defined(PETSC_HAVE_P4EST)
202:     DM dmConv = NULL;

204:     DMPlexSetRefinementUniform(*dm, PETSC_TRUE);
205:     DMPlexRefineSimplexToTensor(*dm, &dmConv);
206:     if (dmConv) {
207:       DMDestroy(dm);
208:       *dm  = dmConv;
209:     }
210:     user->cellSimplex = PETSC_FALSE;

212:     DMConvert(*dm,dim == 2 ? DMP4EST : DMP8EST,&dmConv);
213:     if (dmConv) {
214:       DMSetFromOptions(dmConv);
215:       DMDestroy(dm);
216:       *dm  = dmConv;
217:     }
218:     DMSetUp(*dm);
219:     DMViewFromOptions(*dm, NULL, "-conv_seq_1_dm_view");
220:     DMConvert(*dm,DMPLEX,&dmConv);
221:     if (dmConv) {
222:       DMDestroy(dm);
223:       *dm  = dmConv;
224:     }
225:     DMViewFromOptions(*dm, NULL, "-conv_seq_2_dm_view");
226: #else
227:     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with --download-p4est");
228: #endif
229:   }

231:   PetscLogStagePop();
232:   if (!testp4est_seq) {
233:     DM refinedMesh     = NULL;
234:     DM distributedMesh = NULL;

236:     if (user->testPartition) {
237:       const PetscInt  *sizes = NULL;
238:       const PetscInt  *points = NULL;
239:       PetscPartitioner part;

241:       if (!rank) {
242:         if (dim == 2 && cellSimplex && size == 2) {
243:            sizes = triSizes_n2; points = triPoints_n2;
244:         } else if (dim == 2 && cellSimplex && size == 8) {
245:           sizes = triSizes_n8; points = triPoints_n8;
246:         } else if (dim == 2 && !cellSimplex && size == 2) {
247:           sizes = quadSizes; points = quadPoints;
248:         } else if (dim == 2 && size == 3) {
249:           PetscInt Nc;

251:           DMPlexGetHeightStratum(*dm, 0, NULL, &Nc);
252:           if (Nc == 42) { /* Gmsh 3 & 4 */
253:             sizes = gmshSizes_n3; points = gmshPoints_n3;
254:           } else if (Nc == 150) { /* Fluent 1 */
255:             sizes = fluentSizes_n3; points = fluentPoints_n3;
256:           } else if (Nc == 42) { /* Med 1 */
257:           } else if (Nc == 161) { /* Med 3 */
258:           }
259:         }
260:       }
261:       DMPlexGetPartitioner(*dm, &part);
262:       PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);
263:       PetscPartitionerShellSetPartition(part, size, sizes, points);
264:     } else {
265:       PetscPartitioner part;

267:       DMPlexGetPartitioner(*dm,&part);
268:       PetscPartitionerSetFromOptions(part);
269:     }
270:     /* Distribute mesh over processes */
271:     PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]);
272:     DMPlexDistribute(*dm, 0, NULL, &distributedMesh);
273:     if (distributedMesh) {
274:       DMDestroy(dm);
275:       *dm  = distributedMesh;
276:     }
277:     PetscLogStagePop();
278:     DMViewFromOptions(*dm, NULL, "-distributed_dm_view");
279:     /* Refine mesh using a volume constraint */
280:     PetscLogStagePush(user->stages[STAGE_REFINE]);
281:     DMPlexSetRefinementUniform(*dm, PETSC_FALSE);
282:     DMPlexSetRefinementLimit(*dm, refinementLimit);
283:     DMRefine(*dm, comm, &refinedMesh);
284:     if (refinedMesh) {
285:       DMDestroy(dm);
286:       *dm  = refinedMesh;
287:     }
288:     PetscLogStagePop();
289:   }
290:   PetscLogStagePush(user->stages[STAGE_REFINE]);
291:   DMSetFromOptions(*dm);
292:   PetscLogStagePop();

294:   if (testp4est_par) {
295: #if defined(PETSC_HAVE_P4EST)
296:     DM dmConv = NULL;

298:     DMPlexSetRefinementUniform(*dm, PETSC_TRUE);
299:     DMPlexRefineSimplexToTensor(*dm, &dmConv);
300:     if (dmConv) {
301:       DMDestroy(dm);
302:       *dm  = dmConv;
303:     }
304:     user->cellSimplex = PETSC_FALSE;

306:     DMConvert(*dm,dim == 2 ? DMP4EST : DMP8EST,&dmConv);
307:     PetscObjectSetOptionsPrefix((PetscObject) dmConv, "conv_par_1_");
308:     if (dmConv) {
309:       DMSetFromOptions(dmConv);
310:       DMDestroy(dm);
311:       *dm  = dmConv;
312:     }
313:     DMSetUp(*dm);
314:     DMViewFromOptions(*dm, NULL, "-dm_view");
315:     DMConvert(*dm,DMPLEX,&dmConv);
316:     PetscObjectSetOptionsPrefix((PetscObject) dmConv, "conv_par_2_");
317:     if (dmConv) {
318:       DMSetFromOptions(dmConv);
319:       DMDestroy(dm);
320:       *dm  = dmConv;
321:     }
322:     DMViewFromOptions(*dm, NULL, "-dm_view");
323: #else
324:     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with --download-p4est");
325: #endif
326:   }

328:   if (user->overlap) {
329:     DM overlapMesh = NULL;

331:     /* Add the level-1 overlap to refined mesh */
332:     PetscLogStagePush(user->stages[STAGE_OVERLAP]);
333:     DMPlexDistributeOverlap(*dm, 1, NULL, &overlapMesh);
334:     if (overlapMesh) {
335:       DMView(overlapMesh, PETSC_VIEWER_STDOUT_WORLD);
336:       DMDestroy(dm);
337:       *dm = overlapMesh;
338:     }
339:     PetscLogStagePop();
340:   }

342:   if (simplex2tensor) {
343:     DM rdm = NULL;
344:     DMPlexSetRefinementUniform(*dm, PETSC_TRUE);
345:     DMPlexRefineSimplexToTensor(*dm, &rdm);
346:     if (rdm) {
347:       DMDestroy(dm);
348:       *dm  = rdm;
349:     }
350:     user->cellSimplex = PETSC_FALSE;
351:   }
352:   PetscObjectSetName((PetscObject) *dm, "Simplicial Mesh");
353:   DMViewFromOptions(*dm, NULL, "-dm_view");
354:   PetscLogEventEnd(user->createMeshEvent,0,0,0,0);
355:   user->dm = *dm;
356:   return(0);
357: }

359: int main(int argc, char **argv)
360: {
361:   AppCtx         user;                 /* user-defined work context */

364:   PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
365:   ProcessOptions(PETSC_COMM_WORLD, &user);
366:   CreateMesh(PETSC_COMM_WORLD, &user, &user.dm);
367:   if (user.testShape) {DMPlexCheckCellShape(user.dm, PETSC_TRUE);}
368:   DMDestroy(&user.dm);
369:   PetscFinalize();
370:   return ierr;
371: }

373: /*TEST

375:   # CTetGen 0-1
376:   test:
377:     suffix: 0
378:     requires: ctetgen
379:     args: -dim 3 -ctetgen_verbose 4 -dm_view ascii::ascii_info_detail -info -info_exclude null
380:   test:
381:     suffix: 1
382:     requires: ctetgen
383:     args: -dim 3 -ctetgen_verbose 4 -refinement_limit 0.0625 -dm_view ascii::ascii_info_detail -info -info_exclude null

385:   # 2D LaTex and ASCII output 2-9
386:   test:
387:     suffix: 2
388:     requires: triangle
389:     args: -dim 2 -dm_view ascii::ascii_latex
390:   test:
391:     suffix: 3
392:     requires: triangle
393:     args: -dim 2 -dm_refine 1 -interpolate 1 -dm_view ascii::ascii_info_detail
394:   test:
395:     suffix: 4
396:     requires: triangle
397:     nsize: 2
398:     args: -dim 2 -dm_refine 1 -interpolate 1 -test_partition -dm_view ascii::ascii_info_detail
399:   test:
400:     suffix: 5
401:     requires: triangle
402:     nsize: 2
403:     args: -dim 2 -dm_refine 1 -interpolate 1 -test_partition -dm_view ascii::ascii_latex
404:   test:
405:     suffix: 6
406:     args: -dim 2 -cell_simplex 0 -interpolate -dm_view ascii::ascii_info_detail
407:   test:
408:     suffix: 7
409:     args: -dim 2 -cell_simplex 0 -interpolate -dm_refine 1 -dm_view ascii::ascii_info_detail
410:   test:
411:     suffix: 8
412:     nsize: 2
413:     args: -dim 2 -cell_simplex 0 -interpolate -dm_refine 1 -interpolate 1 -test_partition -dm_view ascii::ascii_latex

415:   # 1D ASCII output
416:   test:
417:     suffix: 1d_0
418:     args: -dim 1 -domain_shape box -dm_view ascii::ascii_info_detail
419:   test:
420:     suffix: 1d_1
421:     args: -dim 1 -domain_shape box -dm_refine 2 -dm_view ascii::ascii_info_detail
422:   test:
423:     suffix: 1d_2
424:     args: -dim 1 -domain_box_sizes 5 -x_periodicity periodic -dm_view ascii::ascii_info_detail -test_shape

426:   # Parallel refinement tests with overlap
427:   test:
428:     suffix: 1d_refine_overlap_0
429:     nsize: 2
430:     args: -dim 1 -domain_box_sizes 4 -dm_refine 1 -overlap 0 -petscpartitioner_type simple -dm_view ascii::ascii_info_detail
431:   test:
432:     suffix: 1d_refine_overlap_1
433:     nsize: 2
434:     args: -dim 1 -domain_box_sizes 4 -dm_refine 1 -overlap 1 -petscpartitioner_type simple -dm_view ascii::ascii_info_detail
435:   test:
436:     suffix: refine_overlap_0
437:     requires: triangle
438:     nsize: 2
439:     requires: triangle
440:     args: -dim 2 -cell_simplex 1 -dm_refine 1 -interpolate 1 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail
441:   test:
442:     suffix: refine_overlap_1
443:     requires: triangle
444:     nsize: 8
445:     args: -dim 2 -cell_simplex 1 -dm_refine 1 -interpolate 1 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail

447:   # Parallel simple partitioner tests
448:   test:
449:     suffix: part_simple_0
450:     requires: triangle
451:     nsize: 2
452:     args: -dim 2 -cell_simplex 1 -dm_refine 0 -interpolate 0 -petscpartitioner_type simple -partition_view -dm_view ascii::ascii_info_detail
453:   test:
454:     suffix: part_simple_1
455:     requires: triangle
456:     nsize: 8
457:     args: -dim 2 -cell_simplex 1 -dm_refine 1 -interpolate 1 -petscpartitioner_type simple -partition_view -dm_view ascii::ascii_info_detail

459:   test:
460:     suffix: part_parmetis_0
461:     requires: parmetis
462:     nsize: 2
463:     args: -dim 2 -cell_simplex 0 -dm_refine 1 -interpolate 1 -petscpartitioner_type parmetis -dm_view -petscpartitioner_view
464:   # Parallel ptscotch partitioner tests
465:   test:
466:     suffix: part_ptscotch_0
467:     requires: ptscotch
468:     nsize: 2
469:     args: -dim 2 -cell_simplex 0 -dm_refine 0 -interpolate 0 -petscpartitioner_type ptscotch -petscpartitioner_view -petscpartitioner_ptscotch_strategy quality
470:   test:
471:     suffix: part_ptscotch_1
472:     requires: ptscotch
473:     nsize: 8
474:     args: -dim 2 -cell_simplex 0 -dm_refine 1 -interpolate 1 -petscpartitioner_type ptscotch -petscpartitioner_view -petscpartitioner_ptscotch_imbalance 0.1

476:   # CGNS reader tests 10-11 (need to find smaller test meshes)
477:   test:
478:     suffix: cgns_0
479:     requires: cgns
480:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/tut21.cgns -interpolate 1 -dm_view

482:   # Gmsh mesh reader tests
483:   test:
484:     suffix: gmsh_0
485:     requires: !single
486:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -interpolate 1 -dm_view
487:   test:
488:     suffix: gmsh_1
489:     requires: !single
490:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -interpolate 1 -dm_view
491:   test:
492:     suffix: gmsh_2
493:     requires: !single
494:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin.msh -interpolate 1 -dm_view
495:   test:
496:     suffix: gmsh_3
497:     nsize: 3
498:     requires: !single
499:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -test_partition -interpolate 1 -dm_view
500:   test:
501:     suffix: gmsh_4
502:     nsize: 3
503:     requires: !single
504:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin.msh -test_partition -interpolate 1 -dm_view
505:   test:
506:     suffix: gmsh_5
507:     requires: !single
508:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_quad.msh -interpolate 1 -dm_view
509:   test:
510:     suffix: gmsh_6
511:     requires: !single
512:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin_physnames.msh -interpolate 1 -dm_view
513:   test:
514:     suffix: gmsh_7
515:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere_bin.msh -dm_plex_gmsh_periodic -dm_view ::ascii_info_detail -interpolate -test_shape
516:   test:
517:     suffix: gmsh_8
518:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere.msh -dm_plex_gmsh_periodic -dm_view ::ascii_info_detail -interpolate -test_shape
519:   test:
520:     suffix: gmsh_9
521:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic_bin.msh -dm_plex_gmsh_periodic -dm_view ::ascii_info_detail -interpolate -test_shape
522:   test:
523:     suffix: gmsh_10
524:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic -dm_view ::ascii_info_detail -interpolate -test_shape
525:   test:
526:     suffix: gmsh_11
527:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic -dm_view ::ascii_info_detail -interpolate -test_shape -dm_refine 1
528:   test:
529:     suffix: gmsh_12
530:     nsize: 4
531:     requires: !single mpiio
532:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin_physnames.msh -viewer_binary_mpiio -petscpartitioner_type simple -interpolate 1 -dm_view
533:   test:
534:     suffix: gmsh_13_hybs2t
535:     nsize: 4
536:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_triquad.msh -petscpartitioner_type simple -interpolate 1 -dm_view -test_shape -simplex2tensor -dm_plex_gmsh_hybrid -dm_plex_check_faces unknown -dm_plex_check_skeleton unknown -dm_plex_check_symmetry
537:   test:
538:     suffix: gmsh_14_ext
539:     requires: !single
540:     args: -ext_layers 2 -ext_thickness 1.5 -ext_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin.msh -dm_view -dm_plex_check_symmetry -dm_plex_check_skeleton unknown
541:   test:
542:     suffix: gmsh_14_ext_s2t
543:     requires: !single
544:     args: -ext_layers 2 -ext_thickness 1.5 -ext_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_bin.msh -dm_view -interpolate -dm_plex_check_faces unknown -dm_plex_check_symmetry -dm_plex_check_skeleton unknown -simplex2tensor -test_shape
545:   test:
546:     suffix: gmsh_15_hyb3d
547:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_tetwedge.msh -dm_view -interpolate -dm_plex_check_faces unknown -dm_plex_check_symmetry -dm_plex_check_skeleton unknown -dm_plex_gmsh_hybrid
548:   test:
549:     suffix: gmsh_15_hyb3d_vtk
550:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_tetwedge.msh -dm_view vtk:
551:   test:
552:     suffix: gmsh_15_hyb3d_s2t
553:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_tetwedge.msh -dm_view -interpolate -dm_plex_check_faces unknown -dm_plex_check_symmetry -dm_plex_check_skeleton unknown -dm_plex_gmsh_hybrid -simplex2tensor -test_shape
554:   test:
555:     suffix: gmsh_16_spheresurface
556:     nsize : 4
557:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3 -dm_plex_check_symmetry -dm_plex_check_faces unknown -dm_plex_check_skeleton unknown -dm_view -interpolate -test_shape -petscpartitioner_type simple
558:   test:
559:     suffix: gmsh_16_spheresurface_s2t
560:     nsize : 4
561:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3 -simplex2tensor -dm_plex_check_symmetry -dm_plex_check_faces unknown -dm_plex_check_skeleton unknown -dm_view -interpolate -test_shape -petscpartitioner_type simple
562:   test:
563:     suffix: gmsh_16_spheresurface_extruded
564:     nsize : 4
565:     args: -ext_layers 3 -ext_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_hybrid -dm_plex_gmsh_spacedim 3 -dm_plex_check_symmetry -dm_plex_check_faces unknown -dm_plex_check_skeleton unknown -dm_view -interpolate -petscpartitioner_type simple
566:   test:
567:     suffix: gmsh_16_spheresurface_extruded_s2t
568:     nsize : 4
569:     args: -ext_layers 3 -ext_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3 -simplex2tensor -dm_plex_check_symmetry -dm_plex_check_faces unknown -dm_plex_check_skeleton unknown -dm_view -interpolate -test_shape -petscpartitioner_type simple

571:   # Legacy Gmsh v22/v40 ascii/binary reader tests
572:   testset:
573:     output_file: output/ex1_gmsh_3d_legacy.out
574:     args: -dm_plex_gmsh_periodic -dm_view ::ascii_info_detail -interpolate -dm_plex_check_symmetry -dm_plex_check_faces unknown -test_shape
575:     test:
576:       suffix: gmsh_3d_ascii_v22
577:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii.msh2
578:     test:
579:       suffix: gmsh_3d_ascii_v40
580:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii.msh4
581:     test:
582:       suffix: gmsh_3d_binary_v22
583:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary.msh2
584:     test:
585:       suffix: gmsh_3d_binary_v40
586:       requires: long64
587:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary.msh4

589:   # Gmsh v41 ascii/binary reader tests
590:   testset: # 32bit mesh, sequential
591:     args: -dm_plex_gmsh_periodic -dm_view ::ascii_info_detail -interpolate -dm_plex_check_symmetry -dm_plex_check_faces unknown -test_shape
592:     output_file: output/ex1_gmsh_3d_32.out
593:     test:
594:       suffix: gmsh_3d_ascii_v41_32
595:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii-32.msh
596:     test:
597:       suffix: gmsh_3d_binary_v41_32
598:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-32.msh
599:     test:
600:       suffix: gmsh_3d_binary_v41_32_mpiio
601:       requires: define(PETSC_HAVE_MPIIO) !define(PETSC_HAVE_LIBMSMPI)
602:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-32.msh -viewer_binary_mpiio
603:   testset:  # 32bit mesh, parallel
604:     args:  -petscpartitioner_type simple -dm_plex_gmsh_periodic -dm_view ::ascii_info_detail -interpolate -dm_plex_check_symmetry -dm_plex_check_faces unknown -test_shape
605:     nsize: 2
606:     output_file: output/ex1_gmsh_3d_32_np2.out
607:     test:
608:       suffix: gmsh_3d_ascii_v41_32_np2
609:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii-32.msh
610:     test:
611:       suffix: gmsh_3d_binary_v41_32_np2
612:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-32.msh
613:     test:
614:       suffix: gmsh_3d_binary_v41_32_np2_mpiio
615:       requires: define(PETSC_HAVE_MPIIO) !define(PETSC_HAVE_LIBMSMPI)
616:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-32.msh -viewer_binary_mpiio
617:   testset: # 64bit mesh, sequential
618:     args: -dm_plex_gmsh_periodic -dm_view ::ascii_info_detail -interpolate -dm_plex_check_symmetry -dm_plex_check_faces unknown -test_shape
619:     output_file: output/ex1_gmsh_3d_64.out
620:     test:
621:       suffix: gmsh_3d_ascii_v41_64
622:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii-64.msh
623:     test:
624:       suffix: gmsh_3d_binary_v41_64
625:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-64.msh
626:     test:
627:       suffix: gmsh_3d_binary_v41_64_mpiio
628:       requires: define(PETSC_HAVE_MPIIO) !define(PETSC_HAVE_LIBMSMPI)
629:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-64.msh -viewer_binary_mpiio
630:   testset:  # 64bit mesh, parallel
631:     args:  -petscpartitioner_type simple -dm_plex_gmsh_periodic -dm_view ::ascii_info_detail -interpolate -dm_plex_check_symmetry -dm_plex_check_faces unknown -test_shape
632:     nsize: 2
633:     output_file: output/ex1_gmsh_3d_64_np2.out
634:     test:
635:       suffix: gmsh_3d_ascii_v41_64_np2
636:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-ascii-64.msh
637:     test:
638:       suffix: gmsh_3d_binary_v41_64_np2
639:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-64.msh
640:     test:
641:       suffix: gmsh_3d_binary_v41_64_np2_mpiio
642:       requires: define(PETSC_HAVE_MPIIO) !define(PETSC_HAVE_LIBMSMPI)
643:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/gmsh-3d-binary-64.msh -viewer_binary_mpiio

645:   # Fluent mesh reader tests
646:   test:
647:     suffix: fluent_0
648:     requires: !complex
649:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.cas -interpolate 1 -dm_view
650:   test:
651:     suffix: fluent_1
652:     nsize: 3
653:     requires: !complex
654:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.cas -interpolate 1 -test_partition -dm_view
655:   test:
656:     suffix: fluent_2
657:     requires: !complex
658:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/cube_5tets_ascii.cas -interpolate 1 -dm_view
659:   test:
660:     suffix: fluent_3
661:     requires: !complex
662:     TODO: broken
663:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/cube_5tets.cas -interpolate 1 -dm_view

665:   # Med mesh reader tests, including parallel file reads
666:   test:
667:     suffix: med_0
668:     requires: med
669:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.med -interpolate 1 -dm_view
670:   test:
671:     suffix: med_1
672:     requires: med
673:     nsize: 3
674:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.med -interpolate 1 -petscpartitioner_type simple -dm_view
675:   test:
676:     suffix: med_2
677:     requires: med
678:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/cylinder.med -interpolate 1 -dm_view
679:   test:
680:     suffix: med_3
681:     requires: med
682:     nsize: 3
683:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/cylinder.med -interpolate 1 -petscpartitioner_type simple -dm_view

685:   # Test shape quality
686:   test:
687:     suffix: test_shape
688:     requires: ctetgen
689:     args: -dim 3 -interpolate -dm_refine_hierarchy 3 -test_shape

691:   # Test simplex to tensor conversion
692:   test:
693:     suffix: s2t2
694:     requires: triangle
695:     args: -dim 2 -simplex2tensor -refinement_limit 0.0625 -dm_view ascii::ascii_info_detail

697:   test:
698:     suffix: s2t3
699:     requires: ctetgen
700:     args: -dim 3 -simplex2tensor -refinement_limit 0.0625 -dm_view ascii::ascii_info_detail

702:   # Test domain shapes
703:   test:
704:     suffix: cylinder
705:     args: -dim 3 -cell_simplex 0 -interpolate -domain_shape cylinder -test_shape -dm_view

707:   test:
708:     suffix: cylinder_per
709:     args: -dim 3 -cell_simplex 0 -interpolate -domain_shape cylinder -z_periodicity periodic -test_shape -dm_view

711:   test:
712:     suffix: cylinder_wedge
713:     args: -dim 3 -cell_simplex 0 -interpolate 0 -cell_wedge -domain_shape cylinder -dm_view vtk: -dm_plex_check_symmetry -dm_plex_check_faces tensor -dm_plex_check_skeleton tensor

715:   test:
716:     suffix: cylinder_wedge_int
717:     output_file: output/ex1_cylinder_wedge.out
718:     args: -dim 3 -cell_simplex 0 -interpolate -cell_wedge -domain_shape cylinder -dm_view vtk: -dm_plex_check_symmetry -dm_plex_check_faces tensor -dm_plex_check_skeleton tensor

720:   test:
721:     suffix: box_2d
722:     args: -dim 2 -cell_simplex 0 -interpolate -domain_shape box -dm_refine 2 -test_shape -dm_view

724:   test:
725:     suffix: box_2d_per
726:     args: -dim 2 -cell_simplex 0 -interpolate -domain_shape box -dm_refine 2 -test_shape -dm_view

728:   test:
729:     suffix: box_2d_per_unint
730:     args: -dim 2 -cell_simplex 0 -interpolate 0 -domain_shape box -domain_box_sizes 3,3 -test_shape -dm_view ::ascii_info_detail

732:   test:
733:     suffix: box_3d
734:     args: -dim 3 -cell_simplex 0 -interpolate -domain_shape box -dm_refine 3 -test_shape -dm_view

736:   test:
737:     requires: triangle
738:     suffix: box_wedge
739:     args: -dim 3 -cell_simplex 0 -interpolate -cell_wedge -domain_shape box -dm_view vtk: -dm_plex_check_symmetry -dm_plex_check_faces unknown -dm_plex_check_skeleton unknown

741:   testset:
742:     requires: triangle
743:     args: -dim 3 -cell_simplex 0 -interpolate -cell_wedge -domain_shape box -domain_box_sizes 2,3,1 -dm_view -dm_plex_check_symmetry -dm_plex_check_faces unknown -dm_plex_check_skeleton unknown -simplex2tensor -test_shape
744:     test:
745:       suffix: box_wedge_s2t
746:     test:
747:       nsize: 3
748:       args: -petscpartitioner_type simple
749:       suffix: box_wedge_s2t_parallel

751:   # Test GLVis output
752:   test:
753:     suffix: glvis_2d_tet
754:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_view glvis:

756:   test:
757:     suffix: glvis_2d_tet_per
758:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary 0

760:   test:
761:     suffix: glvis_2d_tet_per_mfem
762:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic -viewer_glvis_dm_plex_enable_boundary -viewer_glvis_dm_plex_enable_mfem -dm_view glvis: -interpolate

764:   test:
765:     suffix: glvis_2d_quad
766:     args: -dim 2 -cell_simplex 0 -interpolate -domain_shape box -domain_box_sizes 3,3 -dm_view glvis:

768:   test:
769:     suffix: glvis_2d_quad_per
770:     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

772:   test:
773:     suffix: glvis_2d_quad_per_mfem
774:     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

776:   test:
777:     suffix: glvis_3d_tet
778:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere_bin.msh -dm_view glvis:

780:   test:
781:     suffix: glvis_3d_tet_per
782:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere_bin.msh -dm_plex_gmsh_periodic -dm_view glvis: -interpolate -viewer_glvis_dm_plex_enable_boundary

784:   test:
785:     suffix: glvis_3d_tet_per_mfem
786:     TODO: broken
787:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere_bin.msh -dm_plex_gmsh_periodic -viewer_glvis_dm_plex_enable_mfem -dm_view glvis: -interpolate

789:   test:
790:     suffix: glvis_3d_hex
791:     args: -dim 3 -cell_simplex 0 -domain_shape box -domain_box_sizes 3,3,3 -dm_view glvis:

793:   test:
794:     suffix: glvis_3d_hex_per
795:     args: -dim 3 -cell_simplex 0 -domain_shape box -domain_box_sizes 3,3,3 -x_periodicity periodic -y_periodicity periodic -z_periodicity periodic -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary 0

797:   test:
798:     suffix: glvis_3d_hex_per_mfem
799:     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

801:   # Test P4EST
802:   testset:
803:     requires: p4est
804:     args: -interpolate -dm_view -test_p4est_seq -test_shape -dm_plex_check_symmetry -dm_plex_check_skeleton unknown -dm_plex_check_faces unknown -dm_forest_minimum_refinement 1
805:     test:
806:       suffix: p4est_periodic
807:       args: -dim 2 -domain_shape box -cell_simplex 0 -x_periodicity periodic -y_periodicity periodic -domain_box_sizes 3,5 -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 2 -dm_p4est_refine_pattern hash
808:     test:
809:       suffix: p4est_periodic_3d
810:       args: -dim 3 -domain_shape box -cell_simplex 0 -x_periodicity periodic -y_periodicity periodic -z_periodicity -domain_box_sizes 3,5,4 -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 2 -dm_p4est_refine_pattern hash
811:     test:
812:       suffix: p4est_gmsh_periodic
813:       args: -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 1 -dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic
814:     test:
815:       suffix: p4est_gmsh_surface
816:       args: -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 1 -dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/surfacesphere_bin.msh -dm_plex_gmsh_spacedim 3
817:     test:
818:       suffix: p4est_hyb_2d
819:       args: -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 1 -dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_triquad.msh -dm_plex_gmsh_hybrid
820:     test:
821:       suffix: p4est_hyb_3d
822:       args: -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 1 -dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_tetwedge.msh -dm_plex_gmsh_hybrid
823:     test:
824:       requires: ctetgen
825:       suffix: p4est_s2t_bugfaces_3d
826:       args: -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 0 -dim 3 -domain_box_sizes 1,1 -cell_simplex
827:     test:
828:       suffix: p4est_bug_overlapsf
829:       nsize: 3
830:       args: -dim 3 -cell_simplex 0 -domain_box_sizes 2,2,1 -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 1 -dm_p4est_refine_pattern hash  -petscpartitioner_type simple
831:     test:
832:       suffix: p4est_gmsh_s2t_3d
833:       args: -dm_forest_initial_refinement 1 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh
834:     test:
835:       suffix: p4est_gmsh_s2t_3d_hash
836:       args: -dm_forest_initial_refinement 1 -dm_forest_maximum_refinement 2 -dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh
837:     test:
838:       requires: long_runtime
839:       suffix: p4est_gmsh_periodic_3d
840:       args: -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 1 -dm_p4est_refine_pattern hash -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/mesh-3d-box-innersphere.msh -dm_plex_gmsh_periodic

842:   testset:
843:     requires: p4est
844:     nsize: 6
845:     args: -interpolate -test_p4est_par -test_shape -conv_par_2_dm_plex_check_symmetry -conv_par_2_dm_plex_check_skeleton tensor -conv_par_2_dm_plex_check_faces tensor -conv_par_1_dm_forest_minimum_refinement 1 -conv_par_1_dm_forest_partition_overlap 0
846:     test:
847:       suffix: p4est_par_periodic
848:       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
849:     test:
850:       suffix: p4est_par_periodic_3d
851:       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
852:     test:
853:       suffix: p4est_par_gmsh_periodic
854:       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 -dm_plex_gmsh_periodic
855:     test:
856:       suffix: p4est_par_gmsh_surface
857:       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
858:     test:
859:       suffix: p4est_par_gmsh_s2t_3d
860:       args: -conv_par_1_dm_forest_initial_refinement 1 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh
861:     test:
862:       suffix: p4est_par_gmsh_s2t_3d_hash
863:       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
864:     test:
865:       requires: long_runtime
866:       suffix: p4est_par_gmsh_periodic_3d
867:       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 -dm_plex_gmsh_periodic

869:   testset:
870:     requires: p4est
871:     nsize: 6
872:     args: -interpolate -test_p4est_par -test_shape -conv_par_2_dm_plex_check_symmetry -conv_par_2_dm_plex_check_skeleton tensor -conv_par_2_dm_plex_check_faces tensor -conv_par_1_dm_forest_minimum_refinement 1 -conv_par_1_dm_forest_partition_overlap 1 -petscpartitioner_type simple
873:     test:
874:       suffix: p4est_par_ovl_periodic
875:       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
876:     test:
877:       suffix: p4est_par_ovl_periodic_3d
878:       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
879:     test:
880:       suffix: p4est_par_ovl_gmsh_periodic
881:       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 -dm_plex_gmsh_periodic
882:     test:
883:       suffix: p4est_par_ovl_gmsh_surface
884:       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
885:     test:
886:       suffix: p4est_par_ovl_gmsh_s2t_3d
887:       args: -conv_par_1_dm_forest_initial_refinement 1 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh
888:     test:
889:       suffix: p4est_par_ovl_gmsh_s2t_3d_hash
890:       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
891:     test:
892:       requires: long_runtime
893:       suffix: p4est_par_ovl_gmsh_periodic_3d
894:       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 -dm_plex_gmsh_periodic
895:     test:
896:       suffix: p4est_par_ovl_hyb_2d
897:       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 -dm_plex_gmsh_hybrid
898:     test:
899:       suffix: p4est_par_ovl_hyb_3d
900:       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 -dm_plex_gmsh_hybrid

902:   test:
903:     TODO: broken
904:     requires: p4est
905:     nsize: 2
906:     suffix: p4est_bug_labels_noovl
907:     args: -interpolate -test_p4est_seq -test_shape -dm_plex_check_symmetry -dm_plex_check_skeleton unknown -dm_plex_check_faces unknown -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

909:   test:
910:     TODO: broken #probably this can be drop once we fix this problem
911:     requires: p4est
912:     nsize: 21
913:     suffix: p4est_bug_labels
914:     args: -dim 2 -interpolate -test_p4est_seq -test_shape -dm_plex_check_symmetry -dm_plex_check_skeleton unknown -dm_plex_check_faces unknown -cell_simplex 0 -domain_box_sizes 5,4,5 -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash  -petscpartitioner_type simple -dm_forest_print_label_error

916:   test:
917:     TODO: broken
918:     requires: p4est
919:     nsize: 2
920:     suffix: p4est_bug_distribute_overlap
921:     args: -interpolate -test_p4est_seq -test_shape -dm_plex_check_symmetry -dm_plex_check_skeleton unknown -dm_plex_check_faces unknown -dm_forest_minimum_refinement 0 -dm_forest_partition_overlap 0 -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 -overlap 1

923:   test:
924:     suffix: glvis_2d_hyb
925:     args: -dim 2 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_triquad.msh -dm_plex_gmsh_hybrid -interpolate -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary -petscpartitioner_type simple

927:   test:
928:     suffix: glvis_3d_hyb
929:     args: -dim 3 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/hybrid_tetwedge.msh -dm_plex_gmsh_hybrid -interpolate -dm_view glvis: -viewer_glvis_dm_plex_enable_boundary -petscpartitioner_type simple

931: TEST*/