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