Actual source code: ex2.c

  1: static char help[] = "Create a mesh, refine and coarsen simultaneously, and transfer a field\n\n";

  3: #include <petscds.h>
  4: #include <petscdmplex.h>
  5: #include <petscdmforest.h>
  6: #include <petscoptions.h>

  8: static PetscErrorCode AddIdentityLabel(DM dm)
  9: {
 10:   PetscInt       pStart,pEnd,p;

 14:   DMCreateLabel(dm, "identity");
 15:   DMPlexGetChart(dm, &pStart, &pEnd);
 16:   for (p = pStart; p < pEnd; p++) {DMSetLabelValue(dm, "identity", p, p);}
 17:   return(0);
 18: }

 20: static PetscErrorCode CreateAdaptivityLabel(DM forest,DMLabel *adaptLabel)
 21: {
 22:   DMLabel        identLabel;
 23:   PetscInt       cStart, cEnd, c;

 27:   DMLabelCreate(PETSC_COMM_SELF,"adapt",adaptLabel);
 28:   DMLabelSetDefaultValue(*adaptLabel,DM_ADAPT_COARSEN);
 29:   DMGetLabel(forest,"identity",&identLabel);
 30:   DMForestGetCellChart(forest,&cStart,&cEnd);
 31:   for (c = cStart; c < cEnd; c++) {
 32:     PetscInt basePoint;

 34:     DMLabelGetValue(identLabel,c,&basePoint);
 35:     if (!basePoint) {DMLabelSetValue(*adaptLabel,c,DM_ADAPT_REFINE);}
 36:   }
 37:   return(0);
 38: }

 40: static PetscErrorCode LinearFunction(PetscInt dim,PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
 41: {
 43:   u[0] = (x[0] * 2.0 + 1.) + (x[1] * 20.0 + 10.) + ((dim == 3) ? (x[2] * 200.0 + 100.) : 0.);
 44:   return(0);
 45: }

 47: static PetscErrorCode MultiaffineFunction(PetscInt dim,PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
 48: {
 50:   u[0] = (x[0] * 1.0 + 2.0) * (x[1] * 3.0 - 4.0) * ((dim == 3) ? (x[2] * 5.0 + 6.0) : 1.);
 51:   return(0);
 52: }

 54: static PetscErrorCode CoordsFunction(PetscInt dim,PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
 55: {
 56:   PetscInt f;

 59:   for (f=0;f<Nf;f++) u[f] = x[f];
 60:   return(0);
 61: }

 63: typedef struct _bc_func_ctx
 64: {
 65:   PetscErrorCode (*func) (PetscInt,PetscReal,const PetscReal [], PetscInt, PetscScalar [], void *);
 66:   PetscInt dim;
 67:   PetscInt Nf;
 68:   void *ctx;
 69: }
 70: bc_func_ctx;

 72: static PetscErrorCode bc_func_fv (PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
 73: {
 74:   bc_func_ctx    *bcCtx;

 78:   bcCtx = (bc_func_ctx *) ctx;
 79:   (bcCtx->func)(bcCtx->dim,time,c,bcCtx->Nf,xG,bcCtx->ctx);
 80:   return(0);
 81: }

 83: static PetscErrorCode IdentifyBadPoints (DM dm, Vec vec, PetscReal tol)
 84: {
 85:   DM             dmplex;
 86:   PetscInt       p, pStart, pEnd, maxDof;
 87:   Vec            vecLocal;
 88:   DMLabel        depthLabel;
 89:   PetscSection   section;

 93:   DMCreateLocalVector(dm, &vecLocal);
 94:   DMGlobalToLocalBegin(dm, vec, INSERT_VALUES, vecLocal);
 95:   DMGlobalToLocalEnd(dm, vec, INSERT_VALUES, vecLocal);
 96:   DMConvert(dm ,DMPLEX, &dmplex);
 97:   DMPlexGetChart(dmplex, &pStart, &pEnd);
 98:   DMPlexGetDepthLabel(dmplex, &depthLabel);
 99:   DMGetLocalSection(dmplex, &section);
100:   PetscSectionGetMaxDof(section, &maxDof);
101:   for (p = pStart; p < pEnd; p++) {
102:     PetscInt     s, c, cSize, parent, childID, numChildren;
103:     PetscInt     cl, closureSize, *closure = NULL;
104:     PetscScalar *values = NULL;
105:     PetscBool    bad = PETSC_FALSE;

107:     VecGetValuesSection(vecLocal, section, p, &values);
108:     PetscSectionGetDof(section, p, &cSize);
109:     for (c = 0; c < cSize; c++) {
110:       PetscReal absDiff = PetscAbsScalar(values[c]);
111:       if (absDiff > tol) {bad = PETSC_TRUE; break;}
112:     }
113:     if (!bad) continue;
114:     PetscPrintf(PETSC_COMM_SELF, "Bad point %D\n", p);
115:     DMLabelGetValue(depthLabel, p, &s);
116:     PetscPrintf(PETSC_COMM_SELF, "  Depth %D\n", s);
117:     DMPlexGetTransitiveClosure(dmplex, p, PETSC_TRUE, &closureSize, &closure);
118:     for (cl = 0; cl < closureSize; cl++) {
119:       PetscInt cp = closure[2 * cl];
120:       DMPlexGetTreeParent(dmplex, cp, &parent, &childID);
121:       if (parent != cp) {
122:         PetscPrintf(PETSC_COMM_SELF, "  Closure point %D (%D) child of %D (ID %D)\n", cl, cp, parent, childID);
123:       }
124:       DMPlexGetTreeChildren(dmplex, cp, &numChildren, NULL);
125:       if (numChildren) {
126:         PetscPrintf(PETSC_COMM_SELF, "  Closure point %D (%D) is parent\n", cl, cp);
127:       }
128:     }
129:     DMPlexRestoreTransitiveClosure(dmplex, p, PETSC_TRUE, &closureSize, &closure);
130:     for (c = 0; c < cSize; c++) {
131:       PetscReal absDiff = PetscAbsScalar(values[c]);
132:       if (absDiff > tol) {
133:         PetscPrintf(PETSC_COMM_SELF, "  Bad dof %D\n", c);
134:       }
135:     }
136:   }
137:   DMDestroy(&dmplex);
138:   VecDestroy(&vecLocal);
139:   return(0);
140: }

142: int main(int argc, char **argv)
143: {
144:   MPI_Comm       comm;
145:   DM             base, preForest, postForest;
146:   PetscInt       dim, Nf = 1;
147:   PetscInt       step, adaptSteps = 1;
148:   PetscInt       preCount, postCount;
149:   Vec            preVec, postVecTransfer, postVecExact;
150:   PetscErrorCode (*funcs[1]) (PetscInt,PetscReal,const PetscReal [],PetscInt,PetscScalar [], void *) = {MultiaffineFunction};
151:   void           *ctxs[1] = {NULL};
152:   PetscReal      diff, tol = PETSC_SMALL;
153:   PetscBool      linear = PETSC_FALSE;
154:   PetscBool      coords = PETSC_FALSE;
155:   PetscBool      useFV = PETSC_FALSE;
156:   PetscBool      conv = PETSC_FALSE;
157:   PetscBool      transfer_from_base[2] = {PETSC_TRUE,PETSC_FALSE};
158:   PetscBool      use_bcs = PETSC_TRUE;
159:   bc_func_ctx    bcCtx;
160:   DMLabel        adaptLabel;

163:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
164:   comm = PETSC_COMM_WORLD;
165:   PetscOptionsBegin(comm, "", "DMForestTransferVec() Test Options", "DMFOREST");
166:   PetscOptionsBool("-linear","Transfer a simple linear function", "ex2.c", linear, &linear, NULL);
167:   PetscOptionsBool("-coords","Transfer a simple coordinate function", "ex2.c", coords, &coords, NULL);
168:   PetscOptionsBool("-use_fv","Use a finite volume approximation", "ex2.c", useFV, &useFV, NULL);
169:   PetscOptionsBool("-test_convert","Test conversion to DMPLEX",NULL,conv,&conv,NULL);
170:   PetscOptionsBool("-transfer_from_base","Transfer a vector from base DM to DMForest", "ex2.c", transfer_from_base[0], &transfer_from_base[0], NULL);
171:   transfer_from_base[1] = transfer_from_base[0];
172:   PetscOptionsBool("-transfer_from_base_steps","Transfer a vector from base DM to the latest DMForest after the adaptivity steps", "ex2.c", transfer_from_base[1], &transfer_from_base[1], NULL);
173:   PetscOptionsBool("-use_bcs","Use dirichlet boundary conditions", "ex2.c", use_bcs, &use_bcs, NULL);
174:   PetscOptionsBoundedInt("-adapt_steps","Number of adaptivity steps", "ex2.c", adaptSteps, &adaptSteps, NULL,0);
175:   PetscOptionsEnd();

177:   tol = PetscMax(1.e-10,tol); /* XXX fix for quadruple precision -> why do I need to do this? */

179:   /* the base mesh */
180:   DMCreate(comm, &base);
181:   DMSetType(base, DMPLEX);
182:   DMSetFromOptions(base);

184:   AddIdentityLabel(base);
185:   DMGetDimension(base, &dim);

187:   if (linear) {
188:     funcs[0] = LinearFunction;
189:   }
190:   if (coords) {
191:     funcs[0] = CoordsFunction;
192:     Nf = dim;
193:   }

195:   bcCtx.func = funcs[0];
196:   bcCtx.dim  = dim;
197:   bcCtx.Nf   = Nf;
198:   bcCtx.ctx  = NULL;

200:   if (useFV) {
201:     PetscFV      fv;
202:     PetscLimiter limiter;
203:     DM           baseFV;

205:     DMPlexConstructGhostCells(base,NULL,NULL,&baseFV);
206:     DMViewFromOptions(baseFV, NULL, "-fv_dm_view");
207:     DMDestroy(&base);
208:     base = baseFV;
209:     PetscFVCreate(comm, &fv);
210:     PetscFVSetSpatialDimension(fv,dim);
211:     PetscFVSetType(fv,PETSCFVLEASTSQUARES);
212:     PetscFVSetNumComponents(fv,Nf);
213:     PetscLimiterCreate(comm,&limiter);
214:     PetscLimiterSetType(limiter,PETSCLIMITERNONE);
215:     PetscFVSetLimiter(fv,limiter);
216:     PetscLimiterDestroy(&limiter);
217:     PetscFVSetFromOptions(fv);
218:     DMSetField(base,0,NULL,(PetscObject)fv);
219:     PetscFVDestroy(&fv);
220:   } else {
221:     PetscFE fe;

223:     PetscFECreateDefault(comm,dim,Nf,PETSC_FALSE,NULL,PETSC_DEFAULT,&fe);
224:     DMSetField(base,0,NULL,(PetscObject)fe);
225:     PetscFEDestroy(&fe);
226:   }
227:   DMCreateDS(base);

229:   if (use_bcs) {
230:     PetscInt ids[] = {1, 2, 3, 4, 5, 6};
231:     DMLabel  label;

233:     DMGetLabel(base, "marker", &label);
234:     DMAddBoundary(base,DM_BC_ESSENTIAL, "bc", label, 2 * dim, ids, 0, 0, NULL, useFV ? (void(*)(void)) bc_func_fv : (void(*)(void)) funcs[0], NULL, useFV ? (void *) &bcCtx : NULL, NULL);
235:   }
236:   DMViewFromOptions(base,NULL,"-dm_base_view");

238:   /* the pre adaptivity forest */
239:   DMCreate(comm,&preForest);
240:   DMSetType(preForest,(dim == 2) ? DMP4EST : DMP8EST);
241:   DMCopyDisc(base,preForest);
242:   DMForestSetBaseDM(preForest,base);
243:   DMForestSetMinimumRefinement(preForest,0);
244:   DMForestSetInitialRefinement(preForest,1);
245:   DMSetFromOptions(preForest);
246:   DMSetUp(preForest);
247:   DMViewFromOptions(preForest,NULL,"-dm_pre_view");

249:   /* the pre adaptivity field */
250:   DMCreateGlobalVector(preForest,&preVec);
251:   DMProjectFunction(preForest,0.,funcs,ctxs,INSERT_VALUES,preVec);
252:   VecViewFromOptions(preVec,NULL,"-vec_pre_view");

254:   /* communicate between base and pre adaptivity forest */
255:   if (transfer_from_base[0]) {
256:     Vec baseVec, baseVecMapped;

258:     DMGetGlobalVector(base,&baseVec);
259:     DMProjectFunction(base,0.,funcs,ctxs,INSERT_VALUES,baseVec);
260:     PetscObjectSetName((PetscObject)baseVec,"Function Base");
261:     VecViewFromOptions(baseVec,NULL,"-vec_base_view");

263:     DMGetGlobalVector(preForest,&baseVecMapped);
264:     DMForestTransferVecFromBase(preForest,baseVec,baseVecMapped);
265:     VecViewFromOptions(baseVecMapped,NULL,"-vec_map_base_view");

267:     /* compare */
268:     VecAXPY(baseVecMapped,-1.,preVec);
269:     VecViewFromOptions(baseVecMapped,NULL,"-vec_map_diff_view");
270:     VecNorm(baseVecMapped,NORM_2,&diff);

272:     /* output */
273:     if (diff < tol) {
274:       PetscPrintf(comm,"DMForestTransferVecFromBase() passes.\n");
275:     } else {
276:       PetscPrintf(comm,"DMForestTransferVecFromBase() fails with error %g and tolerance %g\n",(double)diff,(double)tol);
277:     }

279:     DMRestoreGlobalVector(base,&baseVec);
280:     DMRestoreGlobalVector(preForest,&baseVecMapped);
281:   }

283:   for (step = 0; step < adaptSteps; ++step) {

285:     if (!transfer_from_base[1]) {
286:       PetscObjectGetReference((PetscObject)preForest,&preCount);
287:     }

289:     /* adapt */
290:     CreateAdaptivityLabel(preForest,&adaptLabel);
291:     DMForestTemplate(preForest,comm,&postForest);
292:     if (step) { DMForestSetAdaptivityLabel(postForest,adaptLabel); }
293:     DMLabelDestroy(&adaptLabel);
294:     DMSetUp(postForest);
295:     DMViewFromOptions(postForest,NULL,"-dm_post_view");

297:     /* transfer */
298:     DMCreateGlobalVector(postForest,&postVecTransfer);
299:     DMForestTransferVec(preForest,preVec,postForest,postVecTransfer,PETSC_TRUE,0.0);
300:     VecViewFromOptions(postVecTransfer,NULL,"-vec_post_transfer_view");

302:     /* the exact post adaptivity field */
303:     DMCreateGlobalVector(postForest,&postVecExact);
304:     DMProjectFunction(postForest,0.,funcs,ctxs,INSERT_VALUES,postVecExact);
305:     VecViewFromOptions(postVecExact,NULL,"-vec_post_exact_view");

307:     /* compare */
308:     VecAXPY(postVecExact,-1.,postVecTransfer);
309:     VecViewFromOptions(postVecExact,NULL,"-vec_diff_view");
310:     VecNorm(postVecExact,NORM_2,&diff);

312:     /* output */
313:     if (diff < tol) {
314:       PetscPrintf(comm,"DMForestTransferVec() passes.\n");
315:     } else {
316:       PetscPrintf(comm,"DMForestTransferVec() fails with error %g and tolerance %g\n",(double)diff,(double)tol);
317:       IdentifyBadPoints(postForest, postVecExact, tol);
318:     }
319:     VecDestroy(&postVecExact);

321:     /* disconnect preForest from postForest if we don't test the transfer throughout the entire refinement process */
322:     if (!transfer_from_base[1]) {
323:       DMForestSetAdaptivityForest(postForest,NULL);
324:       PetscObjectGetReference((PetscObject)preForest,&postCount);
325:       if (postCount != preCount) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Adaptation not memory neutral: reference count increase from %d to %d\n",preCount,postCount);
326:     }

328:     if (conv) {
329:       DM dmConv;

331:       DMConvert(postForest,DMPLEX,&dmConv);
332:       DMViewFromOptions(dmConv,NULL,"-dm_conv_view");
333:       DMPlexCheckCellShape(dmConv,PETSC_TRUE,PETSC_DETERMINE);
334:       DMDestroy(&dmConv);
335:     }

337:     VecDestroy(&preVec);
338:     DMDestroy(&preForest);

340:     preVec    = postVecTransfer;
341:     preForest = postForest;
342:   }

344:   if (transfer_from_base[1]) {
345:     Vec baseVec, baseVecMapped;

347:     /* communicate between base and last adapted forest */
348:     DMGetGlobalVector(base,&baseVec);
349:     DMProjectFunction(base,0.,funcs,ctxs,INSERT_VALUES,baseVec);
350:     PetscObjectSetName((PetscObject)baseVec,"Function Base");
351:     VecViewFromOptions(baseVec,NULL,"-vec_base_view");

353:     DMGetGlobalVector(preForest,&baseVecMapped);
354:     DMForestTransferVecFromBase(preForest,baseVec,baseVecMapped);
355:     VecViewFromOptions(baseVecMapped,NULL,"-vec_map_base_view");

357:     /* compare */
358:     VecAXPY(baseVecMapped,-1.,preVec);
359:     VecViewFromOptions(baseVecMapped,NULL,"-vec_map_diff_view");
360:     VecNorm(baseVecMapped,NORM_2,&diff);

362:     /* output */
363:     if (diff < tol) {
364:       PetscPrintf(comm,"DMForestTransferVecFromBase() passes.\n");
365:     } else {
366:       PetscPrintf(comm,"DMForestTransferVecFromBase() fails with error %g and tolerance %g\n",(double)diff,(double)tol);
367:     }

369:     DMRestoreGlobalVector(base,&baseVec);
370:     DMRestoreGlobalVector(preForest,&baseVecMapped);
371:   }

373:   /* cleanup */
374:   VecDestroy(&preVec);
375:   DMDestroy(&preForest);
376:   DMDestroy(&base);
377:   PetscFinalize();
378:   return ierr;
379: }

381: /*TEST
382:   testset:
383:     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3 -petscspace_type tensor

385:     test:
386:       output_file: output/ex2_2d.out
387:       suffix: p4est_2d
388:       args: -petscspace_degree 2
389:       nsize: 3
390:       requires: p4est !single

392:     test:
393:       output_file: output/ex2_2d.out
394:       suffix: p4est_2d_deg4
395:       args: -petscspace_degree 4
396:       requires: p4est !single

398:     test:
399:       output_file: output/ex2_2d.out
400:       suffix: p4est_2d_deg8
401:       args: -petscspace_degree 8
402:       requires: p4est !single

404:     test:
405:       output_file: output/ex2_steps2.out
406:       suffix: p4est_2d_deg2_steps2
407:       args: -petscspace_degree 2 -coords -adapt_steps 2
408:       nsize: 3
409:       requires: p4est !single

411:     test:
412:       output_file: output/ex2_steps3.out
413:       suffix: p4est_2d_deg3_steps3
414:       args: -petscspace_degree 3 -coords -adapt_steps 3 -petscdualspace_lagrange_node_type equispaced -petscdualspace_lagrange_node_endpoints 1
415:       nsize: 3
416:       requires: p4est !single

418:     test:
419:       output_file: output/ex2_steps3.out
420:       suffix: p4est_2d_deg3_steps3_L2_periodic
421:       args: -petscspace_degree 3 -petscdualspace_lagrange_continuity 0 -coords -adapt_steps 3 -dm_plex_box_bd periodic,periodic -use_bcs 0 -petscdualspace_lagrange_node_type equispaced
422:       nsize: 3
423:       requires: p4est !single

425:     test:
426:       output_file: output/ex2_steps3.out
427:       suffix: p4est_3d_deg2_steps3_L2_periodic
428:       args: -dm_plex_dim 3 -petscspace_degree 2 -petscdualspace_lagrange_continuity 0 -coords -adapt_steps 3 -dm_plex_box_bd periodic,periodic,periodic -use_bcs 0
429:       nsize: 3
430:       requires: p4est !single

432:     test:
433:       output_file: output/ex2_steps2.out
434:       suffix: p4est_3d_deg2_steps2
435:       args: -dm_plex_dim 3 -petscspace_degree 2 -coords -adapt_steps 2
436:       nsize: 3
437:       requires: p4est !single

439:     test:
440:       output_file: output/ex2_steps3.out
441:       suffix: p4est_3d_deg3_steps3
442:       args: -dm_plex_dim 3 -petscspace_degree 3 -coords -adapt_steps 3 -petscdualspace_lagrange_node_type equispaced -petscdualspace_lagrange_node_endpoints 1
443:       nsize: 3
444:       requires: p4est !single

446:     test:
447:       output_file: output/ex2_3d.out
448:       suffix: p4est_3d
449:       args: -dm_plex_dim 3 -petscspace_degree 1
450:       nsize: 3
451:       requires: p4est !single

453:     test:
454:       output_file: output/ex2_3d.out
455:       suffix: p4est_3d_deg3
456:       args: -dm_plex_dim 3 -petscspace_degree 3
457:       nsize: 3
458:       requires: p4est !single

460:     test:
461:       output_file: output/ex2_2d.out
462:       suffix: p4est_2d_deg2_coords
463:       args: -petscspace_degree 2 -coords
464:       nsize: 3
465:       requires: p4est !single

467:     test:
468:       output_file: output/ex2_3d.out
469:       suffix: p4est_3d_deg2_coords
470:       args: -dm_plex_dim 3 -petscspace_degree 2 -coords
471:       nsize: 3
472:       requires: p4est !single

474:     test:
475:       suffix: p4est_3d_nans
476:       args: -dm_plex_dim 3 -dm_forest_partition_overlap 1 -test_convert -petscspace_degree 1
477:       nsize: 2
478:       requires: p4est !single

480:     test:
481:       TODO: not broken, but the 3D case below is broken, so I do not trust this one
482:       output_file: output/ex2_steps2.out
483:       suffix: p4est_2d_tfb_distributed_nc
484:       args: -petscspace_degree 3 -dm_forest_maximum_refinement 2 -dm_p4est_refine_pattern hash -use_bcs 0 -coords -adapt_steps 2 -dm_distribute -petscpartitioner_type shell -petscpartitioner_shell_random
485:       nsize: 3
486:       requires: p4est !single

488:     test:
489:       TODO: broken
490:       output_file: output/ex2_steps2.out
491:       suffix: p4est_3d_tfb_distributed_nc
492:       args: -dm_plex_dim 3 -petscspace_degree 2 -dm_forest_maximum_refinement 2 -dm_p4est_refine_pattern hash -use_bcs 0 -coords -adapt_steps 2 -dm_distribute -petscpartitioner_type shell -petscpartitioner_shell_random
493:       nsize: 3
494:       requires: p4est !single

496:   testset:
497:     args: -petscspace_type tensor -dm_coord_space 0 -dm_plex_transform_type refine_tobox

499:     test:
500:       TODO: broken
501:       output_file: output/ex2_3d.out
502:       suffix: p4est_3d_transfer_fails
503:       args: -petscspace_degree 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -adapt_steps 1 -dm_forest_initial_refinement 1 -use_bcs 0 -dm_refine
504:       requires: p4est !single

506:     test:
507:       TODO: broken
508:       output_file: output/ex2_steps2_notfb.out
509:       suffix: p4est_3d_transfer_fails_2
510:       args: -petscspace_degree 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -adapt_steps 2 -dm_forest_initial_refinement 0 -transfer_from_base 0 -use_bcs 0 -dm_refine
511:       requires: p4est !single

513:     test:
514:       output_file: output/ex2_steps2.out
515:       suffix: p4est_3d_multi_transfer_s2t
516:       args: -petscspace_degree 3 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -adapt_steps 2 -dm_forest_initial_refinement 1 -petscdualspace_lagrange_continuity 0 -use_bcs 0 -dm_refine 1
517:       requires: p4est !single

519:     test:
520:       output_file: output/ex2_steps2.out
521:       suffix: p4est_3d_coords_transfer_s2t
522:       args: -petscspace_degree 3 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -adapt_steps 2 -dm_forest_initial_refinement 1 -petscdualspace_lagrange_continuity 0 -coords -use_bcs 0 -dm_refine 1
523:       requires: p4est !single

525:   testset:
526:     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3

528:     test:
529:       output_file: output/ex2_2d_fv.out
530:       suffix: p4est_2d_fv
531:       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1
532:       nsize: 3
533:       requires: p4est !single

535:     test:
536:       TODO: broken (codimension adjacency)
537:       output_file: output/ex2_2d_fv.out
538:       suffix: p4est_2d_fv_adjcodim
539:       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1 -dm_forest_adjacency_codimension 1
540:       nsize: 2
541:       requires: p4est !single

543:     test:
544:       TODO: broken (dimension adjacency)
545:       output_file: output/ex2_2d_fv.out
546:       suffix: p4est_2d_fv_adjdim
547:       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1 -dm_forest_adjacency_dimension 1
548:       nsize: 2
549:       requires: p4est !single

551:     test:
552:       output_file: output/ex2_2d_fv.out
553:       suffix: p4est_2d_fv_zerocells
554:       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1
555:       nsize: 10
556:       requires: p4est !single

558:     test:
559:       output_file: output/ex2_3d_fv.out
560:       suffix: p4est_3d_fv
561:       args: -dm_plex_dim 3 -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1
562:       nsize: 3
563:       requires: p4est !single

565: TEST*/