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