xref: /petsc/src/ts/tutorials/ex11.c (revision d8a51d2aac522cef49bfbd84fd4cf473dc3fa6ee)
1 static char help[] = "Second Order TVD Finite Volume Example.\n";
2 /*F
3 
4 We use a second order TVD finite volume method to evolve a system of PDEs. Our simple upwinded residual evaluation loops
5 over all mesh faces and uses a Riemann solver to produce the flux given the face geometry and cell values,
6 \begin{equation}
7   f_i = \mathrm{riemann}(\mathrm{phys}, p_\mathrm{centroid}, \hat n, x^L, x^R)
8 \end{equation}
9 and then update the cell values given the cell volume.
10 \begin{eqnarray}
11     f^L_i &-=& \frac{f_i}{vol^L} \\
12     f^R_i &+=& \frac{f_i}{vol^R}
13 \end{eqnarray}
14 
15 As an example, we can consider the shallow water wave equation,
16 \begin{eqnarray}
17      h_t + \nabla\cdot \left( uh                              \right) &=& 0 \\
18   (uh)_t + \nabla\cdot \left( u\otimes uh + \frac{g h^2}{2} I \right) &=& 0
19 \end{eqnarray}
20 where $h$ is wave height, $u$ is wave velocity, and $g$ is the acceleration due to gravity.
21 
22 A representative Riemann solver for the shallow water equations is given in the PhysicsRiemann_SW() function,
23 \begin{eqnarray}
24   f^{L,R}_h    &=& uh^{L,R} \cdot \hat n \\
25   f^{L,R}_{uh} &=& \frac{f^{L,R}_h}{h^{L,R}} uh^{L,R} + g (h^{L,R})^2 \hat n \\
26   c^{L,R}      &=& \sqrt{g h^{L,R}} \\
27   s            &=& \max\left( \left|\frac{uh^L \cdot \hat n}{h^L}\right| + c^L, \left|\frac{uh^R \cdot \hat n}{h^R}\right| + c^R \right) \\
28   f_i          &=& \frac{A_\mathrm{face}}{2} \left( f^L_i + f^R_i + s \left( x^L_i - x^R_i \right) \right)
29 \end{eqnarray}
30 where $c$ is the local gravity wave speed and $f_i$ is a Rusanov flux.
31 
32 The more sophisticated residual evaluation in RHSFunctionLocal_LS() uses a least-squares fit to a quadratic polynomial
33 over a neighborhood of the given element.
34 
35 The mesh is read in from an ExodusII file, usually generated by Cubit.
36 
37 The example also shows how to handle AMR in a time-dependent TS solver.
38 F*/
39 #include <petscdmplex.h>
40 #include <petscdmforest.h>
41 #include <petscds.h>
42 #include <petscts.h>
43 
44 #include "ex11.h"
45 
46 static PetscFunctionList PhysicsList, PhysicsRiemannList_SW, PhysicsRiemannList_Euler;
47 
48 /* 'User' implements a discretization of a continuous model. */
49 typedef struct _n_User *User;
50 typedef PetscErrorCode (*SolutionFunction)(Model, PetscReal, const PetscReal *, PetscScalar *, void *);
51 typedef PetscErrorCode (*SetUpBCFunction)(DM, PetscDS, Physics);
52 typedef PetscErrorCode (*FunctionalFunction)(Model, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *);
53 typedef PetscErrorCode (*SetupFields)(Physics, PetscSection);
54 static PetscErrorCode ModelSolutionSetDefault(Model, SolutionFunction, void *);
55 static PetscErrorCode ModelFunctionalRegister(Model, const char *, PetscInt *, FunctionalFunction, void *);
56 static PetscErrorCode OutputVTK(DM, const char *, PetscViewer *);
57 
58 typedef struct _n_FunctionalLink *FunctionalLink;
59 struct _n_FunctionalLink {
60   char              *name;
61   FunctionalFunction func;
62   void              *ctx;
63   PetscInt           offset;
64   FunctionalLink     next;
65 };
66 
67 struct _n_Model {
68   MPI_Comm         comm; /* Does not do collective communication, but some error conditions can be collective */
69   Physics          physics;
70   FunctionalLink   functionalRegistry;
71   PetscInt         maxComputed;
72   PetscInt         numMonitored;
73   FunctionalLink  *functionalMonitored;
74   PetscInt         numCall;
75   FunctionalLink  *functionalCall;
76   SolutionFunction solution;
77   SetUpBCFunction  setupbc;
78   void            *solutionctx;
79   PetscReal        maxspeed; /* estimate of global maximum speed (for CFL calculation) */
80   PetscReal        bounds[2 * DIM];
81   PetscErrorCode (*errorIndicator)(PetscInt, PetscReal, PetscInt, const PetscScalar[], const PetscScalar[], PetscReal *, void *);
82   void *errorCtx;
83   PetscErrorCode (*setupCEED)(DM, Physics);
84 };
85 
86 struct _n_User {
87   PetscInt  vtkInterval;                        /* For monitor */
88   char      outputBasename[PETSC_MAX_PATH_LEN]; /* Basename for output files */
89   PetscInt  monitorStepOffset;
90   Model     model;
91   PetscBool vtkmon;
92 };
93 
94 #ifdef PETSC_HAVE_LIBCEED
95 // Free a plain data context that was allocated using PETSc; returning libCEED error codes
96 static int FreeContextPetsc(void *data)
97 {
98   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
99   return CEED_ERROR_SUCCESS;
100 }
101 #endif
102 
103 /******************* Advect ********************/
104 typedef enum {
105   ADVECT_SOL_TILTED,
106   ADVECT_SOL_BUMP,
107   ADVECT_SOL_BUMP_CAVITY
108 } AdvectSolType;
109 static const char *const AdvectSolTypes[] = {"TILTED", "BUMP", "BUMP_CAVITY", "AdvectSolType", "ADVECT_SOL_", 0};
110 typedef enum {
111   ADVECT_SOL_BUMP_CONE,
112   ADVECT_SOL_BUMP_COS
113 } AdvectSolBumpType;
114 static const char *const AdvectSolBumpTypes[] = {"CONE", "COS", "AdvectSolBumpType", "ADVECT_SOL_BUMP_", 0};
115 
116 typedef struct {
117   PetscReal wind[DIM];
118 } Physics_Advect_Tilted;
119 typedef struct {
120   PetscReal         center[DIM];
121   PetscReal         radius;
122   AdvectSolBumpType type;
123 } Physics_Advect_Bump;
124 
125 typedef struct {
126   PetscReal     inflowState;
127   AdvectSolType soltype;
128   union
129   {
130     Physics_Advect_Tilted tilted;
131     Physics_Advect_Bump   bump;
132   } sol;
133   struct {
134     PetscInt Solution;
135     PetscInt Error;
136   } functional;
137 } Physics_Advect;
138 
139 static const struct FieldDescription PhysicsFields_Advect[] = {
140   {"U",  1},
141   {NULL, 0}
142 };
143 
144 static PetscErrorCode PhysicsBoundary_Advect_Inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
145 {
146   Physics         phys   = (Physics)ctx;
147   Physics_Advect *advect = (Physics_Advect *)phys->data;
148 
149   PetscFunctionBeginUser;
150   xG[0] = advect->inflowState;
151   PetscFunctionReturn(PETSC_SUCCESS);
152 }
153 
154 static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
155 {
156   PetscFunctionBeginUser;
157   xG[0] = xI[0];
158   PetscFunctionReturn(PETSC_SUCCESS);
159 }
160 
161 static void PhysicsRiemann_Advect(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
162 {
163   Physics_Advect *advect = (Physics_Advect *)phys->data;
164   PetscReal       wind[DIM], wn;
165 
166   switch (advect->soltype) {
167   case ADVECT_SOL_TILTED: {
168     Physics_Advect_Tilted *tilted = &advect->sol.tilted;
169     wind[0]                       = tilted->wind[0];
170     wind[1]                       = tilted->wind[1];
171   } break;
172   case ADVECT_SOL_BUMP:
173     wind[0] = -qp[1];
174     wind[1] = qp[0];
175     break;
176   case ADVECT_SOL_BUMP_CAVITY: {
177     PetscInt  i;
178     PetscReal comp2[3] = {0., 0., 0.}, rad2;
179 
180     rad2 = 0.;
181     for (i = 0; i < dim; i++) {
182       comp2[i] = qp[i] * qp[i];
183       rad2 += comp2[i];
184     }
185 
186     wind[0] = -qp[1];
187     wind[1] = qp[0];
188     if (rad2 > 1.) {
189       PetscInt  maxI     = 0;
190       PetscReal maxComp2 = comp2[0];
191 
192       for (i = 1; i < dim; i++) {
193         if (comp2[i] > maxComp2) {
194           maxI     = i;
195           maxComp2 = comp2[i];
196         }
197       }
198       wind[maxI] = 0.;
199     }
200   } break;
201   default: {
202     PetscInt i;
203     for (i = 0; i < DIM; ++i) wind[i] = 0.0;
204   }
205     /* default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for solution type %s",AdvectSolBumpTypes[advect->soltype]); */
206   }
207   wn      = Dot2Real(wind, n);
208   flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn;
209 }
210 
211 static PetscErrorCode PhysicsSolution_Advect(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
212 {
213   Physics         phys   = (Physics)ctx;
214   Physics_Advect *advect = (Physics_Advect *)phys->data;
215 
216   PetscFunctionBeginUser;
217   switch (advect->soltype) {
218   case ADVECT_SOL_TILTED: {
219     PetscReal              x0[DIM];
220     Physics_Advect_Tilted *tilted = &advect->sol.tilted;
221     Waxpy2Real(-time, tilted->wind, x, x0);
222     if (x0[1] > 0) u[0] = 1. * x[0] + 3. * x[1];
223     else u[0] = advect->inflowState;
224   } break;
225   case ADVECT_SOL_BUMP_CAVITY:
226   case ADVECT_SOL_BUMP: {
227     Physics_Advect_Bump *bump = &advect->sol.bump;
228     PetscReal            x0[DIM], v[DIM], r, cost, sint;
229     cost  = PetscCosReal(time);
230     sint  = PetscSinReal(time);
231     x0[0] = cost * x[0] + sint * x[1];
232     x0[1] = -sint * x[0] + cost * x[1];
233     Waxpy2Real(-1, bump->center, x0, v);
234     r = Norm2Real(v);
235     switch (bump->type) {
236     case ADVECT_SOL_BUMP_CONE:
237       u[0] = PetscMax(1 - r / bump->radius, 0);
238       break;
239     case ADVECT_SOL_BUMP_COS:
240       u[0] = 0.5 + 0.5 * PetscCosReal(PetscMin(r / bump->radius, 1) * PETSC_PI);
241       break;
242     }
243   } break;
244   default:
245     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type");
246   }
247   PetscFunctionReturn(PETSC_SUCCESS);
248 }
249 
250 static PetscErrorCode PhysicsFunctional_Advect(Model mod, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx)
251 {
252   Physics         phys      = (Physics)ctx;
253   Physics_Advect *advect    = (Physics_Advect *)phys->data;
254   PetscScalar     yexact[1] = {0.0};
255 
256   PetscFunctionBeginUser;
257   PetscCall(PhysicsSolution_Advect(mod, time, x, yexact, phys));
258   f[advect->functional.Solution] = PetscRealPart(y[0]);
259   f[advect->functional.Error]    = PetscAbsScalar(y[0] - yexact[0]);
260   PetscFunctionReturn(PETSC_SUCCESS);
261 }
262 
263 static PetscErrorCode SetUpBC_Advect(DM dm, PetscDS prob, Physics phys)
264 {
265   const PetscInt inflowids[] = {100, 200, 300}, outflowids[] = {101};
266   DMLabel        label;
267 
268   PetscFunctionBeginUser;
269   /* Register "canned" boundary conditions and defaults for where to apply. */
270   PetscCall(DMGetLabel(dm, "Face Sets", &label));
271   PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "inflow", label, PETSC_STATIC_ARRAY_LENGTH(inflowids), inflowids, 0, 0, NULL, (void (*)(void))PhysicsBoundary_Advect_Inflow, NULL, phys, NULL));
272   PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "outflow", label, PETSC_STATIC_ARRAY_LENGTH(outflowids), outflowids, 0, 0, NULL, (void (*)(void))PhysicsBoundary_Advect_Outflow, NULL, phys, NULL));
273   PetscFunctionReturn(PETSC_SUCCESS);
274 }
275 
276 static PetscErrorCode PhysicsCreate_Advect(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject)
277 {
278   Physics_Advect *advect;
279 
280   PetscFunctionBeginUser;
281   phys->field_desc = PhysicsFields_Advect;
282   phys->riemann    = (PetscRiemannFunc)PhysicsRiemann_Advect;
283   PetscCall(PetscNew(&advect));
284   phys->data   = advect;
285   mod->setupbc = SetUpBC_Advect;
286 
287   PetscOptionsHeadBegin(PetscOptionsObject, "Advect options");
288   {
289     PetscInt two = 2, dof = 1;
290     advect->soltype = ADVECT_SOL_TILTED;
291     PetscCall(PetscOptionsEnum("-advect_sol_type", "solution type", "", AdvectSolTypes, (PetscEnum)advect->soltype, (PetscEnum *)&advect->soltype, NULL));
292     switch (advect->soltype) {
293     case ADVECT_SOL_TILTED: {
294       Physics_Advect_Tilted *tilted = &advect->sol.tilted;
295       two                           = 2;
296       tilted->wind[0]               = 0.0;
297       tilted->wind[1]               = 1.0;
298       PetscCall(PetscOptionsRealArray("-advect_tilted_wind", "background wind vx,vy", "", tilted->wind, &two, NULL));
299       advect->inflowState = -2.0;
300       PetscCall(PetscOptionsRealArray("-advect_tilted_inflow", "Inflow state", "", &advect->inflowState, &dof, NULL));
301       phys->maxspeed = Norm2Real(tilted->wind);
302     } break;
303     case ADVECT_SOL_BUMP_CAVITY:
304     case ADVECT_SOL_BUMP: {
305       Physics_Advect_Bump *bump = &advect->sol.bump;
306       two                       = 2;
307       bump->center[0]           = 2.;
308       bump->center[1]           = 0.;
309       PetscCall(PetscOptionsRealArray("-advect_bump_center", "location of center of bump x,y", "", bump->center, &two, NULL));
310       bump->radius = 0.9;
311       PetscCall(PetscOptionsReal("-advect_bump_radius", "radius of bump", "", bump->radius, &bump->radius, NULL));
312       bump->type = ADVECT_SOL_BUMP_CONE;
313       PetscCall(PetscOptionsEnum("-advect_bump_type", "type of bump", "", AdvectSolBumpTypes, (PetscEnum)bump->type, (PetscEnum *)&bump->type, NULL));
314       phys->maxspeed = 3.; /* radius of mesh, kludge */
315     } break;
316     }
317   }
318   PetscOptionsHeadEnd();
319   /* Initial/transient solution with default boundary conditions */
320   PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_Advect, phys));
321   /* Register "canned" functionals */
322   PetscCall(ModelFunctionalRegister(mod, "Solution", &advect->functional.Solution, PhysicsFunctional_Advect, phys));
323   PetscCall(ModelFunctionalRegister(mod, "Error", &advect->functional.Error, PhysicsFunctional_Advect, phys));
324   PetscFunctionReturn(PETSC_SUCCESS);
325 }
326 
327 /******************* Shallow Water ********************/
328 static const struct FieldDescription PhysicsFields_SW[] = {
329   {"Height",   1  },
330   {"Momentum", DIM},
331   {NULL,       0  }
332 };
333 
334 static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
335 {
336   PetscFunctionBeginUser;
337   xG[0] = xI[0];
338   xG[1] = -xI[1];
339   xG[2] = -xI[2];
340   PetscFunctionReturn(PETSC_SUCCESS);
341 }
342 
343 static PetscErrorCode PhysicsSolution_SW(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
344 {
345   PetscReal dx[2], r, sigma;
346 
347   PetscFunctionBeginUser;
348   PetscCheck(time == 0.0, mod->comm, PETSC_ERR_SUP, "No solution known for time %g", (double)time);
349   dx[0] = x[0] - 1.5;
350   dx[1] = x[1] - 1.0;
351   r     = Norm2Real(dx);
352   sigma = 0.5;
353   u[0]  = 1 + 2 * PetscExpReal(-PetscSqr(r) / (2 * PetscSqr(sigma)));
354   u[1]  = 0.0;
355   u[2]  = 0.0;
356   PetscFunctionReturn(PETSC_SUCCESS);
357 }
358 
359 static PetscErrorCode PhysicsFunctional_SW(Model mod, PetscReal time, const PetscReal *coord, const PetscScalar *xx, PetscReal *f, void *ctx)
360 {
361   Physics       phys = (Physics)ctx;
362   Physics_SW   *sw   = (Physics_SW *)phys->data;
363   const SWNode *x    = (const SWNode *)xx;
364   PetscReal     u[2];
365   PetscReal     h;
366 
367   PetscFunctionBeginUser;
368   h = x->h;
369   Scale2Real(1. / x->h, x->uh, u);
370   f[sw->functional.Height] = h;
371   f[sw->functional.Speed]  = Norm2Real(u) + PetscSqrtReal(sw->gravity * h);
372   f[sw->functional.Energy] = 0.5 * (Dot2Real(x->uh, u) + sw->gravity * PetscSqr(h));
373   PetscFunctionReturn(PETSC_SUCCESS);
374 }
375 
376 static PetscErrorCode SetUpBC_SW(DM dm, PetscDS prob, Physics phys)
377 {
378   const PetscInt wallids[] = {100, 101, 200, 300};
379   DMLabel        label;
380 
381   PetscFunctionBeginUser;
382   PetscCall(DMGetLabel(dm, "Face Sets", &label));
383   PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, PETSC_STATIC_ARRAY_LENGTH(wallids), wallids, 0, 0, NULL, (void (*)(void))PhysicsBoundary_SW_Wall, NULL, phys, NULL));
384   PetscFunctionReturn(PETSC_SUCCESS);
385 }
386 
387 #ifdef PETSC_HAVE_LIBCEED
388 static PetscErrorCode CreateQFunctionContext_SW(Physics phys, Ceed ceed, CeedQFunctionContext *qfCtx)
389 {
390   Physics_SW *in = (Physics_SW *)phys->data;
391   Physics_SW *sw;
392 
393   PetscFunctionBeginUser;
394   PetscCall(PetscCalloc1(1, &sw));
395 
396   sw->gravity = in->gravity;
397 
398   PetscCallCEED(CeedQFunctionContextCreate(ceed, qfCtx));
399   PetscCallCEED(CeedQFunctionContextSetData(*qfCtx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*sw), sw));
400   PetscCallCEED(CeedQFunctionContextSetDataDestroy(*qfCtx, CEED_MEM_HOST, FreeContextPetsc));
401   PetscCallCEED(CeedQFunctionContextRegisterDouble(*qfCtx, "gravity", offsetof(Physics_SW, gravity), 1, "Accelaration due to gravity"));
402   PetscFunctionReturn(0);
403 }
404 #endif
405 
406 static PetscErrorCode SetupCEED_SW(DM dm, Physics physics)
407 {
408 #ifdef PETSC_HAVE_LIBCEED
409   Ceed                 ceed;
410   CeedQFunctionContext qfCtx;
411 #endif
412 
413   PetscFunctionBegin;
414 #ifdef PETSC_HAVE_LIBCEED
415   PetscCall(DMGetCeed(dm, &ceed));
416   PetscCall(CreateQFunctionContext_SW(physics, ceed, &qfCtx));
417   PetscCall(DMCeedCreateFVM(dm, PETSC_TRUE, PhysicsRiemann_SW_Rusanov_CEED, PhysicsRiemann_SW_Rusanov_CEED_loc, qfCtx));
418   PetscCallCEED(CeedQFunctionContextDestroy(&qfCtx));
419 #endif
420   PetscFunctionReturn(PETSC_SUCCESS);
421 }
422 
423 static PetscErrorCode PhysicsCreate_SW(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject)
424 {
425   Physics_SW *sw;
426   char        sw_riemann[64] = "rusanov";
427 
428   PetscFunctionBeginUser;
429   phys->field_desc = PhysicsFields_SW;
430   PetscCall(PetscNew(&sw));
431   phys->data     = sw;
432   mod->setupbc   = SetUpBC_SW;
433   mod->setupCEED = SetupCEED_SW;
434 
435   PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov", PhysicsRiemann_SW_Rusanov));
436   PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_SW, "hll", PhysicsRiemann_SW_HLL));
437 #ifdef PETSC_HAVE_LIBCEED
438   PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov_ceed", PhysicsRiemann_SW_Rusanov_CEED));
439 #endif
440 
441   PetscOptionsHeadBegin(PetscOptionsObject, "SW options");
442   {
443     void (*PhysicsRiemann_SW)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics);
444     sw->gravity = 1.0;
445     PetscCall(PetscOptionsReal("-sw_gravity", "Gravitational constant", "", sw->gravity, &sw->gravity, NULL));
446     PetscCall(PetscOptionsFList("-sw_riemann", "Riemann solver", "", PhysicsRiemannList_SW, sw_riemann, sw_riemann, sizeof sw_riemann, NULL));
447     PetscCall(PetscFunctionListFind(PhysicsRiemannList_SW, sw_riemann, &PhysicsRiemann_SW));
448     phys->riemann = (PetscRiemannFunc)PhysicsRiemann_SW;
449   }
450   PetscOptionsHeadEnd();
451   phys->maxspeed = PetscSqrtReal(2.0 * sw->gravity); /* Mach 1 for depth of 2 */
452 
453   PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_SW, phys));
454   PetscCall(ModelFunctionalRegister(mod, "Height", &sw->functional.Height, PhysicsFunctional_SW, phys));
455   PetscCall(ModelFunctionalRegister(mod, "Speed", &sw->functional.Speed, PhysicsFunctional_SW, phys));
456   PetscCall(ModelFunctionalRegister(mod, "Energy", &sw->functional.Energy, PhysicsFunctional_SW, phys));
457 
458   PetscFunctionReturn(PETSC_SUCCESS);
459 }
460 
461 /******************* Euler Density Shock (EULER_IV_SHOCK,EULER_SS_SHOCK) ********************/
462 /* An initial-value and self-similar solutions of the compressible Euler equations */
463 /* Ravi Samtaney and D. I. Pullin */
464 /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
465 static const struct FieldDescription PhysicsFields_Euler[] = {
466   {"Density",  1  },
467   {"Momentum", DIM},
468   {"Energy",   1  },
469   {NULL,       0  }
470 };
471 
472 /* initial condition */
473 int                   initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx);
474 static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
475 {
476   PetscInt       i;
477   Physics        phys = (Physics)ctx;
478   Physics_Euler *eu   = (Physics_Euler *)phys->data;
479   EulerNode     *uu   = (EulerNode *)u;
480   PetscReal      p0, gamma, c = 0.;
481   PetscFunctionBeginUser;
482   PetscCheck(time == 0.0, mod->comm, PETSC_ERR_SUP, "No solution known for time %g", (double)time);
483 
484   for (i = 0; i < DIM; i++) uu->ru[i] = 0.0; /* zero out initial velocity */
485   /* set E and rho */
486   gamma = eu->gamma;
487 
488   if (eu->type == EULER_IV_SHOCK || eu->type == EULER_SS_SHOCK) {
489     /******************* Euler Density Shock ********************/
490     /* On initial-value and self-similar solutions of the compressible Euler equations */
491     /* Ravi Samtaney and D. I. Pullin */
492     /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
493     /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity,  */
494     p0 = 1.;
495     if (x[0] < 0.0 + x[1] * eu->itana) {
496       if (x[0] < mod->bounds[0] * 0.5) { /* left of shock (1) */
497         PetscReal amach, rho, press, gas1, p1;
498         amach     = eu->amach;
499         rho       = 1.;
500         press     = p0;
501         p1        = press * (1.0 + 2.0 * gamma / (gamma + 1.0) * (amach * amach - 1.0));
502         gas1      = (gamma - 1.0) / (gamma + 1.0);
503         uu->r     = rho * (p1 / press + gas1) / (gas1 * p1 / press + 1.0);
504         uu->ru[0] = ((uu->r - rho) * PetscSqrtReal(gamma * press / rho) * amach);
505         uu->E     = p1 / (gamma - 1.0) + .5 / uu->r * uu->ru[0] * uu->ru[0];
506       } else {      /* left of discontinuity (0) */
507         uu->r = 1.; /* rho = 1 */
508         uu->E = p0 / (gamma - 1.0);
509       }
510     } else { /* right of discontinuity (2) */
511       uu->r = eu->rhoR;
512       uu->E = p0 / (gamma - 1.0);
513     }
514   } else if (eu->type == EULER_SHOCK_TUBE) {
515     /* For (x<x0) set (rho,u,p)=(8,0,10) and for (x>x0) set (rho,u,p)=(1,0,1). Choose x0 to the midpoint of the domain in the x-direction. */
516     if (x[0] < 0.0) {
517       uu->r = 8.;
518       uu->E = 10. / (gamma - 1.);
519     } else {
520       uu->r = 1.;
521       uu->E = 1. / (gamma - 1.);
522     }
523   } else if (eu->type == EULER_LINEAR_WAVE) {
524     initLinearWave(uu, gamma, x, mod->bounds[1] - mod->bounds[0]);
525   } else SETERRQ(mod->comm, PETSC_ERR_SUP, "Unknown type %d", eu->type);
526 
527   /* set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main; */
528   PetscCall(SpeedOfSound_PG(gamma, uu, &c));
529   c = (uu->ru[0] / uu->r) + c;
530   if (c > phys->maxspeed) phys->maxspeed = c;
531 
532   PetscFunctionReturn(PETSC_SUCCESS);
533 }
534 
535 /* PetscReal* => EulerNode* conversion */
536 static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, PetscScalar *a_xG, void *ctx)
537 {
538   PetscInt         i;
539   const EulerNode *xI   = (const EulerNode *)a_xI;
540   EulerNode       *xG   = (EulerNode *)a_xG;
541   Physics          phys = (Physics)ctx;
542   Physics_Euler   *eu   = (Physics_Euler *)phys->data;
543   PetscFunctionBeginUser;
544   xG->r = xI->r;                                     /* ghost cell density - same */
545   xG->E = xI->E;                                     /* ghost cell energy - same */
546   if (n[1] != 0.) {                                  /* top and bottom */
547     xG->ru[0] = xI->ru[0];                           /* copy tang to wall */
548     xG->ru[1] = -xI->ru[1];                          /* reflect perp to t/b wall */
549   } else {                                           /* sides */
550     for (i = 0; i < DIM; i++) xG->ru[i] = xI->ru[i]; /* copy */
551   }
552   if (eu->type == EULER_LINEAR_WAVE) { /* debug */
553 #if 0
554     PetscPrintf(PETSC_COMM_WORLD,"%s coord=%g,%g\n",PETSC_FUNCTION_NAME,(double)c[0],(double)c[1]);
555 #endif
556   }
557   PetscFunctionReturn(PETSC_SUCCESS);
558 }
559 
560 static PetscErrorCode PhysicsFunctional_Euler(Model mod, PetscReal time, const PetscReal *coord, const PetscScalar *xx, PetscReal *f, void *ctx)
561 {
562   Physics          phys = (Physics)ctx;
563   Physics_Euler   *eu   = (Physics_Euler *)phys->data;
564   const EulerNode *x    = (const EulerNode *)xx;
565   PetscReal        p;
566 
567   PetscFunctionBeginUser;
568   f[eu->monitor.Density]  = x->r;
569   f[eu->monitor.Momentum] = NormDIM(x->ru);
570   f[eu->monitor.Energy]   = x->E;
571   f[eu->monitor.Speed]    = NormDIM(x->ru) / x->r;
572   PetscCall(Pressure_PG(eu->gamma, x, &p));
573   f[eu->monitor.Pressure] = p;
574   PetscFunctionReturn(PETSC_SUCCESS);
575 }
576 
577 static PetscErrorCode SetUpBC_Euler(DM dm, PetscDS prob, Physics phys)
578 {
579   Physics_Euler *eu = (Physics_Euler *)phys->data;
580   DMLabel        label;
581 
582   PetscFunctionBeginUser;
583   PetscCall(DMGetLabel(dm, "Face Sets", &label));
584   if (eu->type == EULER_LINEAR_WAVE) {
585     const PetscInt wallids[] = {100, 101};
586     PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, PETSC_STATIC_ARRAY_LENGTH(wallids), wallids, 0, 0, NULL, (void (*)(void))PhysicsBoundary_Euler_Wall, NULL, phys, NULL));
587   } else {
588     const PetscInt wallids[] = {100, 101, 200, 300};
589     PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, PETSC_STATIC_ARRAY_LENGTH(wallids), wallids, 0, 0, NULL, (void (*)(void))PhysicsBoundary_Euler_Wall, NULL, phys, NULL));
590   }
591   PetscFunctionReturn(PETSC_SUCCESS);
592 }
593 
594 #ifdef PETSC_HAVE_LIBCEED
595 static PetscErrorCode CreateQFunctionContext_Euler(Physics phys, Ceed ceed, CeedQFunctionContext *qfCtx)
596 {
597   Physics_Euler *in = (Physics_Euler *)phys->data;
598   Physics_Euler *eu;
599 
600   PetscFunctionBeginUser;
601   PetscCall(PetscCalloc1(1, &eu));
602 
603   eu->gamma = in->gamma;
604 
605   PetscCallCEED(CeedQFunctionContextCreate(ceed, qfCtx));
606   PetscCallCEED(CeedQFunctionContextSetData(*qfCtx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*eu), eu));
607   PetscCallCEED(CeedQFunctionContextSetDataDestroy(*qfCtx, CEED_MEM_HOST, FreeContextPetsc));
608   PetscCallCEED(CeedQFunctionContextRegisterDouble(*qfCtx, "gamma", offsetof(Physics_Euler, gamma), 1, "Heat capacity ratio"));
609   PetscFunctionReturn(0);
610 }
611 #endif
612 
613 static PetscErrorCode SetupCEED_Euler(DM dm, Physics physics)
614 {
615 #ifdef PETSC_HAVE_LIBCEED
616   Ceed                 ceed;
617   CeedQFunctionContext qfCtx;
618 #endif
619 
620   PetscFunctionBegin;
621 #ifdef PETSC_HAVE_LIBCEED
622   PetscCall(DMGetCeed(dm, &ceed));
623   PetscCall(CreateQFunctionContext_Euler(physics, ceed, &qfCtx));
624   PetscCall(DMCeedCreateFVM(dm, PETSC_TRUE, PhysicsRiemann_Euler_Godunov_CEED, PhysicsRiemann_Euler_Godunov_CEED_loc, qfCtx));
625   PetscCallCEED(CeedQFunctionContextDestroy(&qfCtx));
626 #endif
627   PetscFunctionReturn(PETSC_SUCCESS);
628 }
629 
630 static PetscErrorCode PhysicsCreate_Euler(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject)
631 {
632   Physics_Euler *eu;
633 
634   PetscFunctionBeginUser;
635   phys->field_desc = PhysicsFields_Euler;
636   phys->riemann    = (PetscRiemannFunc)PhysicsRiemann_Euler_Godunov;
637   PetscCall(PetscNew(&eu));
638   phys->data     = eu;
639   mod->setupbc   = SetUpBC_Euler;
640   mod->setupCEED = SetupCEED_Euler;
641 
642   PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_Euler, "godunov", PhysicsRiemann_Euler_Godunov));
643 #ifdef PETSC_HAVE_LIBCEED
644   PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_Euler, "godunov_ceed", PhysicsRiemann_Euler_Godunov_CEED));
645 #endif
646 
647   PetscOptionsHeadBegin(PetscOptionsObject, "Euler options");
648   {
649     void (*PhysicsRiemann_Euler)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics);
650     PetscReal alpha;
651     char      type[64]       = "linear_wave";
652     char      eu_riemann[64] = "godunov";
653     PetscBool is;
654     eu->gamma = 1.4;
655     eu->amach = 2.02;
656     eu->rhoR  = 3.0;
657     eu->itana = 0.57735026918963; /* angle of Euler self similar (SS) shock */
658     PetscCall(PetscOptionsFList("-eu_riemann", "Riemann solver", "", PhysicsRiemannList_Euler, eu_riemann, eu_riemann, sizeof eu_riemann, NULL));
659     PetscCall(PetscFunctionListFind(PhysicsRiemannList_Euler, eu_riemann, &PhysicsRiemann_Euler));
660     phys->riemann = (PetscRiemannFunc)PhysicsRiemann_Euler;
661     PetscCall(PetscOptionsReal("-eu_gamma", "Heat capacity ratio", "", eu->gamma, &eu->gamma, NULL));
662     PetscCall(PetscOptionsReal("-eu_amach", "Shock speed (Mach)", "", eu->amach, &eu->amach, NULL));
663     PetscCall(PetscOptionsReal("-eu_rho2", "Density right of discontinuity", "", eu->rhoR, &eu->rhoR, NULL));
664     alpha = 60.;
665     PetscCall(PetscOptionsRangeReal("-eu_alpha", "Angle of discontinuity", "", alpha, &alpha, NULL, 0.0, 90.0));
666     eu->itana = 1. / PetscTanReal(alpha * PETSC_PI / 180.0);
667     PetscCall(PetscOptionsString("-eu_type", "Type of Euler test", "", type, type, sizeof(type), NULL));
668     PetscCall(PetscStrcmp(type, "linear_wave", &is));
669     if (is) {
670       /* Remember this should be periodic */
671       eu->type = EULER_LINEAR_WAVE;
672       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "linear_wave"));
673     } else {
674       PetscCheck(DIM == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "DIM must be 2 unless linear wave test %s", type);
675       PetscCall(PetscStrcmp(type, "iv_shock", &is));
676       if (is) {
677         eu->type = EULER_IV_SHOCK;
678         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "iv_shock"));
679       } else {
680         PetscCall(PetscStrcmp(type, "ss_shock", &is));
681         if (is) {
682           eu->type = EULER_SS_SHOCK;
683           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "ss_shock"));
684         } else {
685           PetscCall(PetscStrcmp(type, "shock_tube", &is));
686           if (is) eu->type = EULER_SHOCK_TUBE;
687           else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unknown Euler type %s", type);
688           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "shock_tube"));
689         }
690       }
691     }
692   }
693   PetscOptionsHeadEnd();
694   phys->maxspeed = 0.; /* will get set in solution */
695   PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_Euler, phys));
696   PetscCall(ModelFunctionalRegister(mod, "Speed", &eu->monitor.Speed, PhysicsFunctional_Euler, phys));
697   PetscCall(ModelFunctionalRegister(mod, "Energy", &eu->monitor.Energy, PhysicsFunctional_Euler, phys));
698   PetscCall(ModelFunctionalRegister(mod, "Density", &eu->monitor.Density, PhysicsFunctional_Euler, phys));
699   PetscCall(ModelFunctionalRegister(mod, "Momentum", &eu->monitor.Momentum, PhysicsFunctional_Euler, phys));
700   PetscCall(ModelFunctionalRegister(mod, "Pressure", &eu->monitor.Pressure, PhysicsFunctional_Euler, phys));
701 
702   PetscFunctionReturn(PETSC_SUCCESS);
703 }
704 
705 static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx)
706 {
707   PetscReal err = 0.;
708   PetscInt  i, j;
709 
710   PetscFunctionBeginUser;
711   for (i = 0; i < numComps; i++) {
712     for (j = 0; j < dim; j++) err += PetscSqr(PetscRealPart(grad[i * dim + j]));
713   }
714   *error = volume * err;
715   PetscFunctionReturn(PETSC_SUCCESS);
716 }
717 
718 PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
719 {
720   PetscSF      sfPoint;
721   PetscSection coordSection;
722   Vec          coordinates;
723   PetscSection sectionCell;
724   PetscScalar *part;
725   PetscInt     cStart, cEnd, c;
726   PetscMPIInt  rank;
727 
728   PetscFunctionBeginUser;
729   PetscCall(DMGetCoordinateSection(dm, &coordSection));
730   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
731   PetscCall(DMClone(dm, dmCell));
732   PetscCall(DMGetPointSF(dm, &sfPoint));
733   PetscCall(DMSetPointSF(*dmCell, sfPoint));
734   PetscCall(DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection));
735   PetscCall(DMSetCoordinatesLocal(*dmCell, coordinates));
736   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
737   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionCell));
738   PetscCall(DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd));
739   PetscCall(PetscSectionSetChart(sectionCell, cStart, cEnd));
740   for (c = cStart; c < cEnd; ++c) PetscCall(PetscSectionSetDof(sectionCell, c, 1));
741   PetscCall(PetscSectionSetUp(sectionCell));
742   PetscCall(DMSetLocalSection(*dmCell, sectionCell));
743   PetscCall(PetscSectionDestroy(&sectionCell));
744   PetscCall(DMCreateLocalVector(*dmCell, partition));
745   PetscCall(PetscObjectSetName((PetscObject)*partition, "partition"));
746   PetscCall(VecGetArray(*partition, &part));
747   for (c = cStart; c < cEnd; ++c) {
748     PetscScalar *p;
749 
750     PetscCall(DMPlexPointLocalRef(*dmCell, c, part, &p));
751     p[0] = rank;
752   }
753   PetscCall(VecRestoreArray(*partition, &part));
754   PetscFunctionReturn(PETSC_SUCCESS);
755 }
756 
757 PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
758 {
759   DM                 plex, dmMass, dmFace, dmCell, dmCoord;
760   PetscSection       coordSection;
761   Vec                coordinates, facegeom, cellgeom;
762   PetscSection       sectionMass;
763   PetscScalar       *m;
764   const PetscScalar *fgeom, *cgeom, *coords;
765   PetscInt           vStart, vEnd, v;
766 
767   PetscFunctionBeginUser;
768   PetscCall(DMConvert(dm, DMPLEX, &plex));
769   PetscCall(DMGetCoordinateSection(dm, &coordSection));
770   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
771   PetscCall(DMClone(dm, &dmMass));
772   PetscCall(DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection));
773   PetscCall(DMSetCoordinatesLocal(dmMass, coordinates));
774   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionMass));
775   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
776   PetscCall(PetscSectionSetChart(sectionMass, vStart, vEnd));
777   for (v = vStart; v < vEnd; ++v) {
778     PetscInt numFaces;
779 
780     PetscCall(DMPlexGetSupportSize(dmMass, v, &numFaces));
781     PetscCall(PetscSectionSetDof(sectionMass, v, numFaces * numFaces));
782   }
783   PetscCall(PetscSectionSetUp(sectionMass));
784   PetscCall(DMSetLocalSection(dmMass, sectionMass));
785   PetscCall(PetscSectionDestroy(&sectionMass));
786   PetscCall(DMGetLocalVector(dmMass, massMatrix));
787   PetscCall(VecGetArray(*massMatrix, &m));
788   PetscCall(DMPlexGetGeometryFVM(plex, &facegeom, &cellgeom, NULL));
789   PetscCall(VecGetDM(facegeom, &dmFace));
790   PetscCall(VecGetArrayRead(facegeom, &fgeom));
791   PetscCall(VecGetDM(cellgeom, &dmCell));
792   PetscCall(VecGetArrayRead(cellgeom, &cgeom));
793   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
794   PetscCall(VecGetArrayRead(coordinates, &coords));
795   for (v = vStart; v < vEnd; ++v) {
796     const PetscInt  *faces;
797     PetscFVFaceGeom *fgA, *fgB, *cg;
798     PetscScalar     *vertex;
799     PetscInt         numFaces, sides[2], f, g;
800 
801     PetscCall(DMPlexPointLocalRead(dmCoord, v, coords, &vertex));
802     PetscCall(DMPlexGetSupportSize(dmMass, v, &numFaces));
803     PetscCall(DMPlexGetSupport(dmMass, v, &faces));
804     for (f = 0; f < numFaces; ++f) {
805       sides[0] = faces[f];
806       PetscCall(DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA));
807       for (g = 0; g < numFaces; ++g) {
808         const PetscInt *cells = NULL;
809         PetscReal       area  = 0.0;
810         PetscInt        numCells;
811 
812         sides[1] = faces[g];
813         PetscCall(DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB));
814         PetscCall(DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells));
815         PetscCheck(numCells == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces");
816         PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg));
817         area += PetscAbsScalar((vertex[0] - cg->centroid[0]) * (fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1]) * (fgA->centroid[0] - cg->centroid[0]));
818         area += PetscAbsScalar((vertex[0] - cg->centroid[0]) * (fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1]) * (fgB->centroid[0] - cg->centroid[0]));
819         m[f * numFaces + g] = Dot2Real(fgA->normal, fgB->normal) * area * 0.5;
820         PetscCall(DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells));
821       }
822     }
823   }
824   PetscCall(VecRestoreArrayRead(facegeom, &fgeom));
825   PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
826   PetscCall(VecRestoreArrayRead(coordinates, &coords));
827   PetscCall(VecRestoreArray(*massMatrix, &m));
828   PetscCall(DMDestroy(&dmMass));
829   PetscCall(DMDestroy(&plex));
830   PetscFunctionReturn(PETSC_SUCCESS);
831 }
832 
833 /* Behavior will be different for multi-physics or when using non-default boundary conditions */
834 static PetscErrorCode ModelSolutionSetDefault(Model mod, SolutionFunction func, void *ctx)
835 {
836   PetscFunctionBeginUser;
837   mod->solution    = func;
838   mod->solutionctx = ctx;
839   PetscFunctionReturn(PETSC_SUCCESS);
840 }
841 
842 static PetscErrorCode ModelFunctionalRegister(Model mod, const char *name, PetscInt *offset, FunctionalFunction func, void *ctx)
843 {
844   FunctionalLink link, *ptr;
845   PetscInt       lastoffset = -1;
846 
847   PetscFunctionBeginUser;
848   for (ptr = &mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
849   PetscCall(PetscNew(&link));
850   PetscCall(PetscStrallocpy(name, &link->name));
851   link->offset = lastoffset + 1;
852   link->func   = func;
853   link->ctx    = ctx;
854   link->next   = NULL;
855   *ptr         = link;
856   *offset      = link->offset;
857   PetscFunctionReturn(PETSC_SUCCESS);
858 }
859 
860 static PetscErrorCode ModelFunctionalSetFromOptions(Model mod, PetscOptionItems *PetscOptionsObject)
861 {
862   PetscInt       i, j;
863   FunctionalLink link;
864   char          *names[256];
865 
866   PetscFunctionBeginUser;
867   mod->numMonitored = PETSC_STATIC_ARRAY_LENGTH(names);
868   PetscCall(PetscOptionsStringArray("-monitor", "list of functionals to monitor", "", names, &mod->numMonitored, NULL));
869   /* Create list of functionals that will be computed somehow */
870   PetscCall(PetscMalloc1(mod->numMonitored, &mod->functionalMonitored));
871   /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
872   PetscCall(PetscMalloc1(mod->numMonitored, &mod->functionalCall));
873   mod->numCall = 0;
874   for (i = 0; i < mod->numMonitored; i++) {
875     for (link = mod->functionalRegistry; link; link = link->next) {
876       PetscBool match;
877       PetscCall(PetscStrcasecmp(names[i], link->name, &match));
878       if (match) break;
879     }
880     PetscCheck(link, mod->comm, PETSC_ERR_USER, "No known functional '%s'", names[i]);
881     mod->functionalMonitored[i] = link;
882     for (j = 0; j < i; j++) {
883       if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
884     }
885     mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
886   next_name:
887     PetscCall(PetscFree(names[i]));
888   }
889 
890   /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
891   mod->maxComputed = -1;
892   for (link = mod->functionalRegistry; link; link = link->next) {
893     for (i = 0; i < mod->numCall; i++) {
894       FunctionalLink call = mod->functionalCall[i];
895       if (link->func == call->func && link->ctx == call->ctx) mod->maxComputed = PetscMax(mod->maxComputed, link->offset);
896     }
897   }
898   PetscFunctionReturn(PETSC_SUCCESS);
899 }
900 
901 static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
902 {
903   FunctionalLink l, next;
904 
905   PetscFunctionBeginUser;
906   if (!link) PetscFunctionReturn(PETSC_SUCCESS);
907   l     = *link;
908   *link = NULL;
909   for (; l; l = next) {
910     next = l->next;
911     PetscCall(PetscFree(l->name));
912     PetscCall(PetscFree(l));
913   }
914   PetscFunctionReturn(PETSC_SUCCESS);
915 }
916 
917 /* put the solution callback into a functional callback */
918 static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx)
919 {
920   Model mod;
921   PetscFunctionBeginUser;
922   mod = (Model)modctx;
923   PetscCall((*mod->solution)(mod, time, x, u, mod->solutionctx));
924   PetscFunctionReturn(PETSC_SUCCESS);
925 }
926 
927 PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
928 {
929   PetscErrorCode (*func[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
930   void *ctx[1];
931   Model mod = user->model;
932 
933   PetscFunctionBeginUser;
934   func[0] = SolutionFunctional;
935   ctx[0]  = (void *)mod;
936   PetscCall(DMProjectFunction(dm, 0.0, func, ctx, INSERT_ALL_VALUES, X));
937   PetscFunctionReturn(PETSC_SUCCESS);
938 }
939 
940 static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer)
941 {
942   PetscFunctionBeginUser;
943   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer));
944   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERVTK));
945   PetscCall(PetscViewerFileSetName(*viewer, filename));
946   PetscFunctionReturn(PETSC_SUCCESS);
947 }
948 
949 static PetscErrorCode MonitorVTK(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx)
950 {
951   User        user = (User)ctx;
952   DM          dm, plex;
953   PetscViewer viewer;
954   char        filename[PETSC_MAX_PATH_LEN], *ftable = NULL;
955   PetscReal   xnorm;
956 
957   PetscFunctionBeginUser;
958   PetscCall(PetscObjectSetName((PetscObject)X, "u"));
959   PetscCall(VecGetDM(X, &dm));
960   PetscCall(VecNorm(X, NORM_INFINITY, &xnorm));
961 
962   if (stepnum >= 0) stepnum += user->monitorStepOffset;
963   if (stepnum >= 0) { /* No summary for final time */
964     Model              mod = user->model;
965     Vec                cellgeom;
966     PetscInt           c, cStart, cEnd, fcount, i;
967     size_t             ftableused, ftablealloc;
968     const PetscScalar *cgeom, *x;
969     DM                 dmCell;
970     DMLabel            vtkLabel;
971     PetscReal         *fmin, *fmax, *fintegral, *ftmp;
972 
973     PetscCall(DMConvert(dm, DMPLEX, &plex));
974     PetscCall(DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL));
975     fcount = mod->maxComputed + 1;
976     PetscCall(PetscMalloc4(fcount, &fmin, fcount, &fmax, fcount, &fintegral, fcount, &ftmp));
977     for (i = 0; i < fcount; i++) {
978       fmin[i]      = PETSC_MAX_REAL;
979       fmax[i]      = PETSC_MIN_REAL;
980       fintegral[i] = 0;
981     }
982     PetscCall(VecGetDM(cellgeom, &dmCell));
983     PetscCall(DMPlexGetSimplexOrBoxCells(dmCell, 0, &cStart, &cEnd));
984     PetscCall(VecGetArrayRead(cellgeom, &cgeom));
985     PetscCall(VecGetArrayRead(X, &x));
986     PetscCall(DMGetLabel(dm, "vtk", &vtkLabel));
987     for (c = cStart; c < cEnd; ++c) {
988       PetscFVCellGeom   *cg;
989       const PetscScalar *cx     = NULL;
990       PetscInt           vtkVal = 0;
991 
992       /* not that these two routines as currently implemented work for any dm with a
993        * localSection/globalSection */
994       PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg));
995       PetscCall(DMPlexPointGlobalRead(dm, c, x, &cx));
996       if (vtkLabel) PetscCall(DMLabelGetValue(vtkLabel, c, &vtkVal));
997       if (!vtkVal || !cx) continue; /* ghost, or not a global cell */
998       for (i = 0; i < mod->numCall; i++) {
999         FunctionalLink flink = mod->functionalCall[i];
1000         PetscCall((*flink->func)(mod, time, cg->centroid, cx, ftmp, flink->ctx));
1001       }
1002       for (i = 0; i < fcount; i++) {
1003         fmin[i] = PetscMin(fmin[i], ftmp[i]);
1004         fmax[i] = PetscMax(fmax[i], ftmp[i]);
1005         fintegral[i] += cg->volume * ftmp[i];
1006       }
1007     }
1008     PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
1009     PetscCall(VecRestoreArrayRead(X, &x));
1010     PetscCall(DMDestroy(&plex));
1011     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts)));
1012     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts)));
1013     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fintegral, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts)));
1014 
1015     ftablealloc = fcount * 100;
1016     ftableused  = 0;
1017     PetscCall(PetscMalloc1(ftablealloc, &ftable));
1018     for (i = 0; i < mod->numMonitored; i++) {
1019       size_t         countused;
1020       char           buffer[256], *p;
1021       FunctionalLink flink = mod->functionalMonitored[i];
1022       PetscInt       id    = flink->offset;
1023       if (i % 3) {
1024         PetscCall(PetscArraycpy(buffer, "  ", 2));
1025         p = buffer + 2;
1026       } else if (i) {
1027         char newline[] = "\n";
1028         PetscCall(PetscMemcpy(buffer, newline, sizeof(newline) - 1));
1029         p = buffer + sizeof(newline) - 1;
1030       } else {
1031         p = buffer;
1032       }
1033       PetscCall(PetscSNPrintfCount(p, sizeof buffer - (p - buffer), "%12s [%10.7g,%10.7g] int %10.7g", &countused, flink->name, (double)fmin[id], (double)fmax[id], (double)fintegral[id]));
1034       countused--;
1035       countused += p - buffer;
1036       if (countused > ftablealloc - ftableused - 1) { /* reallocate */
1037         char *ftablenew;
1038         ftablealloc = 2 * ftablealloc + countused;
1039         PetscCall(PetscMalloc(ftablealloc, &ftablenew));
1040         PetscCall(PetscArraycpy(ftablenew, ftable, ftableused));
1041         PetscCall(PetscFree(ftable));
1042         ftable = ftablenew;
1043       }
1044       PetscCall(PetscArraycpy(ftable + ftableused, buffer, countused));
1045       ftableused += countused;
1046       ftable[ftableused] = 0;
1047     }
1048     PetscCall(PetscFree4(fmin, fmax, fintegral, ftmp));
1049 
1050     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "% 3" PetscInt_FMT "  time %8.4g  |x| %8.4g  %s\n", stepnum, (double)time, (double)xnorm, ftable ? ftable : ""));
1051     PetscCall(PetscFree(ftable));
1052   }
1053   if (user->vtkInterval < 1) PetscFunctionReturn(PETSC_SUCCESS);
1054   if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1055     if (stepnum == -1) { /* Final time is not multiple of normal time interval, write it anyway */
1056       PetscCall(TSGetStepNumber(ts, &stepnum));
1057     }
1058     PetscCall(PetscSNPrintf(filename, sizeof filename, "%s-%03" PetscInt_FMT ".vtu", user->outputBasename, stepnum));
1059     PetscCall(OutputVTK(dm, filename, &viewer));
1060     PetscCall(VecView(X, viewer));
1061     PetscCall(PetscViewerDestroy(&viewer));
1062   }
1063   PetscFunctionReturn(PETSC_SUCCESS);
1064 }
1065 
1066 static PetscErrorCode initializeTS(DM dm, User user, TS *ts)
1067 {
1068   PetscBool useCeed;
1069 
1070   PetscFunctionBeginUser;
1071   PetscCall(TSCreate(PetscObjectComm((PetscObject)dm), ts));
1072   PetscCall(TSSetType(*ts, TSSSP));
1073   PetscCall(TSSetDM(*ts, dm));
1074   if (user->vtkmon) PetscCall(TSMonitorSet(*ts, MonitorVTK, user, NULL));
1075   PetscCall(DMPlexGetUseCeed(dm, &useCeed));
1076   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user));
1077   if (useCeed) PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVMCEED, user));
1078   else PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user));
1079   PetscCall(TSSetMaxTime(*ts, 2.0));
1080   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
1081   PetscFunctionReturn(PETSC_SUCCESS);
1082 }
1083 
1084 typedef struct {
1085   PetscFV      fvm;
1086   VecTagger    refineTag;
1087   VecTagger    coarsenTag;
1088   DM           adaptedDM;
1089   User         user;
1090   PetscReal    cfl;
1091   PetscLimiter limiter;
1092   PetscLimiter noneLimiter;
1093 } TransferCtx;
1094 
1095 static PetscErrorCode adaptToleranceFVMSetUp(TS ts, PetscInt nstep, PetscReal time, Vec sol, PetscBool *resize, void *ctx)
1096 {
1097   TransferCtx       *tctx       = (TransferCtx *)ctx;
1098   PetscFV            fvm        = tctx->fvm;
1099   VecTagger          refineTag  = tctx->refineTag;
1100   VecTagger          coarsenTag = tctx->coarsenTag;
1101   User               user       = tctx->user;
1102   DM                 dm, gradDM, plex, cellDM, adaptedDM = NULL;
1103   Vec                cellGeom, faceGeom;
1104   PetscBool          isForest, computeGradient;
1105   Vec                grad, locGrad, locX, errVec;
1106   PetscInt           cStart, cEnd, c, dim, nRefine, nCoarsen;
1107   PetscReal          minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2];
1108   PetscScalar       *errArray;
1109   const PetscScalar *pointVals;
1110   const PetscScalar *pointGrads;
1111   const PetscScalar *pointGeom;
1112   DMLabel            adaptLabel = NULL;
1113   IS                 refineIS, coarsenIS;
1114 
1115   PetscFunctionBeginUser;
1116   *resize = PETSC_FALSE;
1117   PetscCall(VecGetDM(sol, &dm));
1118   PetscCall(DMGetDimension(dm, &dim));
1119   PetscCall(PetscFVSetLimiter(fvm, tctx->noneLimiter));
1120   PetscCall(PetscFVGetComputeGradients(fvm, &computeGradient));
1121   PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE));
1122   PetscCall(DMIsForest(dm, &isForest));
1123   PetscCall(DMConvert(dm, DMPLEX, &plex));
1124   PetscCall(DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM));
1125   PetscCall(DMCreateLocalVector(plex, &locX));
1126   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL));
1127   PetscCall(DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX));
1128   PetscCall(DMGlobalToLocalEnd(plex, sol, INSERT_VALUES, locX));
1129   PetscCall(DMCreateGlobalVector(gradDM, &grad));
1130   PetscCall(DMPlexReconstructGradientsFVM(plex, locX, grad));
1131   PetscCall(DMCreateLocalVector(gradDM, &locGrad));
1132   PetscCall(DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad));
1133   PetscCall(DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad));
1134   PetscCall(VecDestroy(&grad));
1135   PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &cStart, &cEnd));
1136   PetscCall(VecGetArrayRead(locGrad, &pointGrads));
1137   PetscCall(VecGetArrayRead(cellGeom, &pointGeom));
1138   PetscCall(VecGetArrayRead(locX, &pointVals));
1139   PetscCall(VecGetDM(cellGeom, &cellDM));
1140   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
1141   PetscCall(VecCreateFromOptions(PetscObjectComm((PetscObject)plex), NULL, 1, cEnd - cStart, PETSC_DETERMINE, &errVec));
1142   PetscCall(VecSetUp(errVec));
1143   PetscCall(VecGetArray(errVec, &errArray));
1144   for (c = cStart; c < cEnd; c++) {
1145     PetscReal        errInd = 0.;
1146     PetscScalar     *pointGrad;
1147     PetscScalar     *pointVal;
1148     PetscFVCellGeom *cg;
1149 
1150     PetscCall(DMPlexPointLocalRead(gradDM, c, pointGrads, &pointGrad));
1151     PetscCall(DMPlexPointLocalRead(cellDM, c, pointGeom, &cg));
1152     PetscCall(DMPlexPointLocalRead(plex, c, pointVals, &pointVal));
1153 
1154     PetscCall((user->model->errorIndicator)(dim, cg->volume, user->model->physics->dof, pointVal, pointGrad, &errInd, user->model->errorCtx));
1155     errArray[c - cStart] = errInd;
1156     minMaxInd[0]         = PetscMin(minMaxInd[0], errInd);
1157     minMaxInd[1]         = PetscMax(minMaxInd[1], errInd);
1158   }
1159   PetscCall(VecRestoreArray(errVec, &errArray));
1160   PetscCall(VecRestoreArrayRead(locX, &pointVals));
1161   PetscCall(VecRestoreArrayRead(cellGeom, &pointGeom));
1162   PetscCall(VecRestoreArrayRead(locGrad, &pointGrads));
1163   PetscCall(VecDestroy(&locGrad));
1164   PetscCall(VecDestroy(&locX));
1165   PetscCall(DMDestroy(&plex));
1166 
1167   PetscCall(VecTaggerComputeIS(refineTag, errVec, &refineIS, NULL));
1168   PetscCall(VecTaggerComputeIS(coarsenTag, errVec, &coarsenIS, NULL));
1169   PetscCall(ISGetSize(refineIS, &nRefine));
1170   PetscCall(ISGetSize(coarsenIS, &nCoarsen));
1171   if (nRefine) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS));
1172   if (nCoarsen) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_COARSEN, coarsenIS));
1173   PetscCall(ISDestroy(&coarsenIS));
1174   PetscCall(ISDestroy(&refineIS));
1175   PetscCall(VecDestroy(&errVec));
1176 
1177   PetscCall(PetscFVSetComputeGradients(fvm, computeGradient));
1178   PetscCall(PetscFVSetLimiter(fvm, tctx->limiter));
1179   minMaxInd[1] = -minMaxInd[1];
1180   PetscCall(MPIU_Allreduce(minMaxInd, minMaxIndGlobal, 2, MPIU_REAL, MPI_MIN, PetscObjectComm((PetscObject)dm)));
1181   PetscCall(PetscInfo(ts, "error indicator range (%E, %E)\n", (double)minMaxIndGlobal[0], (double)(-minMaxIndGlobal[1])));
1182   if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */
1183     PetscCall(DMAdaptLabel(dm, adaptLabel, &adaptedDM));
1184   }
1185   PetscCall(DMLabelDestroy(&adaptLabel));
1186   if (adaptedDM) {
1187     tctx->adaptedDM = adaptedDM;
1188     PetscCall(PetscInfo(ts, "Step %" PetscInt_FMT " adapted mesh, marking %" PetscInt_FMT " cells for refinement, and %" PetscInt_FMT " cells for coarsening\n", nstep, nRefine, nCoarsen));
1189     *resize = PETSC_TRUE;
1190   } else {
1191     PetscCall(PetscInfo(ts, "Step %" PetscInt_FMT " no adaptation\n", nstep));
1192   }
1193   PetscFunctionReturn(PETSC_SUCCESS);
1194 }
1195 
1196 static PetscErrorCode Transfer(TS ts, PetscInt nv, Vec vecsin[], Vec vecsout[], void *ctx)
1197 {
1198   TransferCtx *tctx = (TransferCtx *)ctx;
1199   DM           dm;
1200   PetscReal    time;
1201 
1202   PetscFunctionBeginUser;
1203   PetscCall(TSGetDM(ts, &dm));
1204   PetscCall(TSGetTime(ts, &time));
1205   PetscCheck(tctx->adaptedDM, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Missing adaptedDM");
1206   for (PetscInt i = 0; i < nv; i++) {
1207     const char *name;
1208 
1209     PetscCall(DMCreateGlobalVector(tctx->adaptedDM, &vecsout[i]));
1210     PetscCall(DMForestTransferVec(dm, vecsin[i], tctx->adaptedDM, vecsout[i], PETSC_TRUE, time));
1211     PetscCall(PetscObjectGetName((PetscObject)vecsin[i], &name));
1212     PetscCall(PetscObjectSetName((PetscObject)vecsout[i], name));
1213   }
1214   PetscCall(DMForestSetAdaptivityForest(tctx->adaptedDM, NULL)); /* clear internal references to the previous dm */
1215 
1216   Model     mod  = tctx->user->model;
1217   Physics   phys = mod->physics;
1218   PetscReal minRadius;
1219 
1220   PetscCall(DMPlexGetGeometryFVM(tctx->adaptedDM, NULL, NULL, &minRadius));
1221   PetscCall(MPIU_Allreduce(&phys->maxspeed, &mod->maxspeed, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts)));
1222   PetscCheck(mod->maxspeed > 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Physics did not set maxspeed");
1223 
1224   PetscReal dt = tctx->cfl * minRadius / mod->maxspeed;
1225   PetscCall(TSSetTimeStep(ts, dt));
1226 
1227   PetscCall(TSSetDM(ts, tctx->adaptedDM));
1228   PetscCall(DMDestroy(&tctx->adaptedDM));
1229 
1230   PetscFunctionReturn(PETSC_SUCCESS);
1231 }
1232 
1233 int main(int argc, char **argv)
1234 {
1235   MPI_Comm          comm;
1236   PetscDS           prob;
1237   PetscFV           fvm;
1238   PetscLimiter      limiter = NULL, noneLimiter = NULL;
1239   User              user;
1240   Model             mod;
1241   Physics           phys;
1242   DM                dm, plex;
1243   PetscReal         ftime, cfl, dt, minRadius;
1244   PetscInt          dim, nsteps;
1245   TS                ts;
1246   TSConvergedReason reason;
1247   Vec               X;
1248   PetscViewer       viewer;
1249   PetscBool         vtkCellGeom, useAMR;
1250   PetscInt          adaptInterval;
1251   char              physname[256] = "advect";
1252   VecTagger         refineTag = NULL, coarsenTag = NULL;
1253   TransferCtx       tctx;
1254 
1255   PetscFunctionBeginUser;
1256   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1257   comm = PETSC_COMM_WORLD;
1258 
1259   PetscCall(PetscNew(&user));
1260   PetscCall(PetscNew(&user->model));
1261   PetscCall(PetscNew(&user->model->physics));
1262   mod           = user->model;
1263   phys          = mod->physics;
1264   mod->comm     = comm;
1265   useAMR        = PETSC_FALSE;
1266   adaptInterval = 1;
1267 
1268   /* Register physical models to be available on the command line */
1269   PetscCall(PetscFunctionListAdd(&PhysicsList, "advect", PhysicsCreate_Advect));
1270   PetscCall(PetscFunctionListAdd(&PhysicsList, "sw", PhysicsCreate_SW));
1271   PetscCall(PetscFunctionListAdd(&PhysicsList, "euler", PhysicsCreate_Euler));
1272 
1273   PetscOptionsBegin(comm, NULL, "Unstructured Finite Volume Mesh Options", "");
1274   {
1275     cfl = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1276     PetscCall(PetscOptionsReal("-ufv_cfl", "CFL number per step", "", cfl, &cfl, NULL));
1277     user->vtkInterval = 1;
1278     PetscCall(PetscOptionsInt("-ufv_vtk_interval", "VTK output interval (0 to disable)", "", user->vtkInterval, &user->vtkInterval, NULL));
1279     user->vtkmon = PETSC_TRUE;
1280     PetscCall(PetscOptionsBool("-ufv_vtk_monitor", "Use VTKMonitor routine", "", user->vtkmon, &user->vtkmon, NULL));
1281     vtkCellGeom = PETSC_FALSE;
1282     PetscCall(PetscStrncpy(user->outputBasename, "ex11", sizeof(user->outputBasename)));
1283     PetscCall(PetscOptionsString("-ufv_vtk_basename", "VTK output basename", "", user->outputBasename, user->outputBasename, sizeof(user->outputBasename), NULL));
1284     PetscCall(PetscOptionsBool("-ufv_vtk_cellgeom", "Write cell geometry (for debugging)", "", vtkCellGeom, &vtkCellGeom, NULL));
1285     PetscCall(PetscOptionsBool("-ufv_use_amr", "use local adaptive mesh refinement", "", useAMR, &useAMR, NULL));
1286     PetscCall(PetscOptionsInt("-ufv_adapt_interval", "time steps between AMR", "", adaptInterval, &adaptInterval, NULL));
1287   }
1288   PetscOptionsEnd();
1289 
1290   if (useAMR) {
1291     VecTaggerBox refineBox, coarsenBox;
1292 
1293     refineBox.min = refineBox.max = PETSC_MAX_REAL;
1294     coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL;
1295 
1296     PetscCall(VecTaggerCreate(comm, &refineTag));
1297     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)refineTag, "refine_"));
1298     PetscCall(VecTaggerSetType(refineTag, VECTAGGERABSOLUTE));
1299     PetscCall(VecTaggerAbsoluteSetBox(refineTag, &refineBox));
1300     PetscCall(VecTaggerSetFromOptions(refineTag));
1301     PetscCall(VecTaggerSetUp(refineTag));
1302     PetscCall(PetscObjectViewFromOptions((PetscObject)refineTag, NULL, "-tag_view"));
1303 
1304     PetscCall(VecTaggerCreate(comm, &coarsenTag));
1305     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)coarsenTag, "coarsen_"));
1306     PetscCall(VecTaggerSetType(coarsenTag, VECTAGGERABSOLUTE));
1307     PetscCall(VecTaggerAbsoluteSetBox(coarsenTag, &coarsenBox));
1308     PetscCall(VecTaggerSetFromOptions(coarsenTag));
1309     PetscCall(VecTaggerSetUp(coarsenTag));
1310     PetscCall(PetscObjectViewFromOptions((PetscObject)coarsenTag, NULL, "-tag_view"));
1311   }
1312 
1313   PetscOptionsBegin(comm, NULL, "Unstructured Finite Volume Physics Options", "");
1314   {
1315     PetscErrorCode (*physcreate)(Model, Physics, PetscOptionItems *);
1316     PetscCall(PetscOptionsFList("-physics", "Physics module to solve", "", PhysicsList, physname, physname, sizeof physname, NULL));
1317     PetscCall(PetscFunctionListFind(PhysicsList, physname, &physcreate));
1318     PetscCall(PetscMemzero(phys, sizeof(struct _n_Physics)));
1319     PetscCall((*physcreate)(mod, phys, PetscOptionsObject));
1320     /* Count number of fields and dofs */
1321     for (phys->nfields = 0, phys->dof = 0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1322     PetscCheck(phys->dof > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set dof", physname);
1323     PetscCall(ModelFunctionalSetFromOptions(mod, PetscOptionsObject));
1324   }
1325   PetscOptionsEnd();
1326 
1327   /* Create mesh */
1328   {
1329     PetscInt i;
1330 
1331     PetscCall(DMCreate(comm, &dm));
1332     PetscCall(DMSetType(dm, DMPLEX));
1333     PetscCall(DMSetFromOptions(dm));
1334     for (i = 0; i < DIM; i++) {
1335       mod->bounds[2 * i]     = 0.;
1336       mod->bounds[2 * i + 1] = 1.;
1337     };
1338     dim = DIM;
1339     { /* a null name means just do a hex box */
1340       PetscInt  cells[3] = {1, 1, 1}, n = 3;
1341       PetscBool flg2, skew              = PETSC_FALSE;
1342       PetscInt  nret2 = 2 * DIM;
1343       PetscOptionsBegin(comm, NULL, "Rectangular mesh options", "");
1344       PetscCall(PetscOptionsRealArray("-grid_bounds", "bounds of the mesh in each direction (i.e., x_min,x_max,y_min,y_max", "", mod->bounds, &nret2, &flg2));
1345       PetscCall(PetscOptionsBool("-grid_skew_60", "Skew grid for 60 degree shock mesh", "", skew, &skew, NULL));
1346       PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", cells, &n, NULL));
1347       PetscOptionsEnd();
1348       /* TODO Rewrite this with Mark, and remove grid_bounds at that time */
1349       if (flg2) {
1350         PetscInt     dimEmbed, i;
1351         PetscInt     nCoords;
1352         PetscScalar *coords;
1353         Vec          coordinates;
1354 
1355         PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1356         PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
1357         PetscCall(VecGetLocalSize(coordinates, &nCoords));
1358         PetscCheck(!(nCoords % dimEmbed), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Coordinate vector the wrong size");
1359         PetscCall(VecGetArray(coordinates, &coords));
1360         for (i = 0; i < nCoords; i += dimEmbed) {
1361           PetscInt j;
1362 
1363           PetscScalar *coord = &coords[i];
1364           for (j = 0; j < dimEmbed; j++) {
1365             coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]);
1366             if (dim == 2 && cells[1] == 1 && j == 0 && skew) {
1367               if (cells[0] == 2 && i == 8) {
1368                 coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */
1369               } else if (cells[0] == 3) {
1370                 if (i == 2 || i == 10) coord[j] = mod->bounds[1] / 4.;
1371                 else if (i == 4) coord[j] = mod->bounds[1] / 2.;
1372                 else if (i == 12) coord[j] = 1.57735026918963 * mod->bounds[1] / 2.;
1373               }
1374             }
1375           }
1376         }
1377         PetscCall(VecRestoreArray(coordinates, &coords));
1378         PetscCall(DMSetCoordinatesLocal(dm, coordinates));
1379       }
1380     }
1381   }
1382   PetscCall(DMViewFromOptions(dm, NULL, "-orig_dm_view"));
1383   PetscCall(DMGetDimension(dm, &dim));
1384 
1385   /* set up BCs, functions, tags */
1386   PetscCall(DMCreateLabel(dm, "Face Sets"));
1387   mod->errorIndicator = ErrorIndicator_Simple;
1388 
1389   {
1390     DM gdm;
1391 
1392     PetscCall(DMPlexConstructGhostCells(dm, NULL, NULL, &gdm));
1393     PetscCall(DMDestroy(&dm));
1394     dm = gdm;
1395     PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1396   }
1397 
1398   PetscCall(PetscFVCreate(comm, &fvm));
1399   PetscCall(PetscFVSetFromOptions(fvm));
1400   PetscCall(PetscFVSetNumComponents(fvm, phys->dof));
1401   PetscCall(PetscFVSetSpatialDimension(fvm, dim));
1402   PetscCall(PetscObjectSetName((PetscObject)fvm, ""));
1403   {
1404     PetscInt f, dof;
1405     for (f = 0, dof = 0; f < phys->nfields; f++) {
1406       PetscInt newDof = phys->field_desc[f].dof;
1407 
1408       if (newDof == 1) {
1409         PetscCall(PetscFVSetComponentName(fvm, dof, phys->field_desc[f].name));
1410       } else {
1411         PetscInt j;
1412 
1413         for (j = 0; j < newDof; j++) {
1414           char compName[256] = "Unknown";
1415 
1416           PetscCall(PetscSNPrintf(compName, sizeof(compName), "%s_%" PetscInt_FMT, phys->field_desc[f].name, j));
1417           PetscCall(PetscFVSetComponentName(fvm, dof + j, compName));
1418         }
1419       }
1420       dof += newDof;
1421     }
1422   }
1423   /* FV is now structured with one field having all physics as components */
1424   PetscCall(DMAddField(dm, NULL, (PetscObject)fvm));
1425   PetscCall(DMCreateDS(dm));
1426   PetscCall(DMGetDS(dm, &prob));
1427   PetscCall(PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann));
1428   PetscCall(PetscDSSetContext(prob, 0, user->model->physics));
1429   PetscCall((*mod->setupbc)(dm, prob, phys));
1430   PetscCall(PetscDSSetFromOptions(prob));
1431   {
1432     char      convType[256];
1433     PetscBool flg;
1434 
1435     PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
1436     PetscCall(PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex12", DMList, DMPLEX, convType, 256, &flg));
1437     PetscOptionsEnd();
1438     if (flg) {
1439       DM dmConv;
1440 
1441       PetscCall(DMConvert(dm, convType, &dmConv));
1442       if (dmConv) {
1443         PetscCall(DMViewFromOptions(dmConv, NULL, "-dm_conv_view"));
1444         PetscCall(DMDestroy(&dm));
1445         dm = dmConv;
1446         PetscCall(DMSetFromOptions(dm));
1447       }
1448     }
1449   }
1450 #ifdef PETSC_HAVE_LIBCEED
1451   {
1452     PetscBool useCeed;
1453     PetscCall(DMPlexGetUseCeed(dm, &useCeed));
1454     if (useCeed) PetscCall((*user->model->setupCEED)(dm, user->model->physics));
1455   }
1456 #endif
1457 
1458   PetscCall(initializeTS(dm, user, &ts));
1459 
1460   PetscCall(DMCreateGlobalVector(dm, &X));
1461   PetscCall(PetscObjectSetName((PetscObject)X, "solution"));
1462   PetscCall(SetInitialCondition(dm, X, user));
1463   if (useAMR) {
1464     PetscInt adaptIter;
1465 
1466     /* use no limiting when reconstructing gradients for adaptivity */
1467     PetscCall(PetscFVGetLimiter(fvm, &limiter));
1468     PetscCall(PetscObjectReference((PetscObject)limiter));
1469     PetscCall(PetscLimiterCreate(PetscObjectComm((PetscObject)fvm), &noneLimiter));
1470     PetscCall(PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE));
1471 
1472     /* Refinement context */
1473     tctx.fvm         = fvm;
1474     tctx.refineTag   = refineTag;
1475     tctx.coarsenTag  = coarsenTag;
1476     tctx.adaptedDM   = NULL;
1477     tctx.user        = user;
1478     tctx.noneLimiter = noneLimiter;
1479     tctx.limiter     = limiter;
1480     tctx.cfl         = cfl;
1481 
1482     /* Do some initial refinement steps */
1483     for (adaptIter = 0;; ++adaptIter) {
1484       PetscLogDouble bytes;
1485       PetscBool      resize;
1486 
1487       PetscCall(PetscMemoryGetCurrentUsage(&bytes));
1488       PetscCall(PetscInfo(ts, "refinement loop %" PetscInt_FMT ": memory used %g\n", adaptIter, (double)bytes));
1489       PetscCall(DMViewFromOptions(dm, NULL, "-initial_dm_view"));
1490       PetscCall(VecViewFromOptions(X, NULL, "-initial_vec_view"));
1491 
1492       PetscCall(adaptToleranceFVMSetUp(ts, -1, 0.0, X, &resize, &tctx));
1493       if (!resize) break;
1494       PetscCall(DMDestroy(&dm));
1495       PetscCall(VecDestroy(&X));
1496       dm             = tctx.adaptedDM;
1497       tctx.adaptedDM = NULL;
1498       PetscCall(TSSetDM(ts, dm));
1499       PetscCall(DMCreateGlobalVector(dm, &X));
1500       PetscCall(PetscObjectSetName((PetscObject)X, "solution"));
1501       PetscCall(SetInitialCondition(dm, X, user));
1502     }
1503   }
1504 
1505   PetscCall(DMConvert(dm, DMPLEX, &plex));
1506   if (vtkCellGeom) {
1507     DM  dmCell;
1508     Vec cellgeom, partition;
1509 
1510     PetscCall(DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL));
1511     PetscCall(OutputVTK(dm, "ex11-cellgeom.vtk", &viewer));
1512     PetscCall(VecView(cellgeom, viewer));
1513     PetscCall(PetscViewerDestroy(&viewer));
1514     PetscCall(CreatePartitionVec(dm, &dmCell, &partition));
1515     PetscCall(OutputVTK(dmCell, "ex11-partition.vtk", &viewer));
1516     PetscCall(VecView(partition, viewer));
1517     PetscCall(PetscViewerDestroy(&viewer));
1518     PetscCall(VecDestroy(&partition));
1519     PetscCall(DMDestroy(&dmCell));
1520   }
1521   /* collect max maxspeed from all processes -- todo */
1522   PetscCall(DMPlexGetGeometryFVM(plex, NULL, NULL, &minRadius));
1523   PetscCall(DMDestroy(&plex));
1524   PetscCall(MPIU_Allreduce(&phys->maxspeed, &mod->maxspeed, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts)));
1525   PetscCheck(mod->maxspeed > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set maxspeed", physname);
1526   dt = cfl * minRadius / mod->maxspeed;
1527   PetscCall(TSSetTimeStep(ts, dt));
1528   PetscCall(TSSetFromOptions(ts));
1529 
1530   /* When using adaptive mesh refinement
1531      specify callbacks to refine the solution
1532      and interpolate data from old to new mesh */
1533   if (useAMR) { PetscCall(TSSetResize(ts, adaptToleranceFVMSetUp, Transfer, &tctx)); }
1534   PetscCall(TSSetSolution(ts, X));
1535   PetscCall(VecDestroy(&X));
1536   PetscCall(TSSolve(ts, NULL));
1537   PetscCall(TSGetSolveTime(ts, &ftime));
1538   PetscCall(TSGetStepNumber(ts, &nsteps));
1539   PetscCall(TSGetConvergedReason(ts, &reason));
1540   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, nsteps));
1541   PetscCall(TSDestroy(&ts));
1542 
1543   PetscCall(VecTaggerDestroy(&refineTag));
1544   PetscCall(VecTaggerDestroy(&coarsenTag));
1545   PetscCall(PetscFunctionListDestroy(&PhysicsList));
1546   PetscCall(PetscFunctionListDestroy(&PhysicsRiemannList_SW));
1547   PetscCall(PetscFunctionListDestroy(&PhysicsRiemannList_Euler));
1548   PetscCall(FunctionalLinkDestroy(&user->model->functionalRegistry));
1549   PetscCall(PetscFree(user->model->functionalMonitored));
1550   PetscCall(PetscFree(user->model->functionalCall));
1551   PetscCall(PetscFree(user->model->physics->data));
1552   PetscCall(PetscFree(user->model->physics));
1553   PetscCall(PetscFree(user->model));
1554   PetscCall(PetscFree(user));
1555   PetscCall(PetscLimiterDestroy(&limiter));
1556   PetscCall(PetscLimiterDestroy(&noneLimiter));
1557   PetscCall(PetscFVDestroy(&fvm));
1558   PetscCall(DMDestroy(&dm));
1559   PetscCall(PetscFinalize());
1560   return 0;
1561 }
1562 
1563 /* Subroutine to set up the initial conditions for the */
1564 /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */
1565 /* ----------------------------------------------------------------------- */
1566 int projecteqstate(PetscReal wc[], const PetscReal ueq[], PetscReal lv[][3])
1567 {
1568   int j, k;
1569   /*      Wc=matmul(lv,Ueq) 3 vars */
1570   for (k = 0; k < 3; ++k) {
1571     wc[k] = 0.;
1572     for (j = 0; j < 3; ++j) wc[k] += lv[k][j] * ueq[j];
1573   }
1574   return 0;
1575 }
1576 /* ----------------------------------------------------------------------- */
1577 int projecttoprim(PetscReal v[], const PetscReal wc[], PetscReal rv[][3])
1578 {
1579   int k, j;
1580   /*      V=matmul(rv,WC) 3 vars */
1581   for (k = 0; k < 3; ++k) {
1582     v[k] = 0.;
1583     for (j = 0; j < 3; ++j) v[k] += rv[k][j] * wc[j];
1584   }
1585   return 0;
1586 }
1587 /* ---------------------------------------------------------------------- */
1588 int eigenvectors(PetscReal rv[][3], PetscReal lv[][3], const PetscReal ueq[], PetscReal gamma)
1589 {
1590   int       j, k;
1591   PetscReal rho, csnd, p0;
1592   /* PetscScalar u; */
1593 
1594   for (k = 0; k < 3; ++k)
1595     for (j = 0; j < 3; ++j) {
1596       lv[k][j] = 0.;
1597       rv[k][j] = 0.;
1598     }
1599   rho = ueq[0];
1600   /* u = ueq[1]; */
1601   p0       = ueq[2];
1602   csnd     = PetscSqrtReal(gamma * p0 / rho);
1603   lv[0][1] = rho * .5;
1604   lv[0][2] = -.5 / csnd;
1605   lv[1][0] = csnd;
1606   lv[1][2] = -1. / csnd;
1607   lv[2][1] = rho * .5;
1608   lv[2][2] = .5 / csnd;
1609   rv[0][0] = -1. / csnd;
1610   rv[1][0] = 1. / rho;
1611   rv[2][0] = -csnd;
1612   rv[0][1] = 1. / csnd;
1613   rv[0][2] = 1. / csnd;
1614   rv[1][2] = 1. / rho;
1615   rv[2][2] = csnd;
1616   return 0;
1617 }
1618 
1619 int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx)
1620 {
1621   PetscReal p0, u0, wcp[3], wc[3];
1622   PetscReal lv[3][3];
1623   PetscReal vp[3];
1624   PetscReal rv[3][3];
1625   PetscReal eps, ueq[3], rho0, twopi;
1626 
1627   /* Function Body */
1628   twopi  = 2. * PETSC_PI;
1629   eps    = 1e-4;    /* perturbation */
1630   rho0   = 1e3;     /* density of water */
1631   p0     = 101325.; /* init pressure of 1 atm (?) */
1632   u0     = 0.;
1633   ueq[0] = rho0;
1634   ueq[1] = u0;
1635   ueq[2] = p0;
1636   /* Project initial state to characteristic variables */
1637   eigenvectors(rv, lv, ueq, gamma);
1638   projecteqstate(wc, ueq, lv);
1639   wcp[0] = wc[0];
1640   wcp[1] = wc[1];
1641   wcp[2] = wc[2] + eps * PetscCosReal(coord[0] * 2. * twopi / Lx);
1642   projecttoprim(vp, wcp, rv);
1643   ux->r     = vp[0];         /* density */
1644   ux->ru[0] = vp[0] * vp[1]; /* x momentum */
1645   ux->ru[1] = 0.;
1646 #if defined DIM > 2
1647   if (dim > 2) ux->ru[2] = 0.;
1648 #endif
1649   /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */
1650   ux->E = vp[2] / (gamma - 1.) + 0.5 * vp[0] * vp[1] * vp[1];
1651   return 0;
1652 }
1653 
1654 /*TEST
1655 
1656   testset:
1657     args: -dm_plex_adj_cone -dm_plex_adj_closure 0
1658 
1659     test:
1660       suffix: adv_2d_tri_0
1661       requires: triangle
1662       TODO: how did this ever get in main when there is no support for this
1663       args: -ufv_vtk_interval 0 -simplex -dm_refine 3 -dm_plex_faces 1,1 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
1664 
1665     test:
1666       suffix: adv_2d_tri_1
1667       requires: triangle
1668       TODO: how did this ever get in main when there is no support for this
1669       args: -ufv_vtk_interval 0 -simplex -dm_refine 5 -dm_plex_faces 1,1 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1
1670 
1671     test:
1672       suffix: tut_1
1673       requires: exodusii
1674       nsize: 1
1675       args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -
1676 
1677     test:
1678       suffix: tut_2
1679       requires: exodusii
1680       nsize: 1
1681       args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw
1682 
1683     test:
1684       suffix: tut_3
1685       requires: exodusii
1686       nsize: 4
1687       args: -dm_distribute_overlap 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -monitor Error -advect_sol_type bump -petscfv_type leastsquares -petsclimiter_type sin
1688 
1689     test:
1690       suffix: tut_4
1691       requires: exodusii
1692       nsize: 4
1693       args: -dm_distribute_overlap 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -physics sw -monitor Height,Energy -petscfv_type leastsquares -petsclimiter_type minmod
1694 
1695   testset:
1696     args: -dm_plex_adj_cone -dm_plex_adj_closure 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1
1697 
1698     # 2D Advection 0-10
1699     test:
1700       suffix: 0
1701       requires: exodusii
1702       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
1703 
1704     test:
1705       suffix: 1
1706       requires: exodusii
1707       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
1708 
1709     test:
1710       suffix: 2
1711       requires: exodusii
1712       nsize: 2
1713       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
1714 
1715     test:
1716       suffix: 3
1717       requires: exodusii
1718       nsize: 2
1719       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
1720 
1721     test:
1722       suffix: 4
1723       requires: exodusii
1724       nsize: 4
1725       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -petscpartitioner_type simple
1726 
1727     test:
1728       suffix: 5
1729       requires: exodusii
1730       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1
1731 
1732     test:
1733       suffix: 7
1734       requires: exodusii
1735       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
1736 
1737     test:
1738       suffix: 8
1739       requires: exodusii
1740       nsize: 2
1741       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
1742 
1743     test:
1744       suffix: 9
1745       requires: exodusii
1746       nsize: 8
1747       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
1748 
1749     test:
1750       suffix: 10
1751       requires: exodusii
1752       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo
1753 
1754   # 2D Shallow water
1755   testset:
1756     args: -physics sw -ufv_vtk_interval 0 -dm_plex_adj_cone -dm_plex_adj_closure 0
1757 
1758     test:
1759       suffix: sw_0
1760       requires: exodusii
1761       args: -bc_wall 100,101 -ufv_cfl 5 -petscfv_type leastsquares -petsclimiter_type sin \
1762             -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo \
1763             -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
1764             -monitor height,energy
1765 
1766     test:
1767       suffix: sw_ceed
1768       requires: exodusii libceed
1769       args: -sw_riemann rusanov_ceed -bc_wall 100,101 -ufv_cfl 5 -petsclimiter_type sin \
1770             -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -dm_plex_use_ceed \
1771             -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
1772             -monitor height,energy
1773 
1774     test:
1775       suffix: sw_ceed_small
1776       requires: exodusii libceed
1777       args: -sw_riemann rusanov_ceed -bc_wall 1,3 -ufv_cfl 5 -petsclimiter_type sin -dm_plex_use_ceed \
1778             -dm_plex_shape annulus -dm_plex_simplex 0 -dm_plex_box_lower 0,1 -dm_plex_box_upper 6.28,3 -dm_plex_box_faces 8,2 \
1779             -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
1780             -monitor height,energy
1781 
1782     test:
1783       suffix: sw_1
1784       nsize: 2
1785       args: -bc_wall 1,3 -ufv_cfl 5 -petsclimiter_type sin \
1786             -dm_plex_shape annulus -dm_plex_simplex 0 -dm_plex_box_faces 24,12 -dm_plex_box_lower 0,1 -dm_plex_box_upper 1,3 -dm_distribute_overlap 1 \
1787             -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
1788             -monitor height,energy
1789 
1790     test:
1791       suffix: sw_hll
1792       args: -sw_riemann hll -bc_wall 1,2,3,4 -ufv_cfl 3 -petscfv_type leastsquares -petsclimiter_type sin \
1793             -grid_bounds 0,5,0,5 -dm_plex_simplex 0 -dm_plex_box_faces 25,25 \
1794             -ts_max_steps 5 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
1795             -monitor height,energy
1796 
1797   # 2D Euler
1798   testset:
1799     args: -physics euler -eu_type linear_wave -eu_gamma 1.4 -dm_plex_adj_cone -dm_plex_adj_closure 0 \
1800           -ufv_vtk_interval 0 -ufv_vtk_basename ${wPETSC_DIR}/ex11 -monitor density,energy
1801 
1802     test:
1803       suffix: euler_0
1804       requires: exodusii !complex
1805       args: -eu_riemann godunov -bc_wall 100,101 -ufv_cfl 5 -petsclimiter_type sin \
1806             -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo \
1807             -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10
1808 
1809     test:
1810       suffix: euler_ceed
1811       requires: exodusii libceed
1812       args: -eu_riemann godunov_ceed -bc_wall 100,101 -ufv_cfl 5 -petsclimiter_type sin \
1813             -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -dm_plex_use_ceed \
1814             -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10
1815 
1816   testset:
1817     args: -dm_plex_adj_cone -dm_plex_adj_closure 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1
1818 
1819     # 2D Advection: p4est
1820     test:
1821       suffix: p4est_advec_2d
1822       requires: p4est
1823       args: -ufv_vtk_interval 0 -dm_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 2 -dm_p4est_refine_pattern hash -dm_forest_maximum_refinement 5
1824 
1825     # Advection in a box
1826     test:
1827       suffix: adv_2d_quad_0
1828       args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
1829 
1830     test:
1831       suffix: adv_2d_quad_1
1832       args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1
1833       timeoutfactor: 3
1834 
1835     test:
1836       suffix: adv_2d_quad_p4est_0
1837       requires: p4est
1838       args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
1839 
1840     test:
1841       suffix: adv_2d_quad_p4est_1
1842       requires: p4est
1843       args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1
1844       timeoutfactor: 3
1845 
1846     test: # broken for quad precision
1847       suffix: adv_2d_quad_p4est_adapt_0
1848       requires: p4est !__float128
1849       args: -ufv_vtk_interval 0 -dm_refine 3 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 -ufv_use_amr -refine_vec_tagger_box 0.005,inf -coarsen_vec_tagger_box 0,1.e-5 -petscfv_type leastsquares -ts_max_time 0.01
1850       timeoutfactor: 3
1851 
1852     test:
1853       suffix: adv_0
1854       requires: exodusii
1855       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -bc_inflow 100,101,200 -bc_outflow 201
1856 
1857     test:
1858       suffix: shock_0
1859       requires: p4est !single !complex
1860       args: -dm_plex_box_faces 2,1 -grid_bounds -1,1.,0.,1 -grid_skew_60 \
1861       -dm_type p4est -dm_forest_partition_overlap 1 -dm_forest_maximum_refinement 6 -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 \
1862       -ufv_use_amr -refine_vec_tagger_box 0.5,inf -coarsen_vec_tagger_box 0,1.e-2 -refine_tag_view -coarsen_tag_view \
1863       -bc_wall 1,2,3,4 -physics euler -eu_type iv_shock -ufv_cfl 10 -eu_alpha 60. -eu_gamma 1.4 -eu_amach 2.02 -eu_rho2 3. \
1864       -petscfv_type leastsquares -petsclimiter_type minmod -petscfv_compute_gradients 0 \
1865       -ts_max_time 0.5 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
1866       -ufv_vtk_basename ${wPETSC_DIR}/ex11 -ufv_vtk_interval 0 -monitor density,energy
1867       timeoutfactor: 3
1868 
1869     # Test GLVis visualization of PetscFV fields
1870     test:
1871       suffix: glvis_adv_2d_tet
1872       args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 \
1873             -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0 \
1874             -ts_monitor_solution glvis: -ts_max_steps 0
1875 
1876     test:
1877       suffix: glvis_adv_2d_quad
1878       args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 -bc_inflow 1,2,4 -bc_outflow 3 \
1879             -dm_refine 5 -dm_plex_separate_marker \
1880             -ts_monitor_solution glvis: -ts_max_steps 0
1881 
1882 TEST*/
1883