Actual source code: ex2.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  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:   DMGetDefaultSection(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 = 2, 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:   char           filename[PETSC_MAX_PATH_LEN];
152:   void           *ctxs[1] = {NULL};
153:   const PetscInt cells[] = {3, 3, 3};
154:   PetscReal      diff, tol = PETSC_SMALL;
155:   DMBoundaryType periodicity[] = {DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE};
156:   PetscBool      linear = PETSC_FALSE;
157:   PetscBool      coords = PETSC_FALSE;
158:   PetscBool      useFV = PETSC_FALSE;
159:   PetscBool      conv = PETSC_FALSE;
160:   PetscBool      distribute_base = PETSC_FALSE;
161:   PetscBool      transfer_from_base[2] = {PETSC_TRUE,PETSC_FALSE};
162:   PetscBool      use_bcs = PETSC_TRUE;
163:   PetscBool      periodic = PETSC_FALSE;
164:   bc_func_ctx    bcCtx;
165:   DMLabel        adaptLabel;
166:   size_t         len;

169:   filename[0] = '\0';
170:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
171:   comm = PETSC_COMM_WORLD;
172:   PetscOptionsBegin(comm, "", "DMForestTransferVec() Test Options", "DMFOREST");
173:   PetscOptionsInt("-dim", "The dimension (2 or 3)", "ex2.c", dim, &dim, NULL);
174:   PetscOptionsBool("-linear","Transfer a simple linear function", "ex2.c", linear, &linear, NULL);
175:   PetscOptionsBool("-coords","Transfer a simple coordinate function", "ex2.c", coords, &coords, NULL);
176:   PetscOptionsBool("-use_fv","Use a finite volume approximation", "ex2.c", useFV, &useFV, NULL);
177:   PetscOptionsBool("-test_convert","Test conversion to DMPLEX",NULL,conv,&conv,NULL);
178:   PetscOptionsString("-filename", "Read the base mesh from file", "ex2.c", filename, filename, PETSC_MAX_PATH_LEN, NULL);
179:   PetscOptionsBool("-distribute_base","Distribute base DM", "ex2.c", distribute_base, &distribute_base, NULL);
180:   PetscOptionsBool("-transfer_from_base","Transfer a vector from base DM to DMForest", "ex2.c", transfer_from_base[0], &transfer_from_base[0], NULL);
181:   transfer_from_base[1] = transfer_from_base[0];
182:   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);
183:   PetscOptionsBool("-use_bcs","Use dirichlet boundary conditions", "ex2.c", use_bcs, &use_bcs, NULL);
184:   PetscOptionsInt("-adapt_steps","Number of adaptivity steps", "ex2.c", adaptSteps, &adaptSteps, NULL);
185:   PetscOptionsBool("-periodic","Use periodic box mesh", "ex2.c", periodic, &periodic, NULL);
186:   PetscOptionsEnd();

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

190:   if (periodic) periodicity[0] = periodicity[1] = periodicity[2] = DM_BOUNDARY_PERIODIC;

192:   if (linear) {
193:     funcs[0] = LinearFunction;
194:   }
195:   if (coords) {
196:     funcs[0] = CoordsFunction;
197:     Nf = dim;
198:   }

200:   bcCtx.func = funcs[0];
201:   bcCtx.dim  = dim;
202:   bcCtx.Nf   = Nf;
203:   bcCtx.ctx  = NULL;

205:   /* the base mesh */
206:   PetscStrlen(filename, &len);
207:   if (!len) {
208:     DMPlexCreateBoxMesh(comm,dim,PETSC_FALSE,cells,NULL,NULL,periodicity,PETSC_TRUE,&base);
209:   } else {
210:     DM tdm = NULL;

212:     DMPlexCreateFromFile(comm,filename,PETSC_TRUE,&base);
213:     DMGetDimension(base,&dim);
214:     DMPlexRefineSimplexToTensor(base,&tdm);
215:     if (tdm) {
216:       DMDestroy(&base);
217:       base = tdm;
218:     }
219:     use_bcs = PETSC_FALSE;
220:   }
221:   AddIdentityLabel(base);
222:   if (distribute_base) {
223:     DM               baseParallel;
224:     PetscPartitioner part;

226:     DMPlexGetPartitioner(base,&part);
227:     PetscPartitionerSetFromOptions(part);
228:     DMPlexDistribute(base,0,NULL,&baseParallel);
229:     if (baseParallel) {
230:       DMDestroy(&base);
231:       base = baseParallel;
232:     }
233:   }
234:   if (useFV) {
235:     PetscFV      fv;
236:     PetscLimiter limiter;
237:     DM           baseFV;

239:     DMPlexConstructGhostCells(base,NULL,NULL,&baseFV);
240:     DMDestroy(&base);
241:     base = baseFV;
242:     PetscFVCreate(comm, &fv);
243:     PetscFVSetSpatialDimension(fv,dim);
244:     PetscFVSetType(fv,PETSCFVLEASTSQUARES);
245:     PetscFVSetNumComponents(fv,Nf);
246:     PetscLimiterCreate(comm,&limiter);
247:     PetscLimiterSetType(limiter,PETSCLIMITERNONE);
248:     PetscFVSetLimiter(fv,limiter);
249:     PetscLimiterDestroy(&limiter);
250:     PetscFVSetFromOptions(fv);
251:     DMSetField(base,0,NULL,(PetscObject)fv);
252:     PetscFVDestroy(&fv);
253:   } else {
254:     PetscFE fe;

256:     PetscFECreateDefault(comm,dim,Nf,PETSC_FALSE,NULL,PETSC_DEFAULT,&fe);
257:     DMSetField(base,0,NULL,(PetscObject)fe);
258:     PetscFEDestroy(&fe);
259:   }
260:   DMCreateDS(base);

262:   if (use_bcs) {
263:     PetscDS  prob;
264:     PetscInt ids[]   = {1, 2, 3, 4, 5, 6};

266:     DMGetDS(base,&prob);
267:     PetscDSAddBoundary(prob,DM_BC_ESSENTIAL, "bc", "marker", 0, 0, NULL, useFV ? (void(*)(void)) bc_func_fv : (void(*)(void)) funcs[0], 2 * dim, ids, useFV ? (void *) &bcCtx : NULL);
268:   }
269:   DMViewFromOptions(base,NULL,"-dm_base_view");

271:   /* the pre adaptivity forest */
272:   DMCreate(comm,&preForest);
273:   DMSetType(preForest,(dim == 2) ? DMP4EST : DMP8EST);
274:   DMCopyDisc(base,preForest);
275:   DMForestSetBaseDM(preForest,base);
276:   DMForestSetMinimumRefinement(preForest,0);
277:   DMForestSetInitialRefinement(preForest,1);
278:   DMSetFromOptions(preForest);
279:   DMSetUp(preForest);
280:   DMViewFromOptions(preForest,NULL,"-dm_pre_view");

282:   /* the pre adaptivity field */
283:   DMCreateGlobalVector(preForest,&preVec);
284:   DMProjectFunction(preForest,0.,funcs,ctxs,INSERT_VALUES,preVec);
285:   VecViewFromOptions(preVec,NULL,"-vec_pre_view");

287:   /* communicate between base and pre adaptivity forest */
288:   if (transfer_from_base[0]) {
289:     Vec baseVec, baseVecMapped;

291:     DMGetGlobalVector(base,&baseVec);
292:     DMProjectFunction(base,0.,funcs,ctxs,INSERT_VALUES,baseVec);
293:     PetscObjectSetName((PetscObject)baseVec,"Function Base");
294:     VecViewFromOptions(baseVec,NULL,"-vec_base_view");

296:     DMGetGlobalVector(preForest,&baseVecMapped);
297:     DMForestTransferVecFromBase(preForest,baseVec,baseVecMapped);
298:     VecViewFromOptions(baseVecMapped,NULL,"-vec_map_base_view");

300:     /* compare */
301:     VecAXPY(baseVecMapped,-1.,preVec);
302:     VecViewFromOptions(baseVecMapped,NULL,"-vec_map_diff_view");
303:     VecNorm(baseVecMapped,NORM_2,&diff);

305:     /* output */
306:     if (diff < tol) {
307:       PetscPrintf(comm,"DMForestTransferVecFromBase() passes.\n");
308:     } else {
309:       PetscPrintf(comm,"DMForestTransferVecFromBase() fails with error %g and tolerance %g\n",(double)diff,(double)tol);
310:     }

312:     DMRestoreGlobalVector(base,&baseVec);
313:     DMRestoreGlobalVector(preForest,&baseVecMapped);
314:   }

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

318:     if (!transfer_from_base[1]) {
319:       PetscObjectGetReference((PetscObject)preForest,&preCount);
320:     }

322:     /* adapt */
323:     CreateAdaptivityLabel(preForest,&adaptLabel);
324:     DMForestTemplate(preForest,comm,&postForest);
325:     if (step) { DMForestSetAdaptivityLabel(postForest,adaptLabel); }
326:     DMLabelDestroy(&adaptLabel);
327:     DMSetUp(postForest);
328:     DMViewFromOptions(postForest,NULL,"-dm_post_view");

330:     /* transfer */
331:     DMCreateGlobalVector(postForest,&postVecTransfer);
332:     DMForestTransferVec(preForest,preVec,postForest,postVecTransfer,PETSC_TRUE,0.0);
333:     VecViewFromOptions(postVecTransfer,NULL,"-vec_post_transfer_view");

335:     /* the exact post adaptivity field */
336:     DMCreateGlobalVector(postForest,&postVecExact);
337:     DMProjectFunction(postForest,0.,funcs,ctxs,INSERT_VALUES,postVecExact);
338:     VecViewFromOptions(postVecExact,NULL,"-vec_post_exact_view");

340:     /* compare */
341:     VecAXPY(postVecExact,-1.,postVecTransfer);
342:     VecViewFromOptions(postVecExact,NULL,"-vec_diff_view");
343:     VecNorm(postVecExact,NORM_2,&diff);

345:     /* output */
346:     if (diff < tol) {
347:       PetscPrintf(comm,"DMForestTransferVec() passes.\n");
348:     } else {
349:       PetscPrintf(comm,"DMForestTransferVec() fails with error %g and tolerance %g\n",(double)diff,(double)tol);
350:       IdentifyBadPoints(postForest, postVecExact, tol);
351:     }
352:     VecDestroy(&postVecExact);

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

361:     if (conv) {
362:       DM dmConv;

364:       DMConvert(postForest,DMPLEX,&dmConv);
365:       DMViewFromOptions(dmConv,NULL,"-dm_conv_view");
366:       DMPlexCheckCellShape(dmConv,PETSC_TRUE);
367:       DMDestroy(&dmConv);
368:     }

370:     VecDestroy(&preVec);
371:     DMDestroy(&preForest);

373:     preVec    = postVecTransfer;
374:     preForest = postForest;
375:   }

377:   if (transfer_from_base[1]) {
378:     Vec baseVec, baseVecMapped;

380:     /* communicate between base and last adapted forest */
381:     DMGetGlobalVector(base,&baseVec);
382:     DMProjectFunction(base,0.,funcs,ctxs,INSERT_VALUES,baseVec);
383:     PetscObjectSetName((PetscObject)baseVec,"Function Base");
384:     VecViewFromOptions(baseVec,NULL,"-vec_base_view");

386:     DMGetGlobalVector(preForest,&baseVecMapped);
387:     DMForestTransferVecFromBase(preForest,baseVec,baseVecMapped);
388:     VecViewFromOptions(baseVecMapped,NULL,"-vec_map_base_view");

390:     /* compare */
391:     VecAXPY(baseVecMapped,-1.,preVec);
392:     VecViewFromOptions(baseVecMapped,NULL,"-vec_map_diff_view");
393:     VecNorm(baseVecMapped,NORM_2,&diff);

395:     /* output */
396:     if (diff < tol) {
397:       PetscPrintf(comm,"DMForestTransferVecFromBase() passes.\n");
398:     } else {
399:       PetscPrintf(comm,"DMForestTransferVecFromBase() fails with error %g and tolerance %g\n",(double)diff,(double)tol);
400:     }

402:     DMRestoreGlobalVector(base,&baseVec);
403:     DMRestoreGlobalVector(preForest,&baseVecMapped);
404:   }

406:   /* cleanup */
407:   VecDestroy(&preVec);
408:   DMDestroy(&preForest);
409:   DMDestroy(&base);
410:   PetscFinalize();
411:   return ierr;
412: }

414: /*TEST

416:      test:
417:        output_file: output/ex2_2d.out
418:        suffix: p4est_2d
419:        args: -petscspace_type tensor -petscspace_degree 2 -dim 2
420:        nsize: 3
421:        requires: p4est !single

423:      test:
424:        output_file: output/ex2_2d.out
425:        suffix: p4est_2d_deg4
426:        args: -petscspace_type tensor -petscspace_degree 4 -dim 2
427:        requires: p4est !single

429:      test:
430:        output_file: output/ex2_2d.out
431:        suffix: p4est_2d_deg8
432:        args: -petscspace_type tensor -petscspace_degree 8 -dim 2
433:        requires: p4est !single

435:      test:
436:        output_file: output/ex2_steps2.out
437:        suffix: p4est_2d_deg2_steps2
438:        args: -petscspace_type tensor -petscspace_degree 2 -dim 2 -coords -adapt_steps 2
439:        nsize: 3
440:        requires: p4est !single

442:      test:
443:        output_file: output/ex2_steps3.out
444:        suffix: p4est_2d_deg3_steps3
445:        args: -petscspace_type tensor -petscspace_degree 3 -dim 2 -coords -adapt_steps 3
446:        nsize: 3
447:        requires: p4est !single

449:      test:
450:        output_file: output/ex2_steps3.out
451:        suffix: p4est_2d_deg3_steps3_L2_periodic
452:        args: -petscspace_type tensor -petscspace_degree 3 -petscdualspace_lagrange_continuity 0 -dim 2 -coords -adapt_steps 3 -periodic -use_bcs 0
453:        nsize: 3
454:        requires: p4est !single

456:      test:
457:        output_file: output/ex2_steps3.out
458:        suffix: p4est_3d_deg2_steps3_L2_periodic
459:        args: -petscspace_type tensor -petscspace_degree 2 -petscdualspace_lagrange_continuity 0 -dim 3 -coords -adapt_steps 3 -periodic -use_bcs 0
460:        nsize: 3
461:        requires: p4est !single

463:      test:
464:        output_file: output/ex2_steps2.out
465:        suffix: p4est_3d_deg2_steps2
466:        args: -petscspace_type tensor -petscspace_degree 2 -dim 3 -coords -adapt_steps 2
467:        nsize: 3
468:        requires: p4est !single

470:      test:
471:        output_file: output/ex2_steps3.out
472:        suffix: p4est_3d_deg3_steps3
473:        args: -petscspace_type tensor -petscspace_degree 3 -dim 3 -coords -adapt_steps 3
474:        nsize: 3
475:        requires: p4est !single

477:      test:
478:        output_file: output/ex2_2d_fv.out
479:        suffix: p4est_2d_fv
480:        args: -transfer_from_base 0 -use_fv -linear -dim 2 -dm_forest_partition_overlap 1
481:        nsize: 3
482:        requires: p4est !single

484:      test:
485:        TODO: broken (codimension adjacency)
486:        output_file: output/ex2_2d_fv.out
487:        suffix: p4est_2d_fv_adjcodim
488:        args: -transfer_from_base 0 -use_fv -linear -dim 2 -dm_forest_partition_overlap 1 -dm_forest_adjacency_codimension 1
489:        nsize: 2
490:        requires: p4est !single

492:      test:
493:        TODO: broken (dimension adjacency)
494:        output_file: output/ex2_2d_fv.out
495:        suffix: p4est_2d_fv_adjdim
496:        args: -transfer_from_base 0 -use_fv -linear -dim 2 -dm_forest_partition_overlap 1 -dm_forest_adjacency_dimension 1
497:        nsize: 2
498:        requires: p4est !single

500:      test:
501:        output_file: output/ex2_2d_fv.out
502:        suffix: p4est_2d_fv_zerocells
503:        args: -transfer_from_base 0 -use_fv -linear -dim 2 -dm_forest_partition_overlap 1
504:        nsize: 10
505:        requires: p4est !single

507:      test:
508:        output_file: output/ex2_3d.out
509:        suffix: p4est_3d
510:        args: -petscspace_type tensor -petscspace_degree 1 -dim 3
511:        nsize: 3
512:        requires: p4est !single

514:      test:
515:        output_file: output/ex2_3d.out
516:        suffix: p4est_3d_deg3
517:        args: -petscspace_type tensor -petscspace_degree 3 -dim 3
518:        nsize: 3
519:        requires: p4est !single

521:      test:
522:        output_file: output/ex2_2d.out
523:        suffix: p4est_2d_deg2_coords
524:        args: -petscspace_type tensor -petscspace_degree 2 -dim 2 -coords
525:        nsize: 3
526:        requires: p4est !single

528:      test:
529:        output_file: output/ex2_3d.out
530:        suffix: p4est_3d_deg2_coords
531:        args: -petscspace_type tensor -petscspace_degree 2 -dim 3 -coords
532:        nsize: 3
533:        requires: p4est !single

535:      test:
536:        output_file: output/ex2_3d_fv.out
537:        suffix: p4est_3d_fv
538:        args: -transfer_from_base 0 -use_fv -linear -dim 3 -dm_forest_partition_overlap 1
539:        nsize: 3
540:        requires: p4est !single

542:      test:
543:        suffix: p4est_3d_nans
544:        args: -dim 3 -dm_forest_partition_overlap 1 -test_convert -petscspace_type tensor -petscspace_degree 1
545:        nsize: 2
546:        requires: p4est !single

548:      test:
549:        TODO: not broken, but the 3D case below is broken, so I do not trust this one
550:        output_file: output/ex2_steps2.out
551:        suffix: p4est_2d_tfb_distributed_nc
552:        args: -petscspace_type tensor -petscspace_degree 3 -dim 2 -dm_forest_maximum_refinement 2 -dm_p4est_refine_pattern hash -use_bcs 0 -coords -adapt_steps 2 -distribute_base -petscpartitioner_type shell -petscpartitioner_shell_random
553:        nsize: 3
554:        requires: p4est !single

556:      test:
557:        TODO: broken
558:        output_file: output/ex2_steps2.out
559:        suffix: p4est_3d_tfb_distributed_nc
560:        args: -petscspace_type tensor -petscspace_degree 2 -dim 3 -dm_forest_maximum_refinement 2 -dm_p4est_refine_pattern hash -use_bcs 0 -coords -adapt_steps 2 -distribute_base -petscpartitioner_type shell -petscpartitioner_shell_random
561:        nsize: 3
562:        requires: p4est !single

564:      test:
565:        TODO: broken
566:        output_file: output/ex2_3d.out
567:        suffix: p4est_3d_transfer_fails
568:        args: -petscspace_type tensor -petscspace_degree 1 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -adapt_steps 1 -dm_forest_initial_refinement 1
569:        requires: p4est !single

571:      test:
572:        TODO: broken
573:        output_file: output/ex2_steps2_notfb.out
574:        suffix: p4est_3d_transfer_fails_2
575:        args: -petscspace_type tensor -petscspace_degree 1 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -adapt_steps 2 -dm_forest_initial_refinement 0 -transfer_from_base 0
576:        requires: p4est !single

578:      test:
579:        output_file: output/ex2_steps2.out
580:        suffix: p4est_3d_multi_transfer_s2t
581:        args: -petscspace_type tensor -petscspace_degree 3 -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
582:        requires: p4est !single

584:      test:
585:        output_file: output/ex2_steps2.out
586:        suffix: p4est_3d_coords_transfer_s2t
587:        args: -petscspace_type tensor -petscspace_degree 3 -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
588:        requires: p4est !single

590: TEST*/