xref: /petsc/src/dm/impls/forest/tests/ex2.c (revision dfd676b1a855b7f967ece75a22ee7f6626d10f89)
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   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
159   comm = PETSC_COMM_WORLD;
160   PetscOptionsBegin(comm, "", "DMForestTransferVec() Test Options", "DMFOREST");
161   PetscCall(PetscOptionsBool("-linear","Transfer a simple linear function", "ex2.c", linear, &linear, NULL));
162   PetscCall(PetscOptionsBool("-coords","Transfer a simple coordinate function", "ex2.c", coords, &coords, NULL));
163   PetscCall(PetscOptionsBool("-use_fv","Use a finite volume approximation", "ex2.c", useFV, &useFV, NULL));
164   PetscCall(PetscOptionsBool("-test_convert","Test conversion to DMPLEX",NULL,conv,&conv,NULL));
165   PetscCall(PetscOptionsBool("-transfer_from_base","Transfer a vector from base DM to DMForest", "ex2.c", transfer_from_base[0], &transfer_from_base[0], NULL));
166   transfer_from_base[1] = transfer_from_base[0];
167   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));
168   PetscCall(PetscOptionsBool("-use_bcs","Use dirichlet boundary conditions", "ex2.c", use_bcs, &use_bcs, NULL));
169   PetscCall(PetscOptionsBoundedInt("-adapt_steps","Number of adaptivity steps", "ex2.c", adaptSteps, &adaptSteps, NULL,0));
170   PetscOptionsEnd();
171 
172   tol = PetscMax(1.e-10,tol); /* XXX fix for quadruple precision -> why do I need to do this? */
173 
174   /* the base mesh */
175   PetscCall(DMCreate(comm, &base));
176   PetscCall(DMSetType(base, DMPLEX));
177   PetscCall(DMSetFromOptions(base));
178 
179   PetscCall(AddIdentityLabel(base));
180   PetscCall(DMGetDimension(base, &dim));
181 
182   if (linear) {
183     funcs[0] = LinearFunction;
184   }
185   if (coords) {
186     funcs[0] = CoordsFunction;
187     Nf = dim;
188   }
189 
190   bcCtx.func = funcs[0];
191   bcCtx.dim  = dim;
192   bcCtx.Nf   = Nf;
193   bcCtx.ctx  = NULL;
194 
195   if (useFV) {
196     PetscFV      fv;
197     PetscLimiter limiter;
198     DM           baseFV;
199 
200     PetscCall(DMPlexConstructGhostCells(base,NULL,NULL,&baseFV));
201     PetscCall(DMViewFromOptions(baseFV, NULL, "-fv_dm_view"));
202     PetscCall(DMDestroy(&base));
203     base = baseFV;
204     PetscCall(PetscFVCreate(comm, &fv));
205     PetscCall(PetscFVSetSpatialDimension(fv,dim));
206     PetscCall(PetscFVSetType(fv,PETSCFVLEASTSQUARES));
207     PetscCall(PetscFVSetNumComponents(fv,Nf));
208     PetscCall(PetscLimiterCreate(comm,&limiter));
209     PetscCall(PetscLimiterSetType(limiter,PETSCLIMITERNONE));
210     PetscCall(PetscFVSetLimiter(fv,limiter));
211     PetscCall(PetscLimiterDestroy(&limiter));
212     PetscCall(PetscFVSetFromOptions(fv));
213     PetscCall(DMSetField(base,0,NULL,(PetscObject)fv));
214     PetscCall(PetscFVDestroy(&fv));
215   } else {
216     PetscFE fe;
217 
218     PetscCall(PetscFECreateDefault(comm,dim,Nf,PETSC_FALSE,NULL,PETSC_DEFAULT,&fe));
219     PetscCall(DMSetField(base,0,NULL,(PetscObject)fe));
220     PetscCall(PetscFEDestroy(&fe));
221   }
222   PetscCall(DMCreateDS(base));
223 
224   if (use_bcs) {
225     PetscInt ids[] = {1, 2, 3, 4, 5, 6};
226     DMLabel  label;
227 
228     PetscCall(DMGetLabel(base, "marker", &label));
229     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));
230   }
231   PetscCall(DMViewFromOptions(base,NULL,"-dm_base_view"));
232 
233   /* the pre adaptivity forest */
234   PetscCall(DMCreate(comm,&preForest));
235   PetscCall(DMSetType(preForest,(dim == 2) ? DMP4EST : DMP8EST));
236   PetscCall(DMCopyDisc(base,preForest));
237   PetscCall(DMForestSetBaseDM(preForest,base));
238   PetscCall(DMForestSetMinimumRefinement(preForest,0));
239   PetscCall(DMForestSetInitialRefinement(preForest,1));
240   PetscCall(DMSetFromOptions(preForest));
241   PetscCall(DMSetUp(preForest));
242   PetscCall(DMViewFromOptions(preForest,NULL,"-dm_pre_view"));
243 
244   /* the pre adaptivity field */
245   PetscCall(DMCreateGlobalVector(preForest,&preVec));
246   PetscCall(DMProjectFunction(preForest,0.,funcs,ctxs,INSERT_VALUES,preVec));
247   PetscCall(VecViewFromOptions(preVec,NULL,"-vec_pre_view"));
248 
249   /* communicate between base and pre adaptivity forest */
250   if (transfer_from_base[0]) {
251     Vec baseVec, baseVecMapped;
252 
253     PetscCall(DMGetGlobalVector(base,&baseVec));
254     PetscCall(DMProjectFunction(base,0.,funcs,ctxs,INSERT_VALUES,baseVec));
255     PetscCall(PetscObjectSetName((PetscObject)baseVec,"Function Base"));
256     PetscCall(VecViewFromOptions(baseVec,NULL,"-vec_base_view"));
257 
258     PetscCall(DMGetGlobalVector(preForest,&baseVecMapped));
259     PetscCall(DMForestTransferVecFromBase(preForest,baseVec,baseVecMapped));
260     PetscCall(VecViewFromOptions(baseVecMapped,NULL,"-vec_map_base_view"));
261 
262     /* compare */
263     PetscCall(VecAXPY(baseVecMapped,-1.,preVec));
264     PetscCall(VecViewFromOptions(baseVecMapped,NULL,"-vec_map_diff_view"));
265     PetscCall(VecNorm(baseVecMapped,NORM_2,&diff));
266 
267     /* output */
268     if (diff < tol) {
269       PetscCall(PetscPrintf(comm,"DMForestTransferVecFromBase() passes.\n"));
270     } else {
271       PetscCall(PetscPrintf(comm,"DMForestTransferVecFromBase() fails with error %g and tolerance %g\n",(double)diff,(double)tol));
272     }
273 
274     PetscCall(DMRestoreGlobalVector(base,&baseVec));
275     PetscCall(DMRestoreGlobalVector(preForest,&baseVecMapped));
276   }
277 
278   for (step = 0; step < adaptSteps; ++step) {
279 
280     if (!transfer_from_base[1]) {
281       PetscCall(PetscObjectGetReference((PetscObject)preForest,&preCount));
282     }
283 
284     /* adapt */
285     PetscCall(CreateAdaptivityLabel(preForest,&adaptLabel));
286     PetscCall(DMForestTemplate(preForest,comm,&postForest));
287     if (step) PetscCall(DMForestSetAdaptivityLabel(postForest,adaptLabel));
288     PetscCall(DMLabelDestroy(&adaptLabel));
289     PetscCall(DMSetUp(postForest));
290     PetscCall(DMViewFromOptions(postForest,NULL,"-dm_post_view"));
291 
292     /* transfer */
293     PetscCall(DMCreateGlobalVector(postForest,&postVecTransfer));
294     PetscCall(DMForestTransferVec(preForest,preVec,postForest,postVecTransfer,PETSC_TRUE,0.0));
295     PetscCall(VecViewFromOptions(postVecTransfer,NULL,"-vec_post_transfer_view"));
296 
297     /* the exact post adaptivity field */
298     PetscCall(DMCreateGlobalVector(postForest,&postVecExact));
299     PetscCall(DMProjectFunction(postForest,0.,funcs,ctxs,INSERT_VALUES,postVecExact));
300     PetscCall(VecViewFromOptions(postVecExact,NULL,"-vec_post_exact_view"));
301 
302     /* compare */
303     PetscCall(VecAXPY(postVecExact,-1.,postVecTransfer));
304     PetscCall(VecViewFromOptions(postVecExact,NULL,"-vec_diff_view"));
305     PetscCall(VecNorm(postVecExact,NORM_2,&diff));
306 
307     /* output */
308     if (diff < tol) {
309       PetscCall(PetscPrintf(comm,"DMForestTransferVec() passes.\n"));
310     } else {
311       PetscCall(PetscPrintf(comm,"DMForestTransferVec() fails with error %g and tolerance %g\n",(double)diff,(double)tol));
312       PetscCall(IdentifyBadPoints(postForest, postVecExact, tol));
313     }
314     PetscCall(VecDestroy(&postVecExact));
315 
316     /* disconnect preForest from postForest if we don't test the transfer throughout the entire refinement process */
317     if (!transfer_from_base[1]) {
318       PetscCall(DMForestSetAdaptivityForest(postForest,NULL));
319       PetscCall(PetscObjectGetReference((PetscObject)preForest,&postCount));
320       PetscCheck(postCount == preCount,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Adaptation not memory neutral: reference count increase from %" PetscInt_FMT " to %" PetscInt_FMT,preCount,postCount);
321     }
322 
323     if (conv) {
324       DM dmConv;
325 
326       PetscCall(DMConvert(postForest,DMPLEX,&dmConv));
327       PetscCall(DMViewFromOptions(dmConv,NULL,"-dm_conv_view"));
328       PetscCall(DMPlexCheckCellShape(dmConv,PETSC_TRUE,PETSC_DETERMINE));
329       PetscCall(DMDestroy(&dmConv));
330     }
331 
332     PetscCall(VecDestroy(&preVec));
333     PetscCall(DMDestroy(&preForest));
334 
335     preVec    = postVecTransfer;
336     preForest = postForest;
337   }
338 
339   if (transfer_from_base[1]) {
340     Vec baseVec, baseVecMapped;
341 
342     /* communicate between base and last adapted forest */
343     PetscCall(DMGetGlobalVector(base,&baseVec));
344     PetscCall(DMProjectFunction(base,0.,funcs,ctxs,INSERT_VALUES,baseVec));
345     PetscCall(PetscObjectSetName((PetscObject)baseVec,"Function Base"));
346     PetscCall(VecViewFromOptions(baseVec,NULL,"-vec_base_view"));
347 
348     PetscCall(DMGetGlobalVector(preForest,&baseVecMapped));
349     PetscCall(DMForestTransferVecFromBase(preForest,baseVec,baseVecMapped));
350     PetscCall(VecViewFromOptions(baseVecMapped,NULL,"-vec_map_base_view"));
351 
352     /* compare */
353     PetscCall(VecAXPY(baseVecMapped,-1.,preVec));
354     PetscCall(VecViewFromOptions(baseVecMapped,NULL,"-vec_map_diff_view"));
355     PetscCall(VecNorm(baseVecMapped,NORM_2,&diff));
356 
357     /* output */
358     if (diff < tol) {
359       PetscCall(PetscPrintf(comm,"DMForestTransferVecFromBase() passes.\n"));
360     } else {
361       PetscCall(PetscPrintf(comm,"DMForestTransferVecFromBase() fails with error %g and tolerance %g\n",(double)diff,(double)tol));
362     }
363 
364     PetscCall(DMRestoreGlobalVector(base,&baseVec));
365     PetscCall(DMRestoreGlobalVector(preForest,&baseVecMapped));
366   }
367 
368   /* cleanup */
369   PetscCall(VecDestroy(&preVec));
370   PetscCall(DMDestroy(&preForest));
371   PetscCall(DMDestroy(&base));
372   PetscCall(PetscFinalize());
373   return 0;
374 }
375 
376 /*TEST
377   testset:
378     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3 -petscspace_type tensor
379 
380     test:
381       output_file: output/ex2_2d.out
382       suffix: p4est_2d
383       args: -petscspace_degree 2
384       nsize: 3
385       requires: p4est !single
386 
387     test:
388       output_file: output/ex2_2d.out
389       suffix: p4est_2d_deg4
390       args: -petscspace_degree 4
391       requires: p4est !single
392 
393     test:
394       output_file: output/ex2_2d.out
395       suffix: p4est_2d_deg8
396       args: -petscspace_degree 8
397       requires: p4est !single
398 
399     test:
400       output_file: output/ex2_steps2.out
401       suffix: p4est_2d_deg2_steps2
402       args: -petscspace_degree 2 -coords -adapt_steps 2
403       nsize: 3
404       requires: p4est !single
405 
406     test:
407       output_file: output/ex2_steps3.out
408       suffix: p4est_2d_deg3_steps3
409       args: -petscspace_degree 3 -coords -adapt_steps 3 -petscdualspace_lagrange_node_type equispaced -petscdualspace_lagrange_node_endpoints 1
410       nsize: 3
411       requires: p4est !single
412 
413     test:
414       output_file: output/ex2_steps3.out
415       suffix: p4est_2d_deg3_steps3_L2_periodic
416       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
417       nsize: 3
418       requires: p4est !single
419 
420     test:
421       output_file: output/ex2_steps3.out
422       suffix: p4est_3d_deg2_steps3_L2_periodic
423       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
424       nsize: 3
425       requires: p4est !single
426 
427     test:
428       output_file: output/ex2_steps2.out
429       suffix: p4est_3d_deg2_steps2
430       args: -dm_plex_dim 3 -petscspace_degree 2 -coords -adapt_steps 2
431       nsize: 3
432       requires: p4est !single
433 
434     test:
435       output_file: output/ex2_steps3.out
436       suffix: p4est_3d_deg3_steps3
437       args: -dm_plex_dim 3 -petscspace_degree 3 -coords -adapt_steps 3 -petscdualspace_lagrange_node_type equispaced -petscdualspace_lagrange_node_endpoints 1
438       nsize: 3
439       requires: p4est !single
440 
441     test:
442       output_file: output/ex2_3d.out
443       suffix: p4est_3d
444       args: -dm_plex_dim 3 -petscspace_degree 1
445       nsize: 3
446       requires: p4est !single
447 
448     test:
449       output_file: output/ex2_3d.out
450       suffix: p4est_3d_deg3
451       args: -dm_plex_dim 3 -petscspace_degree 3
452       nsize: 3
453       requires: p4est !single
454 
455     test:
456       output_file: output/ex2_2d.out
457       suffix: p4est_2d_deg2_coords
458       args: -petscspace_degree 2 -coords
459       nsize: 3
460       requires: p4est !single
461 
462     test:
463       output_file: output/ex2_3d.out
464       suffix: p4est_3d_deg2_coords
465       args: -dm_plex_dim 3 -petscspace_degree 2 -coords
466       nsize: 3
467       requires: p4est !single
468 
469     test:
470       suffix: p4est_3d_nans
471       args: -dm_plex_dim 3 -dm_forest_partition_overlap 1 -test_convert -petscspace_degree 1
472       nsize: 2
473       requires: p4est !single
474 
475     test:
476       TODO: not broken, but the 3D case below is broken, so I do not trust this one
477       output_file: output/ex2_steps2.out
478       suffix: p4est_2d_tfb_distributed_nc
479       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
480       nsize: 3
481       requires: p4est !single
482 
483     test:
484       TODO: broken
485       output_file: output/ex2_steps2.out
486       suffix: p4est_3d_tfb_distributed_nc
487       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
488       nsize: 3
489       requires: p4est !single
490 
491   testset:
492     args: -petscspace_type tensor -dm_coord_space 0 -dm_plex_transform_type refine_tobox
493 
494     test:
495       TODO: broken
496       output_file: output/ex2_3d.out
497       suffix: p4est_3d_transfer_fails
498       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
499       requires: p4est !single
500 
501     test:
502       TODO: broken
503       output_file: output/ex2_steps2_notfb.out
504       suffix: p4est_3d_transfer_fails_2
505       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
506       requires: p4est !single
507 
508     test:
509       output_file: output/ex2_steps2.out
510       suffix: p4est_3d_multi_transfer_s2t
511       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
512       requires: p4est !single
513 
514     test:
515       output_file: output/ex2_steps2.out
516       suffix: p4est_3d_coords_transfer_s2t
517       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
518       requires: p4est !single
519 
520   testset:
521     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3
522 
523     test:
524       output_file: output/ex2_2d_fv.out
525       suffix: p4est_2d_fv
526       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1
527       nsize: 3
528       requires: p4est !single
529 
530     test:
531       TODO: broken (codimension adjacency)
532       output_file: output/ex2_2d_fv.out
533       suffix: p4est_2d_fv_adjcodim
534       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1 -dm_forest_adjacency_codimension 1
535       nsize: 2
536       requires: p4est !single
537 
538     test:
539       TODO: broken (dimension adjacency)
540       output_file: output/ex2_2d_fv.out
541       suffix: p4est_2d_fv_adjdim
542       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1 -dm_forest_adjacency_dimension 1
543       nsize: 2
544       requires: p4est !single
545 
546     test:
547       output_file: output/ex2_2d_fv.out
548       suffix: p4est_2d_fv_zerocells
549       args: -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1
550       nsize: 10
551       requires: p4est !single
552 
553     test:
554       output_file: output/ex2_3d_fv.out
555       suffix: p4est_3d_fv
556       args: -dm_plex_dim 3 -transfer_from_base 0 -use_fv -linear -dm_forest_partition_overlap 1
557       nsize: 3
558       requires: p4est !single
559 
560 TEST*/
561