xref: /petsc/src/dm/impls/forest/tests/ex2.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 static char help[] = "Create a mesh, refine and coarsen simultaneously, and transfer a field\n\n";
2 
3 #include <petscds.h>
4 #include <petscdmplex.h>
5 #include <petscdmforest.h>
6 #include <petscoptions.h>
7 
8 static PetscErrorCode AddIdentityLabel(DM dm)
9 {
10   PetscInt       pStart,pEnd,p;
11 
12   PetscFunctionBegin;
13   PetscCall(DMCreateLabel(dm, "identity"));
14   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
15   for (p = pStart; p < pEnd; p++) PetscCall(DMSetLabelValue(dm, "identity", p, p));
16   PetscFunctionReturn(0);
17 }
18 
19 static PetscErrorCode CreateAdaptivityLabel(DM forest,DMLabel *adaptLabel)
20 {
21   DMLabel        identLabel;
22   PetscInt       cStart, cEnd, c;
23 
24   PetscFunctionBegin;
25   PetscCall(DMLabelCreate(PETSC_COMM_SELF,"adapt",adaptLabel));
26   PetscCall(DMLabelSetDefaultValue(*adaptLabel,DM_ADAPT_COARSEN));
27   PetscCall(DMGetLabel(forest,"identity",&identLabel));
28   PetscCall(DMForestGetCellChart(forest,&cStart,&cEnd));
29   for (c = cStart; c < cEnd; c++) {
30     PetscInt basePoint;
31 
32     PetscCall(DMLabelGetValue(identLabel,c,&basePoint));
33     if (!basePoint) PetscCall(DMLabelSetValue(*adaptLabel,c,DM_ADAPT_REFINE));
34   }
35   PetscFunctionReturn(0);
36 }
37 
38 static PetscErrorCode LinearFunction(PetscInt dim,PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
39 {
40   PetscFunctionBeginUser;
41   u[0] = (x[0] * 2.0 + 1.) + (x[1] * 20.0 + 10.) + ((dim == 3) ? (x[2] * 200.0 + 100.) : 0.);
42   PetscFunctionReturn(0);
43 }
44 
45 static PetscErrorCode MultiaffineFunction(PetscInt dim,PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
46 {
47   PetscFunctionBeginUser;
48   u[0] = (x[0] * 1.0 + 2.0) * (x[1] * 3.0 - 4.0) * ((dim == 3) ? (x[2] * 5.0 + 6.0) : 1.);
49   PetscFunctionReturn(0);
50 }
51 
52 static PetscErrorCode CoordsFunction(PetscInt dim,PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
53 {
54   PetscInt f;
55 
56   PetscFunctionBeginUser;
57   for (f=0;f<Nf;f++) u[f] = x[f];
58   PetscFunctionReturn(0);
59 }
60 
61 typedef struct _bc_func_ctx
62 {
63   PetscErrorCode (*func) (PetscInt,PetscReal,const PetscReal [], PetscInt, PetscScalar [], void *);
64   PetscInt dim;
65   PetscInt Nf;
66   void *ctx;
67 }
68 bc_func_ctx;
69 
70 static PetscErrorCode bc_func_fv (PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
71 {
72   bc_func_ctx    *bcCtx;
73 
74   PetscFunctionBegin;
75   bcCtx = (bc_func_ctx *) ctx;
76   PetscCall((bcCtx->func)(bcCtx->dim,time,c,bcCtx->Nf,xG,bcCtx->ctx));
77   PetscFunctionReturn(0);
78 }
79 
80 static PetscErrorCode IdentifyBadPoints (DM dm, Vec vec, PetscReal tol)
81 {
82   DM             dmplex;
83   PetscInt       p, pStart, pEnd, maxDof;
84   Vec            vecLocal;
85   DMLabel        depthLabel;
86   PetscSection   section;
87 
88   PetscFunctionBegin;
89   PetscCall(DMCreateLocalVector(dm, &vecLocal));
90   PetscCall(DMGlobalToLocalBegin(dm, vec, INSERT_VALUES, vecLocal));
91   PetscCall(DMGlobalToLocalEnd(dm, vec, INSERT_VALUES, vecLocal));
92   PetscCall(DMConvert(dm ,DMPLEX, &dmplex));
93   PetscCall(DMPlexGetChart(dmplex, &pStart, &pEnd));
94   PetscCall(DMPlexGetDepthLabel(dmplex, &depthLabel));
95   PetscCall(DMGetLocalSection(dmplex, &section));
96   PetscCall(PetscSectionGetMaxDof(section, &maxDof));
97   for (p = pStart; p < pEnd; p++) {
98     PetscInt     s, c, cSize, parent, childID, numChildren;
99     PetscInt     cl, closureSize, *closure = NULL;
100     PetscScalar *values = NULL;
101     PetscBool    bad = PETSC_FALSE;
102 
103     PetscCall(VecGetValuesSection(vecLocal, section, p, &values));
104     PetscCall(PetscSectionGetDof(section, p, &cSize));
105     for (c = 0; c < cSize; c++) {
106       PetscReal absDiff = PetscAbsScalar(values[c]);
107       if (absDiff > tol) {bad = PETSC_TRUE; break;}
108     }
109     if (!bad) continue;
110     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Bad point %" PetscInt_FMT "\n", p));
111     PetscCall(DMLabelGetValue(depthLabel, p, &s));
112     PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Depth %" PetscInt_FMT "\n", s));
113     PetscCall(DMPlexGetTransitiveClosure(dmplex, p, PETSC_TRUE, &closureSize, &closure));
114     for (cl = 0; cl < closureSize; cl++) {
115       PetscInt cp = closure[2 * cl];
116       PetscCall(DMPlexGetTreeParent(dmplex, cp, &parent, &childID));
117       if (parent != cp) {
118         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Closure point %" PetscInt_FMT " (%" PetscInt_FMT ") child of %" PetscInt_FMT " (ID %" PetscInt_FMT ")\n", cl, cp, parent, childID));
119       }
120       PetscCall(DMPlexGetTreeChildren(dmplex, cp, &numChildren, NULL));
121       if (numChildren) {
122         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Closure point %" PetscInt_FMT " (%" PetscInt_FMT ") is parent\n", cl, cp));
123       }
124     }
125     PetscCall(DMPlexRestoreTransitiveClosure(dmplex, p, PETSC_TRUE, &closureSize, &closure));
126     for (c = 0; c < cSize; c++) {
127       PetscReal absDiff = PetscAbsScalar(values[c]);
128       if (absDiff > tol) {
129         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Bad dof %" PetscInt_FMT "\n", c));
130       }
131     }
132   }
133   PetscCall(DMDestroy(&dmplex));
134   PetscCall(VecDestroy(&vecLocal));
135   PetscFunctionReturn(0);
136 }
137 
138 int main(int argc, char **argv)
139 {
140   MPI_Comm       comm;
141   DM             base, preForest, postForest;
142   PetscInt       dim, Nf = 1;
143   PetscInt       step, adaptSteps = 1;
144   PetscInt       preCount, postCount;
145   Vec            preVec, postVecTransfer, postVecExact;
146   PetscErrorCode (*funcs[1]) (PetscInt,PetscReal,const PetscReal [],PetscInt,PetscScalar [], void *) = {MultiaffineFunction};
147   void           *ctxs[1] = {NULL};
148   PetscReal      diff, tol = PETSC_SMALL;
149   PetscBool      linear = PETSC_FALSE;
150   PetscBool      coords = PETSC_FALSE;
151   PetscBool      useFV = PETSC_FALSE;
152   PetscBool      conv = PETSC_FALSE;
153   PetscBool      transfer_from_base[2] = {PETSC_TRUE,PETSC_FALSE};
154   PetscBool      use_bcs = PETSC_TRUE;
155   bc_func_ctx    bcCtx;
156   DMLabel        adaptLabel;
157 
158   PetscFunctionBeginUser;
159   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
160   comm = PETSC_COMM_WORLD;
161   PetscOptionsBegin(comm, "", "DMForestTransferVec() Test Options", "DMFOREST");
162   PetscCall(PetscOptionsBool("-linear","Transfer a simple linear function", "ex2.c", linear, &linear, NULL));
163   PetscCall(PetscOptionsBool("-coords","Transfer a simple coordinate function", "ex2.c", coords, &coords, NULL));
164   PetscCall(PetscOptionsBool("-use_fv","Use a finite volume approximation", "ex2.c", useFV, &useFV, NULL));
165   PetscCall(PetscOptionsBool("-test_convert","Test conversion to DMPLEX",NULL,conv,&conv,NULL));
166   PetscCall(PetscOptionsBool("-transfer_from_base","Transfer a vector from base DM to DMForest", "ex2.c", transfer_from_base[0], &transfer_from_base[0], NULL));
167   transfer_from_base[1] = transfer_from_base[0];
168   PetscCall(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));
169   PetscCall(PetscOptionsBool("-use_bcs","Use dirichlet boundary conditions", "ex2.c", use_bcs, &use_bcs, NULL));
170   PetscCall(PetscOptionsBoundedInt("-adapt_steps","Number of adaptivity steps", "ex2.c", adaptSteps, &adaptSteps, NULL,0));
171   PetscOptionsEnd();
172 
173   tol = PetscMax(1.e-10,tol); /* XXX fix for quadruple precision -> why do I need to do this? */
174 
175   /* the base mesh */
176   PetscCall(DMCreate(comm, &base));
177   PetscCall(DMSetType(base, DMPLEX));
178   PetscCall(DMSetFromOptions(base));
179 
180   PetscCall(AddIdentityLabel(base));
181   PetscCall(DMGetDimension(base, &dim));
182 
183   if (linear) {
184     funcs[0] = LinearFunction;
185   }
186   if (coords) {
187     funcs[0] = CoordsFunction;
188     Nf = dim;
189   }
190 
191   bcCtx.func = funcs[0];
192   bcCtx.dim  = dim;
193   bcCtx.Nf   = Nf;
194   bcCtx.ctx  = NULL;
195 
196   if (useFV) {
197     PetscFV      fv;
198     PetscLimiter limiter;
199     DM           baseFV;
200 
201     PetscCall(DMPlexConstructGhostCells(base,NULL,NULL,&baseFV));
202     PetscCall(DMViewFromOptions(baseFV, NULL, "-fv_dm_view"));
203     PetscCall(DMDestroy(&base));
204     base = baseFV;
205     PetscCall(PetscFVCreate(PETSC_COMM_SELF, &fv));
206     PetscCall(PetscFVSetSpatialDimension(fv,dim));
207     PetscCall(PetscFVSetType(fv,PETSCFVLEASTSQUARES));
208     PetscCall(PetscFVSetNumComponents(fv,Nf));
209     PetscCall(PetscLimiterCreate(comm,&limiter));
210     PetscCall(PetscLimiterSetType(limiter,PETSCLIMITERNONE));
211     PetscCall(PetscFVSetLimiter(fv,limiter));
212     PetscCall(PetscLimiterDestroy(&limiter));
213     PetscCall(PetscFVSetFromOptions(fv));
214     PetscCall(DMSetField(base,0,NULL,(PetscObject)fv));
215     PetscCall(PetscFVDestroy(&fv));
216   } else {
217     PetscFE fe;
218 
219     PetscCall(PetscFECreateDefault(comm,dim,Nf,PETSC_FALSE,NULL,PETSC_DEFAULT,&fe));
220     PetscCall(DMSetField(base,0,NULL,(PetscObject)fe));
221     PetscCall(PetscFEDestroy(&fe));
222   }
223   PetscCall(DMCreateDS(base));
224 
225   if (use_bcs) {
226     PetscInt ids[] = {1, 2, 3, 4, 5, 6};
227     DMLabel  label;
228 
229     PetscCall(DMGetLabel(base, "marker", &label));
230     PetscCall(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));
231   }
232   PetscCall(DMViewFromOptions(base,NULL,"-dm_base_view"));
233 
234   /* the pre adaptivity forest */
235   PetscCall(DMCreate(comm,&preForest));
236   PetscCall(DMSetType(preForest,(dim == 2) ? DMP4EST : DMP8EST));
237   PetscCall(DMCopyDisc(base,preForest));
238   PetscCall(DMForestSetBaseDM(preForest,base));
239   PetscCall(DMForestSetMinimumRefinement(preForest,0));
240   PetscCall(DMForestSetInitialRefinement(preForest,1));
241   PetscCall(DMSetFromOptions(preForest));
242   PetscCall(DMSetUp(preForest));
243   PetscCall(DMViewFromOptions(preForest,NULL,"-dm_pre_view"));
244 
245   /* the pre adaptivity field */
246   PetscCall(DMCreateGlobalVector(preForest,&preVec));
247   PetscCall(DMProjectFunction(preForest,0.,funcs,ctxs,INSERT_VALUES,preVec));
248   PetscCall(VecViewFromOptions(preVec,NULL,"-vec_pre_view"));
249 
250   /* communicate between base and pre adaptivity forest */
251   if (transfer_from_base[0]) {
252     Vec baseVec, baseVecMapped;
253 
254     PetscCall(DMGetGlobalVector(base,&baseVec));
255     PetscCall(DMProjectFunction(base,0.,funcs,ctxs,INSERT_VALUES,baseVec));
256     PetscCall(PetscObjectSetName((PetscObject)baseVec,"Function Base"));
257     PetscCall(VecViewFromOptions(baseVec,NULL,"-vec_base_view"));
258 
259     PetscCall(DMGetGlobalVector(preForest,&baseVecMapped));
260     PetscCall(DMForestTransferVecFromBase(preForest,baseVec,baseVecMapped));
261     PetscCall(VecViewFromOptions(baseVecMapped,NULL,"-vec_map_base_view"));
262 
263     /* compare */
264     PetscCall(VecAXPY(baseVecMapped,-1.,preVec));
265     PetscCall(VecViewFromOptions(baseVecMapped,NULL,"-vec_map_diff_view"));
266     PetscCall(VecNorm(baseVecMapped,NORM_2,&diff));
267 
268     /* output */
269     if (diff < tol) {
270       PetscCall(PetscPrintf(comm,"DMForestTransferVecFromBase() passes.\n"));
271     } else {
272       PetscCall(PetscPrintf(comm,"DMForestTransferVecFromBase() fails with error %g and tolerance %g\n",(double)diff,(double)tol));
273     }
274 
275     PetscCall(DMRestoreGlobalVector(base,&baseVec));
276     PetscCall(DMRestoreGlobalVector(preForest,&baseVecMapped));
277   }
278 
279   for (step = 0; step < adaptSteps; ++step) {
280 
281     if (!transfer_from_base[1]) {
282       PetscCall(PetscObjectGetReference((PetscObject)preForest,&preCount));
283     }
284 
285     /* adapt */
286     PetscCall(CreateAdaptivityLabel(preForest,&adaptLabel));
287     PetscCall(DMForestTemplate(preForest,comm,&postForest));
288     if (step) PetscCall(DMForestSetAdaptivityLabel(postForest,adaptLabel));
289     PetscCall(DMLabelDestroy(&adaptLabel));
290     PetscCall(DMSetUp(postForest));
291     PetscCall(DMViewFromOptions(postForest,NULL,"-dm_post_view"));
292 
293     /* transfer */
294     PetscCall(DMCreateGlobalVector(postForest,&postVecTransfer));
295     PetscCall(DMForestTransferVec(preForest,preVec,postForest,postVecTransfer,PETSC_TRUE,0.0));
296     PetscCall(VecViewFromOptions(postVecTransfer,NULL,"-vec_post_transfer_view"));
297 
298     /* the exact post adaptivity field */
299     PetscCall(DMCreateGlobalVector(postForest,&postVecExact));
300     PetscCall(DMProjectFunction(postForest,0.,funcs,ctxs,INSERT_VALUES,postVecExact));
301     PetscCall(VecViewFromOptions(postVecExact,NULL,"-vec_post_exact_view"));
302 
303     /* compare */
304     PetscCall(VecAXPY(postVecExact,-1.,postVecTransfer));
305     PetscCall(VecViewFromOptions(postVecExact,NULL,"-vec_diff_view"));
306     PetscCall(VecNorm(postVecExact,NORM_2,&diff));
307 
308     /* output */
309     if (diff < tol) {
310       PetscCall(PetscPrintf(comm,"DMForestTransferVec() passes.\n"));
311     } else {
312       PetscCall(PetscPrintf(comm,"DMForestTransferVec() fails with error %g and tolerance %g\n",(double)diff,(double)tol));
313       PetscCall(IdentifyBadPoints(postForest, postVecExact, tol));
314     }
315     PetscCall(VecDestroy(&postVecExact));
316 
317     /* disconnect preForest from postForest if we don't test the transfer throughout the entire refinement process */
318     if (!transfer_from_base[1]) {
319       PetscCall(DMForestSetAdaptivityForest(postForest,NULL));
320       PetscCall(PetscObjectGetReference((PetscObject)preForest,&postCount));
321       PetscCheck(postCount == preCount,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Adaptation not memory neutral: reference count increase from %" PetscInt_FMT " to %" PetscInt_FMT,preCount,postCount);
322     }
323 
324     if (conv) {
325       DM dmConv;
326 
327       PetscCall(DMConvert(postForest,DMPLEX,&dmConv));
328       PetscCall(DMViewFromOptions(dmConv,NULL,"-dm_conv_view"));
329       PetscCall(DMPlexCheckCellShape(dmConv,PETSC_TRUE,PETSC_DETERMINE));
330       PetscCall(DMDestroy(&dmConv));
331     }
332 
333     PetscCall(VecDestroy(&preVec));
334     PetscCall(DMDestroy(&preForest));
335 
336     preVec    = postVecTransfer;
337     preForest = postForest;
338   }
339 
340   if (transfer_from_base[1]) {
341     Vec baseVec, baseVecMapped;
342 
343     /* communicate between base and last adapted forest */
344     PetscCall(DMGetGlobalVector(base,&baseVec));
345     PetscCall(DMProjectFunction(base,0.,funcs,ctxs,INSERT_VALUES,baseVec));
346     PetscCall(PetscObjectSetName((PetscObject)baseVec,"Function Base"));
347     PetscCall(VecViewFromOptions(baseVec,NULL,"-vec_base_view"));
348 
349     PetscCall(DMGetGlobalVector(preForest,&baseVecMapped));
350     PetscCall(DMForestTransferVecFromBase(preForest,baseVec,baseVecMapped));
351     PetscCall(VecViewFromOptions(baseVecMapped,NULL,"-vec_map_base_view"));
352 
353     /* compare */
354     PetscCall(VecAXPY(baseVecMapped,-1.,preVec));
355     PetscCall(VecViewFromOptions(baseVecMapped,NULL,"-vec_map_diff_view"));
356     PetscCall(VecNorm(baseVecMapped,NORM_2,&diff));
357 
358     /* output */
359     if (diff < tol) {
360       PetscCall(PetscPrintf(comm,"DMForestTransferVecFromBase() passes.\n"));
361     } else {
362       PetscCall(PetscPrintf(comm,"DMForestTransferVecFromBase() fails with error %g and tolerance %g\n",(double)diff,(double)tol));
363     }
364 
365     PetscCall(DMRestoreGlobalVector(base,&baseVec));
366     PetscCall(DMRestoreGlobalVector(preForest,&baseVecMapped));
367   }
368 
369   /* cleanup */
370   PetscCall(VecDestroy(&preVec));
371   PetscCall(DMDestroy(&preForest));
372   PetscCall(DMDestroy(&base));
373   PetscCall(PetscFinalize());
374   return 0;
375 }
376 
377 /*TEST
378   testset:
379     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3 -petscspace_type tensor
380 
381     test:
382       output_file: output/ex2_2d.out
383       suffix: p4est_2d
384       args: -petscspace_degree 2
385       nsize: 3
386       requires: p4est !single
387 
388     test:
389       output_file: output/ex2_2d.out
390       suffix: p4est_2d_deg4
391       args: -petscspace_degree 4
392       requires: p4est !single
393 
394     test:
395       output_file: output/ex2_2d.out
396       suffix: p4est_2d_deg8
397       args: -petscspace_degree 8
398       requires: p4est !single
399 
400     test:
401       output_file: output/ex2_steps2.out
402       suffix: p4est_2d_deg2_steps2
403       args: -petscspace_degree 2 -coords -adapt_steps 2
404       nsize: 3
405       requires: p4est !single
406 
407     test:
408       output_file: output/ex2_steps3.out
409       suffix: p4est_2d_deg3_steps3
410       args: -petscspace_degree 3 -coords -adapt_steps 3 -petscdualspace_lagrange_node_type equispaced -petscdualspace_lagrange_node_endpoints 1
411       nsize: 3
412       requires: p4est !single
413 
414     test:
415       output_file: output/ex2_steps3.out
416       suffix: p4est_2d_deg3_steps3_L2_periodic
417       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
418       nsize: 3
419       requires: p4est !single
420 
421     test:
422       output_file: output/ex2_steps3.out
423       suffix: p4est_3d_deg2_steps3_L2_periodic
424       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
425       nsize: 3
426       requires: p4est !single
427 
428     test:
429       output_file: output/ex2_steps2.out
430       suffix: p4est_3d_deg2_steps2
431       args: -dm_plex_dim 3 -petscspace_degree 2 -coords -adapt_steps 2
432       nsize: 3
433       requires: p4est !single
434 
435     test:
436       output_file: output/ex2_steps3.out
437       suffix: p4est_3d_deg3_steps3
438       args: -dm_plex_dim 3 -petscspace_degree 3 -coords -adapt_steps 3 -petscdualspace_lagrange_node_type equispaced -petscdualspace_lagrange_node_endpoints 1
439       nsize: 3
440       requires: p4est !single
441 
442     test:
443       output_file: output/ex2_3d.out
444       suffix: p4est_3d
445       args: -dm_plex_dim 3 -petscspace_degree 1
446       nsize: 3
447       requires: p4est !single
448 
449     test:
450       output_file: output/ex2_3d.out
451       suffix: p4est_3d_deg3
452       args: -dm_plex_dim 3 -petscspace_degree 3
453       nsize: 3
454       requires: p4est !single
455 
456     test:
457       output_file: output/ex2_2d.out
458       suffix: p4est_2d_deg2_coords
459       args: -petscspace_degree 2 -coords
460       nsize: 3
461       requires: p4est !single
462 
463     test:
464       output_file: output/ex2_3d.out
465       suffix: p4est_3d_deg2_coords
466       args: -dm_plex_dim 3 -petscspace_degree 2 -coords
467       nsize: 3
468       requires: p4est !single
469 
470     test:
471       suffix: p4est_3d_nans
472       args: -dm_plex_dim 3 -dm_forest_partition_overlap 1 -test_convert -petscspace_degree 1
473       nsize: 2
474       requires: p4est !single
475 
476     test:
477       TODO: not broken, but the 3D case below is broken, so I do not trust this one
478       output_file: output/ex2_steps2.out
479       suffix: p4est_2d_tfb_distributed_nc
480       args: -petscspace_degree 3 -dm_forest_maximum_refinement 2 -dm_p4est_refine_pattern hash -use_bcs 0 -coords -adapt_steps 2 -petscpartitioner_type shell -petscpartitioner_shell_random
481       nsize: 3
482       requires: p4est !single
483 
484     test:
485       TODO: broken
486       output_file: output/ex2_steps2.out
487       suffix: p4est_3d_tfb_distributed_nc
488       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 -petscpartitioner_type shell -petscpartitioner_shell_random
489       nsize: 3
490       requires: p4est !single
491 
492   testset:
493     args: -petscspace_type tensor -dm_coord_space 0 -dm_plex_transform_type refine_tobox
494 
495     test:
496       TODO: broken
497       output_file: output/ex2_3d.out
498       suffix: p4est_3d_transfer_fails
499       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
500       requires: p4est !single
501 
502     test:
503       TODO: broken
504       output_file: output/ex2_steps2_notfb.out
505       suffix: p4est_3d_transfer_fails_2
506       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
507       requires: p4est !single
508 
509     test:
510       output_file: output/ex2_steps2.out
511       suffix: p4est_3d_multi_transfer_s2t
512       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
513       requires: p4est !single
514 
515     test:
516       output_file: output/ex2_steps2.out
517       suffix: p4est_3d_coords_transfer_s2t
518       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
519       requires: p4est !single
520 
521   testset:
522     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3
523 
524     test:
525       output_file: output/ex2_2d_fv.out
526       suffix: p4est_2d_fv
527       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1
528       nsize: 3
529       requires: p4est !single
530 
531     test:
532       TODO: broken (codimension adjacency)
533       output_file: output/ex2_2d_fv.out
534       suffix: p4est_2d_fv_adjcodim
535       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1 -dm_forest_adjacency_codimension 1
536       nsize: 2
537       requires: p4est !single
538 
539     test:
540       TODO: broken (dimension adjacency)
541       output_file: output/ex2_2d_fv.out
542       suffix: p4est_2d_fv_adjdim
543       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1 -dm_forest_adjacency_dimension 1
544       nsize: 2
545       requires: p4est !single
546 
547     test:
548       output_file: output/ex2_2d_fv.out
549       suffix: p4est_2d_fv_zerocells
550       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1
551       nsize: 10
552       requires: p4est !single
553 
554     test:
555       output_file: output/ex2_3d_fv.out
556       suffix: p4est_3d_fv
557       args: -dm_plex_dim 3 -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1
558       nsize: 3
559       requires: p4est !single
560 
561 TEST*/
562