xref: /petsc/src/ts/tutorials/ex11.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
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 F*/
37 #include <petscdmplex.h>
38 #include <petscdmforest.h>
39 #include <petscds.h>
40 #include <petscts.h>
41 
42 #define DIM 2 /* Geometric dimension */
43 
44 static PetscFunctionList PhysicsList, PhysicsRiemannList_SW;
45 
46 /* Represents continuum physical equations. */
47 typedef struct _n_Physics *Physics;
48 
49 /* Physical model includes boundary conditions, initial conditions, and functionals of interest. It is
50  * discretization-independent, but its members depend on the scenario being solved. */
51 typedef struct _n_Model *Model;
52 
53 /* 'User' implements a discretization of a continuous model. */
54 typedef struct _n_User *User;
55 typedef PetscErrorCode (*SolutionFunction)(Model, PetscReal, const PetscReal *, PetscScalar *, void *);
56 typedef PetscErrorCode (*SetUpBCFunction)(DM, PetscDS, Physics);
57 typedef PetscErrorCode (*FunctionalFunction)(Model, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *);
58 typedef PetscErrorCode (*SetupFields)(Physics, PetscSection);
59 static PetscErrorCode ModelSolutionSetDefault(Model, SolutionFunction, void *);
60 static PetscErrorCode ModelFunctionalRegister(Model, const char *, PetscInt *, FunctionalFunction, void *);
61 static PetscErrorCode OutputVTK(DM, const char *, PetscViewer *);
62 
63 struct FieldDescription {
64   const char *name;
65   PetscInt    dof;
66 };
67 
68 typedef struct _n_FunctionalLink *FunctionalLink;
69 struct _n_FunctionalLink {
70   char              *name;
71   FunctionalFunction func;
72   void              *ctx;
73   PetscInt           offset;
74   FunctionalLink     next;
75 };
76 
77 struct _n_Physics {
78   PetscRiemannFunc               riemann;
79   PetscInt                       dof;      /* number of degrees of freedom per cell */
80   PetscReal                      maxspeed; /* kludge to pick initial time step, need to add monitoring and step control */
81   void                          *data;
82   PetscInt                       nfields;
83   const struct FieldDescription *field_desc;
84 };
85 
86 struct _n_Model {
87   MPI_Comm         comm; /* Does not do collective communicaton, but some error conditions can be collective */
88   Physics          physics;
89   FunctionalLink   functionalRegistry;
90   PetscInt         maxComputed;
91   PetscInt         numMonitored;
92   FunctionalLink  *functionalMonitored;
93   PetscInt         numCall;
94   FunctionalLink  *functionalCall;
95   SolutionFunction solution;
96   SetUpBCFunction  setupbc;
97   void            *solutionctx;
98   PetscReal        maxspeed; /* estimate of global maximum speed (for CFL calculation) */
99   PetscReal        bounds[2 * DIM];
100   PetscErrorCode (*errorIndicator)(PetscInt, PetscReal, PetscInt, const PetscScalar[], const PetscScalar[], PetscReal *, void *);
101   void *errorCtx;
102 };
103 
104 struct _n_User {
105   PetscInt  vtkInterval;                        /* For monitor */
106   char      outputBasename[PETSC_MAX_PATH_LEN]; /* Basename for output files */
107   PetscInt  monitorStepOffset;
108   Model     model;
109   PetscBool vtkmon;
110 };
111 
112 static inline PetscReal DotDIMReal(const PetscReal *x, const PetscReal *y) {
113   PetscInt  i;
114   PetscReal prod = 0.0;
115 
116   for (i = 0; i < DIM; i++) prod += x[i] * y[i];
117   return prod;
118 }
119 static inline PetscReal NormDIM(const PetscReal *x) {
120   return PetscSqrtReal(PetscAbsReal(DotDIMReal(x, x)));
121 }
122 
123 static inline PetscReal Dot2Real(const PetscReal *x, const PetscReal *y) {
124   return x[0] * y[0] + x[1] * y[1];
125 }
126 static inline PetscReal Norm2Real(const PetscReal *x) {
127   return PetscSqrtReal(PetscAbsReal(Dot2Real(x, x)));
128 }
129 static inline void Normalize2Real(PetscReal *x) {
130   PetscReal a = 1. / Norm2Real(x);
131   x[0] *= a;
132   x[1] *= a;
133 }
134 static inline void Waxpy2Real(PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w) {
135   w[0] = a * x[0] + y[0];
136   w[1] = a * x[1] + y[1];
137 }
138 static inline void Scale2Real(PetscReal a, const PetscReal *x, PetscReal *y) {
139   y[0] = a * x[0];
140   y[1] = a * x[1];
141 }
142 
143 /******************* Advect ********************/
144 typedef enum {
145   ADVECT_SOL_TILTED,
146   ADVECT_SOL_BUMP,
147   ADVECT_SOL_BUMP_CAVITY
148 } AdvectSolType;
149 static const char *const AdvectSolTypes[] = {"TILTED", "BUMP", "BUMP_CAVITY", "AdvectSolType", "ADVECT_SOL_", 0};
150 typedef enum {
151   ADVECT_SOL_BUMP_CONE,
152   ADVECT_SOL_BUMP_COS
153 } AdvectSolBumpType;
154 static const char *const AdvectSolBumpTypes[] = {"CONE", "COS", "AdvectSolBumpType", "ADVECT_SOL_BUMP_", 0};
155 
156 typedef struct {
157   PetscReal wind[DIM];
158 } Physics_Advect_Tilted;
159 typedef struct {
160   PetscReal         center[DIM];
161   PetscReal         radius;
162   AdvectSolBumpType type;
163 } Physics_Advect_Bump;
164 
165 typedef struct {
166   PetscReal     inflowState;
167   AdvectSolType soltype;
168   union
169   {
170     Physics_Advect_Tilted tilted;
171     Physics_Advect_Bump   bump;
172   } sol;
173   struct {
174     PetscInt Solution;
175     PetscInt Error;
176   } functional;
177 } Physics_Advect;
178 
179 static const struct FieldDescription PhysicsFields_Advect[] = {
180   {"U",  1},
181   {NULL, 0}
182 };
183 
184 static PetscErrorCode PhysicsBoundary_Advect_Inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) {
185   Physics         phys   = (Physics)ctx;
186   Physics_Advect *advect = (Physics_Advect *)phys->data;
187 
188   PetscFunctionBeginUser;
189   xG[0] = advect->inflowState;
190   PetscFunctionReturn(0);
191 }
192 
193 static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) {
194   PetscFunctionBeginUser;
195   xG[0] = xI[0];
196   PetscFunctionReturn(0);
197 }
198 
199 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) {
200   Physics_Advect *advect = (Physics_Advect *)phys->data;
201   PetscReal       wind[DIM], wn;
202 
203   switch (advect->soltype) {
204   case ADVECT_SOL_TILTED: {
205     Physics_Advect_Tilted *tilted = &advect->sol.tilted;
206     wind[0]                       = tilted->wind[0];
207     wind[1]                       = tilted->wind[1];
208   } break;
209   case ADVECT_SOL_BUMP:
210     wind[0] = -qp[1];
211     wind[1] = qp[0];
212     break;
213   case ADVECT_SOL_BUMP_CAVITY: {
214     PetscInt  i;
215     PetscReal comp2[3] = {0., 0., 0.}, rad2;
216 
217     rad2 = 0.;
218     for (i = 0; i < dim; i++) {
219       comp2[i] = qp[i] * qp[i];
220       rad2 += comp2[i];
221     }
222 
223     wind[0] = -qp[1];
224     wind[1] = qp[0];
225     if (rad2 > 1.) {
226       PetscInt  maxI     = 0;
227       PetscReal maxComp2 = comp2[0];
228 
229       for (i = 1; i < dim; i++) {
230         if (comp2[i] > maxComp2) {
231           maxI     = i;
232           maxComp2 = comp2[i];
233         }
234       }
235       wind[maxI] = 0.;
236     }
237   } break;
238   default: {
239     PetscInt i;
240     for (i = 0; i < DIM; ++i) wind[i] = 0.0;
241   }
242     /* default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for solution type %s",AdvectSolBumpTypes[advect->soltype]); */
243   }
244   wn      = Dot2Real(wind, n);
245   flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn;
246 }
247 
248 static PetscErrorCode PhysicsSolution_Advect(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) {
249   Physics         phys   = (Physics)ctx;
250   Physics_Advect *advect = (Physics_Advect *)phys->data;
251 
252   PetscFunctionBeginUser;
253   switch (advect->soltype) {
254   case ADVECT_SOL_TILTED: {
255     PetscReal              x0[DIM];
256     Physics_Advect_Tilted *tilted = &advect->sol.tilted;
257     Waxpy2Real(-time, tilted->wind, x, x0);
258     if (x0[1] > 0) u[0] = 1. * x[0] + 3. * x[1];
259     else u[0] = advect->inflowState;
260   } break;
261   case ADVECT_SOL_BUMP_CAVITY:
262   case ADVECT_SOL_BUMP: {
263     Physics_Advect_Bump *bump = &advect->sol.bump;
264     PetscReal            x0[DIM], v[DIM], r, cost, sint;
265     cost  = PetscCosReal(time);
266     sint  = PetscSinReal(time);
267     x0[0] = cost * x[0] + sint * x[1];
268     x0[1] = -sint * x[0] + cost * x[1];
269     Waxpy2Real(-1, bump->center, x0, v);
270     r = Norm2Real(v);
271     switch (bump->type) {
272     case ADVECT_SOL_BUMP_CONE: u[0] = PetscMax(1 - r / bump->radius, 0); break;
273     case ADVECT_SOL_BUMP_COS: u[0] = 0.5 + 0.5 * PetscCosReal(PetscMin(r / bump->radius, 1) * PETSC_PI); break;
274     }
275   } break;
276   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type");
277   }
278   PetscFunctionReturn(0);
279 }
280 
281 static PetscErrorCode PhysicsFunctional_Advect(Model mod, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx) {
282   Physics         phys      = (Physics)ctx;
283   Physics_Advect *advect    = (Physics_Advect *)phys->data;
284   PetscScalar     yexact[1] = {0.0};
285 
286   PetscFunctionBeginUser;
287   PetscCall(PhysicsSolution_Advect(mod, time, x, yexact, phys));
288   f[advect->functional.Solution] = PetscRealPart(y[0]);
289   f[advect->functional.Error]    = PetscAbsScalar(y[0] - yexact[0]);
290   PetscFunctionReturn(0);
291 }
292 
293 static PetscErrorCode SetUpBC_Advect(DM dm, PetscDS prob, Physics phys) {
294   const PetscInt inflowids[] = {100, 200, 300}, outflowids[] = {101};
295   DMLabel        label;
296 
297   PetscFunctionBeginUser;
298   /* Register "canned" boundary conditions and defaults for where to apply. */
299   PetscCall(DMGetLabel(dm, "Face Sets", &label));
300   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));
301   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));
302   PetscFunctionReturn(0);
303 }
304 
305 static PetscErrorCode PhysicsCreate_Advect(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject) {
306   Physics_Advect *advect;
307 
308   PetscFunctionBeginUser;
309   phys->field_desc = PhysicsFields_Advect;
310   phys->riemann    = (PetscRiemannFunc)PhysicsRiemann_Advect;
311   PetscCall(PetscNew(&advect));
312   phys->data   = advect;
313   mod->setupbc = SetUpBC_Advect;
314 
315   PetscOptionsHeadBegin(PetscOptionsObject, "Advect options");
316   {
317     PetscInt two = 2, dof = 1;
318     advect->soltype = ADVECT_SOL_TILTED;
319     PetscCall(PetscOptionsEnum("-advect_sol_type", "solution type", "", AdvectSolTypes, (PetscEnum)advect->soltype, (PetscEnum *)&advect->soltype, NULL));
320     switch (advect->soltype) {
321     case ADVECT_SOL_TILTED: {
322       Physics_Advect_Tilted *tilted = &advect->sol.tilted;
323       two                           = 2;
324       tilted->wind[0]               = 0.0;
325       tilted->wind[1]               = 1.0;
326       PetscCall(PetscOptionsRealArray("-advect_tilted_wind", "background wind vx,vy", "", tilted->wind, &two, NULL));
327       advect->inflowState = -2.0;
328       PetscCall(PetscOptionsRealArray("-advect_tilted_inflow", "Inflow state", "", &advect->inflowState, &dof, NULL));
329       phys->maxspeed = Norm2Real(tilted->wind);
330     } break;
331     case ADVECT_SOL_BUMP_CAVITY:
332     case ADVECT_SOL_BUMP: {
333       Physics_Advect_Bump *bump = &advect->sol.bump;
334       two                       = 2;
335       bump->center[0]           = 2.;
336       bump->center[1]           = 0.;
337       PetscCall(PetscOptionsRealArray("-advect_bump_center", "location of center of bump x,y", "", bump->center, &two, NULL));
338       bump->radius = 0.9;
339       PetscCall(PetscOptionsReal("-advect_bump_radius", "radius of bump", "", bump->radius, &bump->radius, NULL));
340       bump->type = ADVECT_SOL_BUMP_CONE;
341       PetscCall(PetscOptionsEnum("-advect_bump_type", "type of bump", "", AdvectSolBumpTypes, (PetscEnum)bump->type, (PetscEnum *)&bump->type, NULL));
342       phys->maxspeed = 3.; /* radius of mesh, kludge */
343     } break;
344     }
345   }
346   PetscOptionsHeadEnd();
347   /* Initial/transient solution with default boundary conditions */
348   PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_Advect, phys));
349   /* Register "canned" functionals */
350   PetscCall(ModelFunctionalRegister(mod, "Solution", &advect->functional.Solution, PhysicsFunctional_Advect, phys));
351   PetscCall(ModelFunctionalRegister(mod, "Error", &advect->functional.Error, PhysicsFunctional_Advect, phys));
352   PetscFunctionReturn(0);
353 }
354 
355 /******************* Shallow Water ********************/
356 typedef struct {
357   PetscReal gravity;
358   PetscReal boundaryHeight;
359   struct {
360     PetscInt Height;
361     PetscInt Speed;
362     PetscInt Energy;
363   } functional;
364 } Physics_SW;
365 typedef struct {
366   PetscReal h;
367   PetscReal uh[DIM];
368 } SWNode;
369 typedef union
370 {
371   SWNode    swnode;
372   PetscReal vals[DIM + 1];
373 } SWNodeUnion;
374 
375 static const struct FieldDescription PhysicsFields_SW[] = {
376   {"Height",   1  },
377   {"Momentum", DIM},
378   {NULL,       0  }
379 };
380 
381 /*
382  * h_t + div(uh) = 0
383  * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0
384  *
385  * */
386 static PetscErrorCode SWFlux(Physics phys, const PetscReal *n, const SWNode *x, SWNode *f) {
387   Physics_SW *sw = (Physics_SW *)phys->data;
388   PetscReal   uhn, u[DIM];
389   PetscInt    i;
390 
391   PetscFunctionBeginUser;
392   Scale2Real(1. / x->h, x->uh, u);
393   uhn  = x->uh[0] * n[0] + x->uh[1] * n[1];
394   f->h = uhn;
395   for (i = 0; i < DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
396   PetscFunctionReturn(0);
397 }
398 
399 static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) {
400   PetscFunctionBeginUser;
401   xG[0] = xI[0];
402   xG[1] = -xI[1];
403   xG[2] = -xI[2];
404   PetscFunctionReturn(0);
405 }
406 
407 static void PhysicsRiemann_SW_HLL(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) {
408   Physics_SW *sw = (Physics_SW *)phys->data;
409   PetscReal   aL, aR;
410   PetscReal   nn[DIM];
411 #if !defined(PETSC_USE_COMPLEX)
412   const SWNode *uL = (const SWNode *)xL, *uR = (const SWNode *)xR;
413 #else
414   SWNodeUnion   uLreal, uRreal;
415   const SWNode *uL = &uLreal.swnode;
416   const SWNode *uR = &uRreal.swnode;
417 #endif
418   SWNodeUnion fL, fR;
419   PetscInt    i;
420   PetscReal   zero = 0.;
421 
422 #if defined(PETSC_USE_COMPLEX)
423   uLreal.swnode.h = 0;
424   uRreal.swnode.h = 0;
425   for (i = 0; i < 1 + dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
426   for (i = 0; i < 1 + dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
427 #endif
428   if (uL->h <= 0 || uR->h <= 0) {
429     for (i = 0; i < 1 + dim; i++) flux[i] = zero;
430     return;
431   } /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */
432   nn[0] = n[0];
433   nn[1] = n[1];
434   Normalize2Real(nn);
435   SWFlux(phys, nn, uL, &(fL.swnode));
436   SWFlux(phys, nn, uR, &(fR.swnode));
437   /* gravity wave speed */
438   aL = PetscSqrtReal(sw->gravity * uL->h);
439   aR = PetscSqrtReal(sw->gravity * uR->h);
440   // Defining u_tilda and v_tilda as u and v
441   PetscReal u_L, u_R;
442   u_L = Dot2Real(uL->uh, nn) / uL->h;
443   u_R = Dot2Real(uR->uh, nn) / uR->h;
444   PetscReal sL, sR;
445   sL = PetscMin(u_L - aL, u_R - aR);
446   sR = PetscMax(u_L + aL, u_R + aR);
447   if (sL > zero) {
448     for (i = 0; i < dim + 1; i++) flux[i] = fL.vals[i] * Norm2Real(n);
449   } else if (sR < zero) {
450     for (i = 0; i < dim + 1; i++) flux[i] = fR.vals[i] * Norm2Real(n);
451   } else {
452     for (i = 0; i < dim + 1; i++) flux[i] = ((sR * fL.vals[i] - sL * fR.vals[i] + sR * sL * (xR[i] - xL[i])) / (sR - sL)) * Norm2Real(n);
453   }
454 }
455 
456 static void PhysicsRiemann_SW_Rusanov(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) {
457   Physics_SW *sw = (Physics_SW *)phys->data;
458   PetscReal   cL, cR, speed;
459   PetscReal   nn[DIM];
460 #if !defined(PETSC_USE_COMPLEX)
461   const SWNode *uL = (const SWNode *)xL, *uR = (const SWNode *)xR;
462 #else
463   SWNodeUnion   uLreal, uRreal;
464   const SWNode *uL = &uLreal.swnode;
465   const SWNode *uR = &uRreal.swnode;
466 #endif
467   SWNodeUnion fL, fR;
468   PetscInt    i;
469   PetscReal   zero = 0.;
470 
471 #if defined(PETSC_USE_COMPLEX)
472   uLreal.swnode.h = 0;
473   uRreal.swnode.h = 0;
474   for (i = 0; i < 1 + dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
475   for (i = 0; i < 1 + dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
476 #endif
477   if (uL->h < 0 || uR->h < 0) {
478     for (i = 0; i < 1 + dim; i++) flux[i] = zero / zero;
479     return;
480   } /* reconstructed thickness is negative */
481   nn[0] = n[0];
482   nn[1] = n[1];
483   Normalize2Real(nn);
484   SWFlux(phys, nn, uL, &(fL.swnode));
485   SWFlux(phys, nn, uR, &(fR.swnode));
486   cL    = PetscSqrtReal(sw->gravity * uL->h);
487   cR    = PetscSqrtReal(sw->gravity * uR->h); /* gravity wave speed */
488   speed = PetscMax(PetscAbsReal(Dot2Real(uL->uh, nn) / uL->h) + cL, PetscAbsReal(Dot2Real(uR->uh, nn) / uR->h) + cR);
489   for (i = 0; i < 1 + dim; i++) flux[i] = (0.5 * (fL.vals[i] + fR.vals[i]) + 0.5 * speed * (xL[i] - xR[i])) * Norm2Real(n);
490 }
491 
492 static PetscErrorCode PhysicsSolution_SW(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) {
493   PetscReal dx[2], r, sigma;
494 
495   PetscFunctionBeginUser;
496   PetscCheck(time == 0.0, mod->comm, PETSC_ERR_SUP, "No solution known for time %g", (double)time);
497   dx[0] = x[0] - 1.5;
498   dx[1] = x[1] - 1.0;
499   r     = Norm2Real(dx);
500   sigma = 0.5;
501   u[0]  = 1 + 2 * PetscExpReal(-PetscSqr(r) / (2 * PetscSqr(sigma)));
502   u[1]  = 0.0;
503   u[2]  = 0.0;
504   PetscFunctionReturn(0);
505 }
506 
507 static PetscErrorCode PhysicsFunctional_SW(Model mod, PetscReal time, const PetscReal *coord, const PetscScalar *xx, PetscReal *f, void *ctx) {
508   Physics       phys = (Physics)ctx;
509   Physics_SW   *sw   = (Physics_SW *)phys->data;
510   const SWNode *x    = (const SWNode *)xx;
511   PetscReal     u[2];
512   PetscReal     h;
513 
514   PetscFunctionBeginUser;
515   h = x->h;
516   Scale2Real(1. / x->h, x->uh, u);
517   f[sw->functional.Height] = h;
518   f[sw->functional.Speed]  = Norm2Real(u) + PetscSqrtReal(sw->gravity * h);
519   f[sw->functional.Energy] = 0.5 * (Dot2Real(x->uh, u) + sw->gravity * PetscSqr(h));
520   PetscFunctionReturn(0);
521 }
522 
523 static PetscErrorCode SetUpBC_SW(DM dm, PetscDS prob, Physics phys) {
524   const PetscInt wallids[] = {100, 101, 200, 300};
525   DMLabel        label;
526 
527   PetscFunctionBeginUser;
528   PetscCall(DMGetLabel(dm, "Face Sets", &label));
529   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));
530   PetscFunctionReturn(0);
531 }
532 
533 static PetscErrorCode PhysicsCreate_SW(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject) {
534   Physics_SW *sw;
535   char        sw_riemann[64] = "rusanov";
536 
537   PetscFunctionBeginUser;
538   phys->field_desc = PhysicsFields_SW;
539   PetscCall(PetscNew(&sw));
540   phys->data   = sw;
541   mod->setupbc = SetUpBC_SW;
542 
543   PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov", PhysicsRiemann_SW_Rusanov);
544   PetscFunctionListAdd(&PhysicsRiemannList_SW, "hll", PhysicsRiemann_SW_HLL);
545 
546   PetscOptionsHeadBegin(PetscOptionsObject, "SW options");
547   {
548     void (*PhysicsRiemann_SW)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics);
549     sw->gravity = 1.0;
550     PetscCall(PetscOptionsReal("-sw_gravity", "Gravitational constant", "", sw->gravity, &sw->gravity, NULL));
551     PetscCall(PetscOptionsFList("-sw_riemann", "Riemann solver", "", PhysicsRiemannList_SW, sw_riemann, sw_riemann, sizeof sw_riemann, NULL));
552     PetscCall(PetscFunctionListFind(PhysicsRiemannList_SW, sw_riemann, &PhysicsRiemann_SW));
553     phys->riemann = (PetscRiemannFunc)PhysicsRiemann_SW;
554   }
555   PetscOptionsHeadEnd();
556   phys->maxspeed = PetscSqrtReal(2.0 * sw->gravity); /* Mach 1 for depth of 2 */
557 
558   PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_SW, phys));
559   PetscCall(ModelFunctionalRegister(mod, "Height", &sw->functional.Height, PhysicsFunctional_SW, phys));
560   PetscCall(ModelFunctionalRegister(mod, "Speed", &sw->functional.Speed, PhysicsFunctional_SW, phys));
561   PetscCall(ModelFunctionalRegister(mod, "Energy", &sw->functional.Energy, PhysicsFunctional_SW, phys));
562 
563   PetscFunctionReturn(0);
564 }
565 
566 /******************* Euler Density Shock (EULER_IV_SHOCK,EULER_SS_SHOCK) ********************/
567 /* An initial-value and self-similar solutions of the compressible Euler equations */
568 /* Ravi Samtaney and D. I. Pullin */
569 /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
570 typedef enum {
571   EULER_PAR_GAMMA,
572   EULER_PAR_RHOR,
573   EULER_PAR_AMACH,
574   EULER_PAR_ITANA,
575   EULER_PAR_SIZE
576 } EulerParamIdx;
577 typedef enum {
578   EULER_IV_SHOCK,
579   EULER_SS_SHOCK,
580   EULER_SHOCK_TUBE,
581   EULER_LINEAR_WAVE
582 } EulerType;
583 typedef struct {
584   PetscReal r;
585   PetscReal ru[DIM];
586   PetscReal E;
587 } EulerNode;
588 typedef union
589 {
590   EulerNode eulernode;
591   PetscReal vals[DIM + 2];
592 } EulerNodeUnion;
593 typedef PetscErrorCode (*EquationOfState)(const PetscReal *, const EulerNode *, PetscReal *);
594 typedef struct {
595   EulerType       type;
596   PetscReal       pars[EULER_PAR_SIZE];
597   EquationOfState sound;
598   struct {
599     PetscInt Density;
600     PetscInt Momentum;
601     PetscInt Energy;
602     PetscInt Pressure;
603     PetscInt Speed;
604   } monitor;
605 } Physics_Euler;
606 
607 static const struct FieldDescription PhysicsFields_Euler[] = {
608   {"Density",  1  },
609   {"Momentum", DIM},
610   {"Energy",   1  },
611   {NULL,       0  }
612 };
613 
614 /* initial condition */
615 int                   initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx);
616 static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) {
617   PetscInt       i;
618   Physics        phys = (Physics)ctx;
619   Physics_Euler *eu   = (Physics_Euler *)phys->data;
620   EulerNode     *uu   = (EulerNode *)u;
621   PetscReal      p0, gamma, c;
622   PetscFunctionBeginUser;
623   PetscCheck(time == 0.0, mod->comm, PETSC_ERR_SUP, "No solution known for time %g", (double)time);
624 
625   for (i = 0; i < DIM; i++) uu->ru[i] = 0.0; /* zero out initial velocity */
626   /* set E and rho */
627   gamma = eu->pars[EULER_PAR_GAMMA];
628 
629   if (eu->type == EULER_IV_SHOCK || eu->type == EULER_SS_SHOCK) {
630     /******************* Euler Density Shock ********************/
631     /* On initial-value and self-similar solutions of the compressible Euler equations */
632     /* Ravi Samtaney and D. I. Pullin */
633     /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
634     /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity,  */
635     p0 = 1.;
636     if (x[0] < 0.0 + x[1] * eu->pars[EULER_PAR_ITANA]) {
637       if (x[0] < mod->bounds[0] * 0.5) { /* left of shock (1) */
638         PetscReal amach, rho, press, gas1, p1;
639         amach     = eu->pars[EULER_PAR_AMACH];
640         rho       = 1.;
641         press     = p0;
642         p1        = press * (1.0 + 2.0 * gamma / (gamma + 1.0) * (amach * amach - 1.0));
643         gas1      = (gamma - 1.0) / (gamma + 1.0);
644         uu->r     = rho * (p1 / press + gas1) / (gas1 * p1 / press + 1.0);
645         uu->ru[0] = ((uu->r - rho) * PetscSqrtReal(gamma * press / rho) * amach);
646         uu->E     = p1 / (gamma - 1.0) + .5 / uu->r * uu->ru[0] * uu->ru[0];
647       } else {      /* left of discontinuity (0) */
648         uu->r = 1.; /* rho = 1 */
649         uu->E = p0 / (gamma - 1.0);
650       }
651     } else { /* right of discontinuity (2) */
652       uu->r = eu->pars[EULER_PAR_RHOR];
653       uu->E = p0 / (gamma - 1.0);
654     }
655   } else if (eu->type == EULER_SHOCK_TUBE) {
656     /* 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. */
657     if (x[0] < 0.0) {
658       uu->r = 8.;
659       uu->E = 10. / (gamma - 1.);
660     } else {
661       uu->r = 1.;
662       uu->E = 1. / (gamma - 1.);
663     }
664   } else if (eu->type == EULER_LINEAR_WAVE) {
665     initLinearWave(uu, gamma, x, mod->bounds[1] - mod->bounds[0]);
666   } else SETERRQ(mod->comm, PETSC_ERR_SUP, "Unknown type %d", eu->type);
667 
668   /* set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main; */
669   eu->sound(&gamma, uu, &c);
670   c = (uu->ru[0] / uu->r) + c;
671   if (c > phys->maxspeed) phys->maxspeed = c;
672 
673   PetscFunctionReturn(0);
674 }
675 
676 static PetscErrorCode Pressure_PG(const PetscReal gamma, const EulerNode *x, PetscReal *p) {
677   PetscReal ru2;
678 
679   PetscFunctionBeginUser;
680   ru2  = DotDIMReal(x->ru, x->ru);
681   (*p) = (x->E - 0.5 * ru2 / x->r) * (gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */
682   PetscFunctionReturn(0);
683 }
684 
685 static PetscErrorCode SpeedOfSound_PG(const PetscReal *gamma, const EulerNode *x, PetscReal *c) {
686   PetscReal p;
687 
688   PetscFunctionBeginUser;
689   Pressure_PG(*gamma, x, &p);
690   PetscCheck(p >= 0., PETSC_COMM_WORLD, PETSC_ERR_SUP, "negative pressure time %g -- NEED TO FIX!!!!!!", (double)p);
691   /* pars[EULER_PAR_GAMMA] = heat capacity ratio */
692   (*c) = PetscSqrtReal(*gamma * p / x->r);
693   PetscFunctionReturn(0);
694 }
695 
696 /*
697  * x = (rho,rho*(u_1),...,rho*e)^T
698  * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
699  *
700  * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
701  *
702  */
703 static PetscErrorCode EulerFlux(Physics phys, const PetscReal *n, const EulerNode *x, EulerNode *f) {
704   Physics_Euler *eu = (Physics_Euler *)phys->data;
705   PetscReal      nu, p;
706   PetscInt       i;
707 
708   PetscFunctionBeginUser;
709   Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p);
710   nu   = DotDIMReal(x->ru, n);
711   f->r = nu;                                                     /* A rho u */
712   nu /= x->r;                                                    /* A u */
713   for (i = 0; i < DIM; i++) f->ru[i] = nu * x->ru[i] + n[i] * p; /* r u^2 + p */
714   f->E = nu * (x->E + p);                                        /* u(e+p) */
715   PetscFunctionReturn(0);
716 }
717 
718 /* PetscReal* => EulerNode* conversion */
719 static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, PetscScalar *a_xG, void *ctx) {
720   PetscInt         i;
721   const EulerNode *xI   = (const EulerNode *)a_xI;
722   EulerNode       *xG   = (EulerNode *)a_xG;
723   Physics          phys = (Physics)ctx;
724   Physics_Euler   *eu   = (Physics_Euler *)phys->data;
725   PetscFunctionBeginUser;
726   xG->r = xI->r;                                     /* ghost cell density - same */
727   xG->E = xI->E;                                     /* ghost cell energy - same */
728   if (n[1] != 0.) {                                  /* top and bottom */
729     xG->ru[0] = xI->ru[0];                           /* copy tang to wall */
730     xG->ru[1] = -xI->ru[1];                          /* reflect perp to t/b wall */
731   } else {                                           /* sides */
732     for (i = 0; i < DIM; i++) xG->ru[i] = xI->ru[i]; /* copy */
733   }
734   if (eu->type == EULER_LINEAR_WAVE) { /* debug */
735 #if 0
736     PetscPrintf(PETSC_COMM_WORLD,"%s coord=%g,%g\n",PETSC_FUNCTION_NAME,(double)c[0],(double)c[1]);
737 #endif
738   }
739   PetscFunctionReturn(0);
740 }
741 int         godunovflux(const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, const int *ndim, const PetscReal *gamma);
742 /* PetscReal* => EulerNode* conversion */
743 static void PhysicsRiemann_Euler_Godunov(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) {
744   Physics_Euler *eu = (Physics_Euler *)phys->data;
745   PetscReal      cL, cR, speed, velL, velR, nn[DIM], s2;
746   PetscInt       i;
747   PetscErrorCode ierr;
748 
749   PetscFunctionBeginUser;
750   for (i = 0, s2 = 0.; i < DIM; i++) {
751     nn[i] = n[i];
752     s2 += nn[i] * nn[i];
753   }
754   s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */
755   for (i = 0.; i < DIM; i++) nn[i] /= s2;
756   if (0) { /* Rusanov */
757     const EulerNode *uL = (const EulerNode *)xL, *uR = (const EulerNode *)xR;
758     EulerNodeUnion   fL, fR;
759     EulerFlux(phys, nn, uL, &(fL.eulernode));
760     EulerFlux(phys, nn, uR, &(fR.eulernode));
761     ierr = eu->sound(&eu->pars[EULER_PAR_GAMMA], uL, &cL);
762     if (ierr) exit(13);
763     ierr = eu->sound(&eu->pars[EULER_PAR_GAMMA], uR, &cR);
764     if (ierr) exit(14);
765     velL  = DotDIMReal(uL->ru, nn) / uL->r;
766     velR  = DotDIMReal(uR->ru, nn) / uR->r;
767     speed = PetscMax(velR + cR, velL + cL);
768     for (i = 0; i < 2 + dim; i++) flux[i] = 0.5 * ((fL.vals[i] + fR.vals[i]) + speed * (xL[i] - xR[i])) * s2;
769   } else {
770     int dim = DIM;
771     /* int iwave =  */
772     godunovflux(xL, xR, flux, nn, &dim, &eu->pars[EULER_PAR_GAMMA]);
773     for (i = 0; i < 2 + dim; i++) flux[i] *= s2;
774   }
775   PetscFunctionReturnVoid();
776 }
777 
778 static PetscErrorCode PhysicsFunctional_Euler(Model mod, PetscReal time, const PetscReal *coord, const PetscScalar *xx, PetscReal *f, void *ctx) {
779   Physics          phys = (Physics)ctx;
780   Physics_Euler   *eu   = (Physics_Euler *)phys->data;
781   const EulerNode *x    = (const EulerNode *)xx;
782   PetscReal        p;
783 
784   PetscFunctionBeginUser;
785   f[eu->monitor.Density]  = x->r;
786   f[eu->monitor.Momentum] = NormDIM(x->ru);
787   f[eu->monitor.Energy]   = x->E;
788   f[eu->monitor.Speed]    = NormDIM(x->ru) / x->r;
789   Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p);
790   f[eu->monitor.Pressure] = p;
791   PetscFunctionReturn(0);
792 }
793 
794 static PetscErrorCode SetUpBC_Euler(DM dm, PetscDS prob, Physics phys) {
795   Physics_Euler *eu = (Physics_Euler *)phys->data;
796   DMLabel        label;
797 
798   PetscFunctionBeginUser;
799   PetscCall(DMGetLabel(dm, "Face Sets", &label));
800   if (eu->type == EULER_LINEAR_WAVE) {
801     const PetscInt wallids[] = {100, 101};
802     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));
803   } else {
804     const PetscInt wallids[] = {100, 101, 200, 300};
805     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));
806   }
807   PetscFunctionReturn(0);
808 }
809 
810 static PetscErrorCode PhysicsCreate_Euler(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject) {
811   Physics_Euler *eu;
812 
813   PetscFunctionBeginUser;
814   phys->field_desc = PhysicsFields_Euler;
815   phys->riemann    = (PetscRiemannFunc)PhysicsRiemann_Euler_Godunov;
816   PetscCall(PetscNew(&eu));
817   phys->data   = eu;
818   mod->setupbc = SetUpBC_Euler;
819   PetscOptionsHeadBegin(PetscOptionsObject, "Euler options");
820   {
821     PetscReal alpha;
822     char      type[64] = "linear_wave";
823     PetscBool is;
824     eu->pars[EULER_PAR_GAMMA] = 1.4;
825     eu->pars[EULER_PAR_AMACH] = 2.02;
826     eu->pars[EULER_PAR_RHOR]  = 3.0;
827     eu->pars[EULER_PAR_ITANA] = 0.57735026918963; /* angle of Euler self similar (SS) shock */
828     PetscCall(PetscOptionsReal("-eu_gamma", "Heat capacity ratio", "", eu->pars[EULER_PAR_GAMMA], &eu->pars[EULER_PAR_GAMMA], NULL));
829     PetscCall(PetscOptionsReal("-eu_amach", "Shock speed (Mach)", "", eu->pars[EULER_PAR_AMACH], &eu->pars[EULER_PAR_AMACH], NULL));
830     PetscCall(PetscOptionsReal("-eu_rho2", "Density right of discontinuity", "", eu->pars[EULER_PAR_RHOR], &eu->pars[EULER_PAR_RHOR], NULL));
831     alpha = 60.;
832     PetscCall(PetscOptionsReal("-eu_alpha", "Angle of discontinuity", "", alpha, &alpha, NULL));
833     PetscCheck(alpha > 0. && alpha <= 90., PETSC_COMM_WORLD, PETSC_ERR_SUP, "Alpha bust be > 0 and <= 90 (%g)", (double)alpha);
834     eu->pars[EULER_PAR_ITANA] = 1. / PetscTanReal(alpha * PETSC_PI / 180.0);
835     PetscCall(PetscOptionsString("-eu_type", "Type of Euler test", "", type, type, sizeof(type), NULL));
836     PetscCall(PetscStrcmp(type, "linear_wave", &is));
837     if (is) {
838       /* Remember this should be periodic */
839       eu->type = EULER_LINEAR_WAVE;
840       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "linear_wave"));
841     } else {
842       PetscCheck(DIM == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "DIM must be 2 unless linear wave test %s", type);
843       PetscCall(PetscStrcmp(type, "iv_shock", &is));
844       if (is) {
845         eu->type = EULER_IV_SHOCK;
846         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "iv_shock"));
847       } else {
848         PetscCall(PetscStrcmp(type, "ss_shock", &is));
849         if (is) {
850           eu->type = EULER_SS_SHOCK;
851           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "ss_shock"));
852         } else {
853           PetscCall(PetscStrcmp(type, "shock_tube", &is));
854           if (is) eu->type = EULER_SHOCK_TUBE;
855           else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unknown Euler type %s", type);
856           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "shock_tube"));
857         }
858       }
859     }
860   }
861   PetscOptionsHeadEnd();
862   eu->sound      = SpeedOfSound_PG;
863   phys->maxspeed = 0.; /* will get set in solution */
864   PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_Euler, phys));
865   PetscCall(ModelFunctionalRegister(mod, "Speed", &eu->monitor.Speed, PhysicsFunctional_Euler, phys));
866   PetscCall(ModelFunctionalRegister(mod, "Energy", &eu->monitor.Energy, PhysicsFunctional_Euler, phys));
867   PetscCall(ModelFunctionalRegister(mod, "Density", &eu->monitor.Density, PhysicsFunctional_Euler, phys));
868   PetscCall(ModelFunctionalRegister(mod, "Momentum", &eu->monitor.Momentum, PhysicsFunctional_Euler, phys));
869   PetscCall(ModelFunctionalRegister(mod, "Pressure", &eu->monitor.Pressure, PhysicsFunctional_Euler, phys));
870 
871   PetscFunctionReturn(0);
872 }
873 
874 static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx) {
875   PetscReal err = 0.;
876   PetscInt  i, j;
877 
878   PetscFunctionBeginUser;
879   for (i = 0; i < numComps; i++) {
880     for (j = 0; j < dim; j++) err += PetscSqr(PetscRealPart(grad[i * dim + j]));
881   }
882   *error = volume * err;
883   PetscFunctionReturn(0);
884 }
885 
886 PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition) {
887   PetscSF      sfPoint;
888   PetscSection coordSection;
889   Vec          coordinates;
890   PetscSection sectionCell;
891   PetscScalar *part;
892   PetscInt     cStart, cEnd, c;
893   PetscMPIInt  rank;
894 
895   PetscFunctionBeginUser;
896   PetscCall(DMGetCoordinateSection(dm, &coordSection));
897   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
898   PetscCall(DMClone(dm, dmCell));
899   PetscCall(DMGetPointSF(dm, &sfPoint));
900   PetscCall(DMSetPointSF(*dmCell, sfPoint));
901   PetscCall(DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection));
902   PetscCall(DMSetCoordinatesLocal(*dmCell, coordinates));
903   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
904   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionCell));
905   PetscCall(DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd));
906   PetscCall(PetscSectionSetChart(sectionCell, cStart, cEnd));
907   for (c = cStart; c < cEnd; ++c) PetscCall(PetscSectionSetDof(sectionCell, c, 1));
908   PetscCall(PetscSectionSetUp(sectionCell));
909   PetscCall(DMSetLocalSection(*dmCell, sectionCell));
910   PetscCall(PetscSectionDestroy(&sectionCell));
911   PetscCall(DMCreateLocalVector(*dmCell, partition));
912   PetscCall(PetscObjectSetName((PetscObject)*partition, "partition"));
913   PetscCall(VecGetArray(*partition, &part));
914   for (c = cStart; c < cEnd; ++c) {
915     PetscScalar *p;
916 
917     PetscCall(DMPlexPointLocalRef(*dmCell, c, part, &p));
918     p[0] = rank;
919   }
920   PetscCall(VecRestoreArray(*partition, &part));
921   PetscFunctionReturn(0);
922 }
923 
924 PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user) {
925   DM                 plex, dmMass, dmFace, dmCell, dmCoord;
926   PetscSection       coordSection;
927   Vec                coordinates, facegeom, cellgeom;
928   PetscSection       sectionMass;
929   PetscScalar       *m;
930   const PetscScalar *fgeom, *cgeom, *coords;
931   PetscInt           vStart, vEnd, v;
932 
933   PetscFunctionBeginUser;
934   PetscCall(DMConvert(dm, DMPLEX, &plex));
935   PetscCall(DMGetCoordinateSection(dm, &coordSection));
936   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
937   PetscCall(DMClone(dm, &dmMass));
938   PetscCall(DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection));
939   PetscCall(DMSetCoordinatesLocal(dmMass, coordinates));
940   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionMass));
941   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
942   PetscCall(PetscSectionSetChart(sectionMass, vStart, vEnd));
943   for (v = vStart; v < vEnd; ++v) {
944     PetscInt numFaces;
945 
946     PetscCall(DMPlexGetSupportSize(dmMass, v, &numFaces));
947     PetscCall(PetscSectionSetDof(sectionMass, v, numFaces * numFaces));
948   }
949   PetscCall(PetscSectionSetUp(sectionMass));
950   PetscCall(DMSetLocalSection(dmMass, sectionMass));
951   PetscCall(PetscSectionDestroy(&sectionMass));
952   PetscCall(DMGetLocalVector(dmMass, massMatrix));
953   PetscCall(VecGetArray(*massMatrix, &m));
954   PetscCall(DMPlexGetGeometryFVM(plex, &facegeom, &cellgeom, NULL));
955   PetscCall(VecGetDM(facegeom, &dmFace));
956   PetscCall(VecGetArrayRead(facegeom, &fgeom));
957   PetscCall(VecGetDM(cellgeom, &dmCell));
958   PetscCall(VecGetArrayRead(cellgeom, &cgeom));
959   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
960   PetscCall(VecGetArrayRead(coordinates, &coords));
961   for (v = vStart; v < vEnd; ++v) {
962     const PetscInt  *faces;
963     PetscFVFaceGeom *fgA, *fgB, *cg;
964     PetscScalar     *vertex;
965     PetscInt         numFaces, sides[2], f, g;
966 
967     PetscCall(DMPlexPointLocalRead(dmCoord, v, coords, &vertex));
968     PetscCall(DMPlexGetSupportSize(dmMass, v, &numFaces));
969     PetscCall(DMPlexGetSupport(dmMass, v, &faces));
970     for (f = 0; f < numFaces; ++f) {
971       sides[0] = faces[f];
972       PetscCall(DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA));
973       for (g = 0; g < numFaces; ++g) {
974         const PetscInt *cells = NULL;
975         PetscReal       area  = 0.0;
976         PetscInt        numCells;
977 
978         sides[1] = faces[g];
979         PetscCall(DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB));
980         PetscCall(DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells));
981         PetscCheck(numCells == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces");
982         PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg));
983         area += PetscAbsScalar((vertex[0] - cg->centroid[0]) * (fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1]) * (fgA->centroid[0] - cg->centroid[0]));
984         area += PetscAbsScalar((vertex[0] - cg->centroid[0]) * (fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1]) * (fgB->centroid[0] - cg->centroid[0]));
985         m[f * numFaces + g] = Dot2Real(fgA->normal, fgB->normal) * area * 0.5;
986         PetscCall(DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells));
987       }
988     }
989   }
990   PetscCall(VecRestoreArrayRead(facegeom, &fgeom));
991   PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
992   PetscCall(VecRestoreArrayRead(coordinates, &coords));
993   PetscCall(VecRestoreArray(*massMatrix, &m));
994   PetscCall(DMDestroy(&dmMass));
995   PetscCall(DMDestroy(&plex));
996   PetscFunctionReturn(0);
997 }
998 
999 /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1000 static PetscErrorCode ModelSolutionSetDefault(Model mod, SolutionFunction func, void *ctx) {
1001   PetscFunctionBeginUser;
1002   mod->solution    = func;
1003   mod->solutionctx = ctx;
1004   PetscFunctionReturn(0);
1005 }
1006 
1007 static PetscErrorCode ModelFunctionalRegister(Model mod, const char *name, PetscInt *offset, FunctionalFunction func, void *ctx) {
1008   FunctionalLink link, *ptr;
1009   PetscInt       lastoffset = -1;
1010 
1011   PetscFunctionBeginUser;
1012   for (ptr = &mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
1013   PetscCall(PetscNew(&link));
1014   PetscCall(PetscStrallocpy(name, &link->name));
1015   link->offset = lastoffset + 1;
1016   link->func   = func;
1017   link->ctx    = ctx;
1018   link->next   = NULL;
1019   *ptr         = link;
1020   *offset      = link->offset;
1021   PetscFunctionReturn(0);
1022 }
1023 
1024 static PetscErrorCode ModelFunctionalSetFromOptions(Model mod, PetscOptionItems *PetscOptionsObject) {
1025   PetscInt       i, j;
1026   FunctionalLink link;
1027   char          *names[256];
1028 
1029   PetscFunctionBeginUser;
1030   mod->numMonitored = PETSC_STATIC_ARRAY_LENGTH(names);
1031   PetscCall(PetscOptionsStringArray("-monitor", "list of functionals to monitor", "", names, &mod->numMonitored, NULL));
1032   /* Create list of functionals that will be computed somehow */
1033   PetscCall(PetscMalloc1(mod->numMonitored, &mod->functionalMonitored));
1034   /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
1035   PetscCall(PetscMalloc1(mod->numMonitored, &mod->functionalCall));
1036   mod->numCall = 0;
1037   for (i = 0; i < mod->numMonitored; i++) {
1038     for (link = mod->functionalRegistry; link; link = link->next) {
1039       PetscBool match;
1040       PetscCall(PetscStrcasecmp(names[i], link->name, &match));
1041       if (match) break;
1042     }
1043     PetscCheck(link, mod->comm, PETSC_ERR_USER, "No known functional '%s'", names[i]);
1044     mod->functionalMonitored[i] = link;
1045     for (j = 0; j < i; j++) {
1046       if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1047     }
1048     mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1049   next_name:
1050     PetscCall(PetscFree(names[i]));
1051   }
1052 
1053   /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1054   mod->maxComputed = -1;
1055   for (link = mod->functionalRegistry; link; link = link->next) {
1056     for (i = 0; i < mod->numCall; i++) {
1057       FunctionalLink call = mod->functionalCall[i];
1058       if (link->func == call->func && link->ctx == call->ctx) mod->maxComputed = PetscMax(mod->maxComputed, link->offset);
1059     }
1060   }
1061   PetscFunctionReturn(0);
1062 }
1063 
1064 static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link) {
1065   FunctionalLink l, next;
1066 
1067   PetscFunctionBeginUser;
1068   if (!link) PetscFunctionReturn(0);
1069   l     = *link;
1070   *link = NULL;
1071   for (; l; l = next) {
1072     next = l->next;
1073     PetscCall(PetscFree(l->name));
1074     PetscCall(PetscFree(l));
1075   }
1076   PetscFunctionReturn(0);
1077 }
1078 
1079 /* put the solution callback into a functional callback */
1080 static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx) {
1081   Model mod;
1082   PetscFunctionBeginUser;
1083   mod = (Model)modctx;
1084   PetscCall((*mod->solution)(mod, time, x, u, mod->solutionctx));
1085   PetscFunctionReturn(0);
1086 }
1087 
1088 PetscErrorCode SetInitialCondition(DM dm, Vec X, User user) {
1089   PetscErrorCode (*func[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1090   void *ctx[1];
1091   Model mod = user->model;
1092 
1093   PetscFunctionBeginUser;
1094   func[0] = SolutionFunctional;
1095   ctx[0]  = (void *)mod;
1096   PetscCall(DMProjectFunction(dm, 0.0, func, ctx, INSERT_ALL_VALUES, X));
1097   PetscFunctionReturn(0);
1098 }
1099 
1100 static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer) {
1101   PetscFunctionBeginUser;
1102   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer));
1103   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERVTK));
1104   PetscCall(PetscViewerFileSetName(*viewer, filename));
1105   PetscFunctionReturn(0);
1106 }
1107 
1108 static PetscErrorCode MonitorVTK(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx) {
1109   User        user = (User)ctx;
1110   DM          dm, plex;
1111   PetscViewer viewer;
1112   char        filename[PETSC_MAX_PATH_LEN], *ftable = NULL;
1113   PetscReal   xnorm;
1114 
1115   PetscFunctionBeginUser;
1116   PetscCall(PetscObjectSetName((PetscObject)X, "u"));
1117   PetscCall(VecGetDM(X, &dm));
1118   PetscCall(VecNorm(X, NORM_INFINITY, &xnorm));
1119 
1120   if (stepnum >= 0) stepnum += user->monitorStepOffset;
1121   if (stepnum >= 0) { /* No summary for final time */
1122     Model              mod = user->model;
1123     Vec                cellgeom;
1124     PetscInt           c, cStart, cEnd, fcount, i;
1125     size_t             ftableused, ftablealloc;
1126     const PetscScalar *cgeom, *x;
1127     DM                 dmCell;
1128     DMLabel            vtkLabel;
1129     PetscReal         *fmin, *fmax, *fintegral, *ftmp;
1130 
1131     PetscCall(DMConvert(dm, DMPLEX, &plex));
1132     PetscCall(DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL));
1133     fcount = mod->maxComputed + 1;
1134     PetscCall(PetscMalloc4(fcount, &fmin, fcount, &fmax, fcount, &fintegral, fcount, &ftmp));
1135     for (i = 0; i < fcount; i++) {
1136       fmin[i]      = PETSC_MAX_REAL;
1137       fmax[i]      = PETSC_MIN_REAL;
1138       fintegral[i] = 0;
1139     }
1140     PetscCall(VecGetDM(cellgeom, &dmCell));
1141     PetscCall(DMPlexGetSimplexOrBoxCells(dmCell, 0, &cStart, &cEnd));
1142     PetscCall(VecGetArrayRead(cellgeom, &cgeom));
1143     PetscCall(VecGetArrayRead(X, &x));
1144     PetscCall(DMGetLabel(dm, "vtk", &vtkLabel));
1145     for (c = cStart; c < cEnd; ++c) {
1146       PetscFVCellGeom   *cg;
1147       const PetscScalar *cx     = NULL;
1148       PetscInt           vtkVal = 0;
1149 
1150       /* not that these two routines as currently implemented work for any dm with a
1151        * localSection/globalSection */
1152       PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg));
1153       PetscCall(DMPlexPointGlobalRead(dm, c, x, &cx));
1154       if (vtkLabel) PetscCall(DMLabelGetValue(vtkLabel, c, &vtkVal));
1155       if (!vtkVal || !cx) continue; /* ghost, or not a global cell */
1156       for (i = 0; i < mod->numCall; i++) {
1157         FunctionalLink flink = mod->functionalCall[i];
1158         PetscCall((*flink->func)(mod, time, cg->centroid, cx, ftmp, flink->ctx));
1159       }
1160       for (i = 0; i < fcount; i++) {
1161         fmin[i] = PetscMin(fmin[i], ftmp[i]);
1162         fmax[i] = PetscMax(fmax[i], ftmp[i]);
1163         fintegral[i] += cg->volume * ftmp[i];
1164       }
1165     }
1166     PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
1167     PetscCall(VecRestoreArrayRead(X, &x));
1168     PetscCall(DMDestroy(&plex));
1169     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts)));
1170     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts)));
1171     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fintegral, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts)));
1172 
1173     ftablealloc = fcount * 100;
1174     ftableused  = 0;
1175     PetscCall(PetscMalloc1(ftablealloc, &ftable));
1176     for (i = 0; i < mod->numMonitored; i++) {
1177       size_t         countused;
1178       char           buffer[256], *p;
1179       FunctionalLink flink = mod->functionalMonitored[i];
1180       PetscInt       id    = flink->offset;
1181       if (i % 3) {
1182         PetscCall(PetscArraycpy(buffer, "  ", 2));
1183         p = buffer + 2;
1184       } else if (i) {
1185         char newline[] = "\n";
1186         PetscCall(PetscMemcpy(buffer, newline, sizeof(newline) - 1));
1187         p = buffer + sizeof(newline) - 1;
1188       } else {
1189         p = buffer;
1190       }
1191       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]));
1192       countused--;
1193       countused += p - buffer;
1194       if (countused > ftablealloc - ftableused - 1) { /* reallocate */
1195         char *ftablenew;
1196         ftablealloc = 2 * ftablealloc + countused;
1197         PetscCall(PetscMalloc(ftablealloc, &ftablenew));
1198         PetscCall(PetscArraycpy(ftablenew, ftable, ftableused));
1199         PetscCall(PetscFree(ftable));
1200         ftable = ftablenew;
1201       }
1202       PetscCall(PetscArraycpy(ftable + ftableused, buffer, countused));
1203       ftableused += countused;
1204       ftable[ftableused] = 0;
1205     }
1206     PetscCall(PetscFree4(fmin, fmax, fintegral, ftmp));
1207 
1208     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "% 3" PetscInt_FMT "  time %8.4g  |x| %8.4g  %s\n", stepnum, (double)time, (double)xnorm, ftable ? ftable : ""));
1209     PetscCall(PetscFree(ftable));
1210   }
1211   if (user->vtkInterval < 1) PetscFunctionReturn(0);
1212   if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1213     if (stepnum == -1) { /* Final time is not multiple of normal time interval, write it anyway */
1214       PetscCall(TSGetStepNumber(ts, &stepnum));
1215     }
1216     PetscCall(PetscSNPrintf(filename, sizeof filename, "%s-%03" PetscInt_FMT ".vtu", user->outputBasename, stepnum));
1217     PetscCall(OutputVTK(dm, filename, &viewer));
1218     PetscCall(VecView(X, viewer));
1219     PetscCall(PetscViewerDestroy(&viewer));
1220   }
1221   PetscFunctionReturn(0);
1222 }
1223 
1224 static PetscErrorCode initializeTS(DM dm, User user, TS *ts) {
1225   PetscFunctionBeginUser;
1226   PetscCall(TSCreate(PetscObjectComm((PetscObject)dm), ts));
1227   PetscCall(TSSetType(*ts, TSSSP));
1228   PetscCall(TSSetDM(*ts, dm));
1229   if (user->vtkmon) PetscCall(TSMonitorSet(*ts, MonitorVTK, user, NULL));
1230   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user));
1231   PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user));
1232   PetscCall(TSSetMaxTime(*ts, 2.0));
1233   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
1234   PetscFunctionReturn(0);
1235 }
1236 
1237 static PetscErrorCode adaptToleranceFVM(PetscFV fvm, TS ts, Vec sol, VecTagger refineTag, VecTagger coarsenTag, User user, TS *tsNew, Vec *solNew) {
1238   DM                 dm, gradDM, plex, cellDM, adaptedDM = NULL;
1239   Vec                cellGeom, faceGeom;
1240   PetscBool          isForest, computeGradient;
1241   Vec                grad, locGrad, locX, errVec;
1242   PetscInt           cStart, cEnd, c, dim, nRefine, nCoarsen;
1243   PetscReal          minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2], minInd, maxInd, time;
1244   PetscScalar       *errArray;
1245   const PetscScalar *pointVals;
1246   const PetscScalar *pointGrads;
1247   const PetscScalar *pointGeom;
1248   DMLabel            adaptLabel = NULL;
1249   IS                 refineIS, coarsenIS;
1250 
1251   PetscFunctionBeginUser;
1252   PetscCall(TSGetTime(ts, &time));
1253   PetscCall(VecGetDM(sol, &dm));
1254   PetscCall(DMGetDimension(dm, &dim));
1255   PetscCall(PetscFVGetComputeGradients(fvm, &computeGradient));
1256   PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE));
1257   PetscCall(DMIsForest(dm, &isForest));
1258   PetscCall(DMConvert(dm, DMPLEX, &plex));
1259   PetscCall(DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM));
1260   PetscCall(DMCreateLocalVector(plex, &locX));
1261   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL));
1262   PetscCall(DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX));
1263   PetscCall(DMGlobalToLocalEnd(plex, sol, INSERT_VALUES, locX));
1264   PetscCall(DMCreateGlobalVector(gradDM, &grad));
1265   PetscCall(DMPlexReconstructGradientsFVM(plex, locX, grad));
1266   PetscCall(DMCreateLocalVector(gradDM, &locGrad));
1267   PetscCall(DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad));
1268   PetscCall(DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad));
1269   PetscCall(VecDestroy(&grad));
1270   PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &cStart, &cEnd));
1271   PetscCall(VecGetArrayRead(locGrad, &pointGrads));
1272   PetscCall(VecGetArrayRead(cellGeom, &pointGeom));
1273   PetscCall(VecGetArrayRead(locX, &pointVals));
1274   PetscCall(VecGetDM(cellGeom, &cellDM));
1275   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
1276   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)plex), cEnd - cStart, PETSC_DETERMINE, &errVec));
1277   PetscCall(VecSetUp(errVec));
1278   PetscCall(VecGetArray(errVec, &errArray));
1279   for (c = cStart; c < cEnd; c++) {
1280     PetscReal        errInd = 0.;
1281     PetscScalar     *pointGrad;
1282     PetscScalar     *pointVal;
1283     PetscFVCellGeom *cg;
1284 
1285     PetscCall(DMPlexPointLocalRead(gradDM, c, pointGrads, &pointGrad));
1286     PetscCall(DMPlexPointLocalRead(cellDM, c, pointGeom, &cg));
1287     PetscCall(DMPlexPointLocalRead(plex, c, pointVals, &pointVal));
1288 
1289     PetscCall((user->model->errorIndicator)(dim, cg->volume, user->model->physics->dof, pointVal, pointGrad, &errInd, user->model->errorCtx));
1290     errArray[c - cStart] = errInd;
1291     minMaxInd[0]         = PetscMin(minMaxInd[0], errInd);
1292     minMaxInd[1]         = PetscMax(minMaxInd[1], errInd);
1293   }
1294   PetscCall(VecRestoreArray(errVec, &errArray));
1295   PetscCall(VecRestoreArrayRead(locX, &pointVals));
1296   PetscCall(VecRestoreArrayRead(cellGeom, &pointGeom));
1297   PetscCall(VecRestoreArrayRead(locGrad, &pointGrads));
1298   PetscCall(VecDestroy(&locGrad));
1299   PetscCall(VecDestroy(&locX));
1300   PetscCall(DMDestroy(&plex));
1301 
1302   PetscCall(VecTaggerComputeIS(refineTag, errVec, &refineIS, NULL));
1303   PetscCall(VecTaggerComputeIS(coarsenTag, errVec, &coarsenIS, NULL));
1304   PetscCall(ISGetSize(refineIS, &nRefine));
1305   PetscCall(ISGetSize(coarsenIS, &nCoarsen));
1306   if (nRefine) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS));
1307   if (nCoarsen) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_COARSEN, coarsenIS));
1308   PetscCall(ISDestroy(&coarsenIS));
1309   PetscCall(ISDestroy(&refineIS));
1310   PetscCall(VecDestroy(&errVec));
1311 
1312   PetscCall(PetscFVSetComputeGradients(fvm, computeGradient));
1313   minMaxInd[1] = -minMaxInd[1];
1314   PetscCallMPI(MPI_Allreduce(minMaxInd, minMaxIndGlobal, 2, MPIU_REAL, MPI_MIN, PetscObjectComm((PetscObject)dm)));
1315   minInd = minMaxIndGlobal[0];
1316   maxInd = -minMaxIndGlobal[1];
1317   PetscCall(PetscInfo(ts, "error indicator range (%E, %E)\n", (double)minInd, (double)maxInd));
1318   if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */
1319     PetscCall(DMAdaptLabel(dm, adaptLabel, &adaptedDM));
1320   }
1321   PetscCall(DMLabelDestroy(&adaptLabel));
1322   if (adaptedDM) {
1323     PetscCall(PetscInfo(ts, "Adapted mesh, marking %" PetscInt_FMT " cells for refinement, and %" PetscInt_FMT " cells for coarsening\n", nRefine, nCoarsen));
1324     if (tsNew) PetscCall(initializeTS(adaptedDM, user, tsNew));
1325     if (solNew) {
1326       PetscCall(DMCreateGlobalVector(adaptedDM, solNew));
1327       PetscCall(PetscObjectSetName((PetscObject)*solNew, "solution"));
1328       PetscCall(DMForestTransferVec(dm, sol, adaptedDM, *solNew, PETSC_TRUE, time));
1329     }
1330     if (isForest) PetscCall(DMForestSetAdaptivityForest(adaptedDM, NULL)); /* clear internal references to the previous dm */
1331     PetscCall(DMDestroy(&adaptedDM));
1332   } else {
1333     if (tsNew) *tsNew = NULL;
1334     if (solNew) *solNew = NULL;
1335   }
1336   PetscFunctionReturn(0);
1337 }
1338 
1339 int main(int argc, char **argv) {
1340   MPI_Comm          comm;
1341   PetscDS           prob;
1342   PetscFV           fvm;
1343   PetscLimiter      limiter = NULL, noneLimiter = NULL;
1344   User              user;
1345   Model             mod;
1346   Physics           phys;
1347   DM                dm, plex;
1348   PetscReal         ftime, cfl, dt, minRadius;
1349   PetscInt          dim, nsteps;
1350   TS                ts;
1351   TSConvergedReason reason;
1352   Vec               X;
1353   PetscViewer       viewer;
1354   PetscBool         vtkCellGeom, useAMR;
1355   PetscInt          adaptInterval;
1356   char              physname[256] = "advect";
1357   VecTagger         refineTag = NULL, coarsenTag = NULL;
1358 
1359   PetscFunctionBeginUser;
1360   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1361   comm = PETSC_COMM_WORLD;
1362 
1363   PetscCall(PetscNew(&user));
1364   PetscCall(PetscNew(&user->model));
1365   PetscCall(PetscNew(&user->model->physics));
1366   mod           = user->model;
1367   phys          = mod->physics;
1368   mod->comm     = comm;
1369   useAMR        = PETSC_FALSE;
1370   adaptInterval = 1;
1371 
1372   /* Register physical models to be available on the command line */
1373   PetscCall(PetscFunctionListAdd(&PhysicsList, "advect", PhysicsCreate_Advect));
1374   PetscCall(PetscFunctionListAdd(&PhysicsList, "sw", PhysicsCreate_SW));
1375   PetscCall(PetscFunctionListAdd(&PhysicsList, "euler", PhysicsCreate_Euler));
1376 
1377   PetscOptionsBegin(comm, NULL, "Unstructured Finite Volume Mesh Options", "");
1378   {
1379     cfl = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1380     PetscCall(PetscOptionsReal("-ufv_cfl", "CFL number per step", "", cfl, &cfl, NULL));
1381     user->vtkInterval = 1;
1382     PetscCall(PetscOptionsInt("-ufv_vtk_interval", "VTK output interval (0 to disable)", "", user->vtkInterval, &user->vtkInterval, NULL));
1383     user->vtkmon = PETSC_TRUE;
1384     PetscCall(PetscOptionsBool("-ufv_vtk_monitor", "Use VTKMonitor routine", "", user->vtkmon, &user->vtkmon, NULL));
1385     vtkCellGeom = PETSC_FALSE;
1386     PetscCall(PetscStrcpy(user->outputBasename, "ex11"));
1387     PetscCall(PetscOptionsString("-ufv_vtk_basename", "VTK output basename", "", user->outputBasename, user->outputBasename, sizeof(user->outputBasename), NULL));
1388     PetscCall(PetscOptionsBool("-ufv_vtk_cellgeom", "Write cell geometry (for debugging)", "", vtkCellGeom, &vtkCellGeom, NULL));
1389     PetscCall(PetscOptionsBool("-ufv_use_amr", "use local adaptive mesh refinement", "", useAMR, &useAMR, NULL));
1390     PetscCall(PetscOptionsInt("-ufv_adapt_interval", "time steps between AMR", "", adaptInterval, &adaptInterval, NULL));
1391   }
1392   PetscOptionsEnd();
1393 
1394   if (useAMR) {
1395     VecTaggerBox refineBox, coarsenBox;
1396 
1397     refineBox.min = refineBox.max = PETSC_MAX_REAL;
1398     coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL;
1399 
1400     PetscCall(VecTaggerCreate(comm, &refineTag));
1401     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)refineTag, "refine_"));
1402     PetscCall(VecTaggerSetType(refineTag, VECTAGGERABSOLUTE));
1403     PetscCall(VecTaggerAbsoluteSetBox(refineTag, &refineBox));
1404     PetscCall(VecTaggerSetFromOptions(refineTag));
1405     PetscCall(VecTaggerSetUp(refineTag));
1406     PetscCall(PetscObjectViewFromOptions((PetscObject)refineTag, NULL, "-tag_view"));
1407 
1408     PetscCall(VecTaggerCreate(comm, &coarsenTag));
1409     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)coarsenTag, "coarsen_"));
1410     PetscCall(VecTaggerSetType(coarsenTag, VECTAGGERABSOLUTE));
1411     PetscCall(VecTaggerAbsoluteSetBox(coarsenTag, &coarsenBox));
1412     PetscCall(VecTaggerSetFromOptions(coarsenTag));
1413     PetscCall(VecTaggerSetUp(coarsenTag));
1414     PetscCall(PetscObjectViewFromOptions((PetscObject)coarsenTag, NULL, "-tag_view"));
1415   }
1416 
1417   PetscOptionsBegin(comm, NULL, "Unstructured Finite Volume Physics Options", "");
1418   {
1419     PetscErrorCode (*physcreate)(Model, Physics, PetscOptionItems *);
1420     PetscCall(PetscOptionsFList("-physics", "Physics module to solve", "", PhysicsList, physname, physname, sizeof physname, NULL));
1421     PetscCall(PetscFunctionListFind(PhysicsList, physname, &physcreate));
1422     PetscCall(PetscMemzero(phys, sizeof(struct _n_Physics)));
1423     PetscCall((*physcreate)(mod, phys, PetscOptionsObject));
1424     /* Count number of fields and dofs */
1425     for (phys->nfields = 0, phys->dof = 0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1426     PetscCheck(phys->dof > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set dof", physname);
1427     PetscCall(ModelFunctionalSetFromOptions(mod, PetscOptionsObject));
1428   }
1429   PetscOptionsEnd();
1430 
1431   /* Create mesh */
1432   {
1433     PetscInt i;
1434 
1435     PetscCall(DMCreate(comm, &dm));
1436     PetscCall(DMSetType(dm, DMPLEX));
1437     PetscCall(DMSetFromOptions(dm));
1438     for (i = 0; i < DIM; i++) {
1439       mod->bounds[2 * i]     = 0.;
1440       mod->bounds[2 * i + 1] = 1.;
1441     };
1442     dim = DIM;
1443     { /* a null name means just do a hex box */
1444       PetscInt  cells[3] = {1, 1, 1}, n = 3;
1445       PetscBool flg2, skew              = PETSC_FALSE;
1446       PetscInt  nret2 = 2 * DIM;
1447       PetscOptionsBegin(comm, NULL, "Rectangular mesh options", "");
1448       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));
1449       PetscCall(PetscOptionsBool("-grid_skew_60", "Skew grid for 60 degree shock mesh", "", skew, &skew, NULL));
1450       PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", cells, &n, NULL));
1451       PetscOptionsEnd();
1452       /* TODO Rewrite this with Mark, and remove grid_bounds at that time */
1453       if (flg2) {
1454         PetscInt     dimEmbed, i;
1455         PetscInt     nCoords;
1456         PetscScalar *coords;
1457         Vec          coordinates;
1458 
1459         PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1460         PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
1461         PetscCall(VecGetLocalSize(coordinates, &nCoords));
1462         PetscCheck(!(nCoords % dimEmbed), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Coordinate vector the wrong size");
1463         PetscCall(VecGetArray(coordinates, &coords));
1464         for (i = 0; i < nCoords; i += dimEmbed) {
1465           PetscInt j;
1466 
1467           PetscScalar *coord = &coords[i];
1468           for (j = 0; j < dimEmbed; j++) {
1469             coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]);
1470             if (dim == 2 && cells[1] == 1 && j == 0 && skew) {
1471               if (cells[0] == 2 && i == 8) {
1472                 coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */
1473               } else if (cells[0] == 3) {
1474                 if (i == 2 || i == 10) coord[j] = mod->bounds[1] / 4.;
1475                 else if (i == 4) coord[j] = mod->bounds[1] / 2.;
1476                 else if (i == 12) coord[j] = 1.57735026918963 * mod->bounds[1] / 2.;
1477               }
1478             }
1479           }
1480         }
1481         PetscCall(VecRestoreArray(coordinates, &coords));
1482         PetscCall(DMSetCoordinatesLocal(dm, coordinates));
1483       }
1484     }
1485   }
1486   PetscCall(DMViewFromOptions(dm, NULL, "-orig_dm_view"));
1487   PetscCall(DMGetDimension(dm, &dim));
1488 
1489   /* set up BCs, functions, tags */
1490   PetscCall(DMCreateLabel(dm, "Face Sets"));
1491   mod->errorIndicator = ErrorIndicator_Simple;
1492 
1493   {
1494     DM gdm;
1495 
1496     PetscCall(DMPlexConstructGhostCells(dm, NULL, NULL, &gdm));
1497     PetscCall(DMDestroy(&dm));
1498     dm = gdm;
1499     PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1500   }
1501 
1502   PetscCall(PetscFVCreate(comm, &fvm));
1503   PetscCall(PetscFVSetFromOptions(fvm));
1504   PetscCall(PetscFVSetNumComponents(fvm, phys->dof));
1505   PetscCall(PetscFVSetSpatialDimension(fvm, dim));
1506   PetscCall(PetscObjectSetName((PetscObject)fvm, ""));
1507   {
1508     PetscInt f, dof;
1509     for (f = 0, dof = 0; f < phys->nfields; f++) {
1510       PetscInt newDof = phys->field_desc[f].dof;
1511 
1512       if (newDof == 1) {
1513         PetscCall(PetscFVSetComponentName(fvm, dof, phys->field_desc[f].name));
1514       } else {
1515         PetscInt j;
1516 
1517         for (j = 0; j < newDof; j++) {
1518           char compName[256] = "Unknown";
1519 
1520           PetscCall(PetscSNPrintf(compName, sizeof(compName), "%s_%" PetscInt_FMT, phys->field_desc[f].name, j));
1521           PetscCall(PetscFVSetComponentName(fvm, dof + j, compName));
1522         }
1523       }
1524       dof += newDof;
1525     }
1526   }
1527   /* FV is now structured with one field having all physics as components */
1528   PetscCall(DMAddField(dm, NULL, (PetscObject)fvm));
1529   PetscCall(DMCreateDS(dm));
1530   PetscCall(DMGetDS(dm, &prob));
1531   PetscCall(PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann));
1532   PetscCall(PetscDSSetContext(prob, 0, user->model->physics));
1533   PetscCall((*mod->setupbc)(dm, prob, phys));
1534   PetscCall(PetscDSSetFromOptions(prob));
1535   {
1536     char      convType[256];
1537     PetscBool flg;
1538 
1539     PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
1540     PetscCall(PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex12", DMList, DMPLEX, convType, 256, &flg));
1541     PetscOptionsEnd();
1542     if (flg) {
1543       DM dmConv;
1544 
1545       PetscCall(DMConvert(dm, convType, &dmConv));
1546       if (dmConv) {
1547         PetscCall(DMViewFromOptions(dmConv, NULL, "-dm_conv_view"));
1548         PetscCall(DMDestroy(&dm));
1549         dm = dmConv;
1550         PetscCall(DMSetFromOptions(dm));
1551       }
1552     }
1553   }
1554 
1555   PetscCall(initializeTS(dm, user, &ts));
1556 
1557   PetscCall(DMCreateGlobalVector(dm, &X));
1558   PetscCall(PetscObjectSetName((PetscObject)X, "solution"));
1559   PetscCall(SetInitialCondition(dm, X, user));
1560   if (useAMR) {
1561     PetscInt adaptIter;
1562 
1563     /* use no limiting when reconstructing gradients for adaptivity */
1564     PetscCall(PetscFVGetLimiter(fvm, &limiter));
1565     PetscCall(PetscObjectReference((PetscObject)limiter));
1566     PetscCall(PetscLimiterCreate(PetscObjectComm((PetscObject)fvm), &noneLimiter));
1567     PetscCall(PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE));
1568 
1569     PetscCall(PetscFVSetLimiter(fvm, noneLimiter));
1570     for (adaptIter = 0;; ++adaptIter) {
1571       PetscLogDouble bytes;
1572       TS             tsNew = NULL;
1573 
1574       PetscCall(PetscMemoryGetCurrentUsage(&bytes));
1575       PetscCall(PetscInfo(ts, "refinement loop %" PetscInt_FMT ": memory used %g\n", adaptIter, (double)bytes));
1576       PetscCall(DMViewFromOptions(dm, NULL, "-initial_dm_view"));
1577       PetscCall(VecViewFromOptions(X, NULL, "-initial_vec_view"));
1578 #if 0
1579       if (viewInitial) {
1580         PetscViewer viewer;
1581         char        buf[256];
1582         PetscBool   isHDF5, isVTK;
1583 
1584         PetscCall(PetscViewerCreate(comm,&viewer));
1585         PetscCall(PetscViewerSetType(viewer,PETSCVIEWERVTK));
1586         PetscCall(PetscViewerSetOptionsPrefix(viewer,"initial_"));
1587         PetscCall(PetscViewerSetFromOptions(viewer));
1588         PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5));
1589         PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK));
1590         if (isHDF5) {
1591           PetscCall(PetscSNPrintf(buf, 256, "ex11-initial-%" PetscInt_FMT ".h5", adaptIter));
1592         } else if (isVTK) {
1593           PetscCall(PetscSNPrintf(buf, 256, "ex11-initial-%" PetscInt_FMT ".vtu", adaptIter));
1594           PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU));
1595         }
1596         PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE));
1597         PetscCall(PetscViewerFileSetName(viewer,buf));
1598         if (isHDF5) {
1599           PetscCall(DMView(dm,viewer));
1600           PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_UPDATE));
1601         }
1602         PetscCall(VecView(X,viewer));
1603         PetscCall(PetscViewerDestroy(&viewer));
1604       }
1605 #endif
1606 
1607       PetscCall(adaptToleranceFVM(fvm, ts, X, refineTag, coarsenTag, user, &tsNew, NULL));
1608       if (!tsNew) {
1609         break;
1610       } else {
1611         PetscCall(DMDestroy(&dm));
1612         PetscCall(VecDestroy(&X));
1613         PetscCall(TSDestroy(&ts));
1614         ts = tsNew;
1615         PetscCall(TSGetDM(ts, &dm));
1616         PetscCall(PetscObjectReference((PetscObject)dm));
1617         PetscCall(DMCreateGlobalVector(dm, &X));
1618         PetscCall(PetscObjectSetName((PetscObject)X, "solution"));
1619         PetscCall(SetInitialCondition(dm, X, user));
1620       }
1621     }
1622     /* restore original limiter */
1623     PetscCall(PetscFVSetLimiter(fvm, limiter));
1624   }
1625 
1626   PetscCall(DMConvert(dm, DMPLEX, &plex));
1627   if (vtkCellGeom) {
1628     DM  dmCell;
1629     Vec cellgeom, partition;
1630 
1631     PetscCall(DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL));
1632     PetscCall(OutputVTK(dm, "ex11-cellgeom.vtk", &viewer));
1633     PetscCall(VecView(cellgeom, viewer));
1634     PetscCall(PetscViewerDestroy(&viewer));
1635     PetscCall(CreatePartitionVec(dm, &dmCell, &partition));
1636     PetscCall(OutputVTK(dmCell, "ex11-partition.vtk", &viewer));
1637     PetscCall(VecView(partition, viewer));
1638     PetscCall(PetscViewerDestroy(&viewer));
1639     PetscCall(VecDestroy(&partition));
1640     PetscCall(DMDestroy(&dmCell));
1641   }
1642   /* collect max maxspeed from all processes -- todo */
1643   PetscCall(DMPlexGetGeometryFVM(plex, NULL, NULL, &minRadius));
1644   PetscCall(DMDestroy(&plex));
1645   PetscCallMPI(MPI_Allreduce(&phys->maxspeed, &mod->maxspeed, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts)));
1646   PetscCheck(mod->maxspeed > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set maxspeed", physname);
1647   dt = cfl * minRadius / mod->maxspeed;
1648   PetscCall(TSSetTimeStep(ts, dt));
1649   PetscCall(TSSetFromOptions(ts));
1650   if (!useAMR) {
1651     PetscCall(TSSolve(ts, X));
1652     PetscCall(TSGetSolveTime(ts, &ftime));
1653     PetscCall(TSGetStepNumber(ts, &nsteps));
1654   } else {
1655     PetscReal finalTime;
1656     PetscInt  adaptIter;
1657     TS        tsNew  = NULL;
1658     Vec       solNew = NULL;
1659 
1660     PetscCall(TSGetMaxTime(ts, &finalTime));
1661     PetscCall(TSSetMaxSteps(ts, adaptInterval));
1662     PetscCall(TSSolve(ts, X));
1663     PetscCall(TSGetSolveTime(ts, &ftime));
1664     PetscCall(TSGetStepNumber(ts, &nsteps));
1665     for (adaptIter = 0; ftime < finalTime; adaptIter++) {
1666       PetscLogDouble bytes;
1667 
1668       PetscCall(PetscMemoryGetCurrentUsage(&bytes));
1669       PetscCall(PetscInfo(ts, "AMR time step loop %" PetscInt_FMT ": memory used %g\n", adaptIter, bytes));
1670       PetscCall(PetscFVSetLimiter(fvm, noneLimiter));
1671       PetscCall(adaptToleranceFVM(fvm, ts, X, refineTag, coarsenTag, user, &tsNew, &solNew));
1672       PetscCall(PetscFVSetLimiter(fvm, limiter));
1673       if (tsNew) {
1674         PetscCall(PetscInfo(ts, "AMR used\n"));
1675         PetscCall(DMDestroy(&dm));
1676         PetscCall(VecDestroy(&X));
1677         PetscCall(TSDestroy(&ts));
1678         ts = tsNew;
1679         X  = solNew;
1680         PetscCall(TSSetFromOptions(ts));
1681         PetscCall(VecGetDM(X, &dm));
1682         PetscCall(PetscObjectReference((PetscObject)dm));
1683         PetscCall(DMConvert(dm, DMPLEX, &plex));
1684         PetscCall(DMPlexGetGeometryFVM(dm, NULL, NULL, &minRadius));
1685         PetscCall(DMDestroy(&plex));
1686         PetscCallMPI(MPI_Allreduce(&phys->maxspeed, &mod->maxspeed, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts)));
1687         PetscCheck(mod->maxspeed > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set maxspeed", physname);
1688         dt = cfl * minRadius / mod->maxspeed;
1689         PetscCall(TSSetStepNumber(ts, nsteps));
1690         PetscCall(TSSetTime(ts, ftime));
1691         PetscCall(TSSetTimeStep(ts, dt));
1692       } else {
1693         PetscCall(PetscInfo(ts, "AMR not used\n"));
1694       }
1695       user->monitorStepOffset = nsteps;
1696       PetscCall(TSSetMaxSteps(ts, nsteps + adaptInterval));
1697       PetscCall(TSSolve(ts, X));
1698       PetscCall(TSGetSolveTime(ts, &ftime));
1699       PetscCall(TSGetStepNumber(ts, &nsteps));
1700     }
1701   }
1702   PetscCall(TSGetConvergedReason(ts, &reason));
1703   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, nsteps));
1704   PetscCall(TSDestroy(&ts));
1705 
1706   PetscCall(VecTaggerDestroy(&refineTag));
1707   PetscCall(VecTaggerDestroy(&coarsenTag));
1708   PetscCall(PetscFunctionListDestroy(&PhysicsList));
1709   PetscCall(PetscFunctionListDestroy(&PhysicsRiemannList_SW));
1710   PetscCall(FunctionalLinkDestroy(&user->model->functionalRegistry));
1711   PetscCall(PetscFree(user->model->functionalMonitored));
1712   PetscCall(PetscFree(user->model->functionalCall));
1713   PetscCall(PetscFree(user->model->physics->data));
1714   PetscCall(PetscFree(user->model->physics));
1715   PetscCall(PetscFree(user->model));
1716   PetscCall(PetscFree(user));
1717   PetscCall(VecDestroy(&X));
1718   PetscCall(PetscLimiterDestroy(&limiter));
1719   PetscCall(PetscLimiterDestroy(&noneLimiter));
1720   PetscCall(PetscFVDestroy(&fvm));
1721   PetscCall(DMDestroy(&dm));
1722   PetscCall(PetscFinalize());
1723   return 0;
1724 }
1725 
1726 /* Godunov fluxs */
1727 PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test) {
1728   /* System generated locals */
1729   PetscScalar ret_val;
1730 
1731   if (PetscRealPart(*test) > 0.) goto L10;
1732   ret_val = *b;
1733   return ret_val;
1734 L10:
1735   ret_val = *a;
1736   return ret_val;
1737 } /* cvmgp_ */
1738 
1739 PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test) {
1740   /* System generated locals */
1741   PetscScalar ret_val;
1742 
1743   if (PetscRealPart(*test) < 0.) goto L10;
1744   ret_val = *b;
1745   return ret_val;
1746 L10:
1747   ret_val = *a;
1748   return ret_val;
1749 } /* cvmgm_ */
1750 
1751 int riem1mdt(PetscScalar *gaml, PetscScalar *gamr, PetscScalar *rl, PetscScalar *pl, PetscScalar *uxl, PetscScalar *rr, PetscScalar *pr, PetscScalar *uxr, PetscScalar *rstarl, PetscScalar *rstarr, PetscScalar *pstar, PetscScalar *ustar) {
1752   /* Initialized data */
1753 
1754   static PetscScalar smallp = 1e-8;
1755 
1756   /* System generated locals */
1757   int         i__1;
1758   PetscScalar d__1, d__2;
1759 
1760   /* Local variables */
1761   static int         i0;
1762   static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2;
1763   static int         iwave;
1764   static PetscScalar gascl4, gascr4, cstarl, dpstar, cstarr;
1765   /* static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascr1, gascr2, gascr3; */
1766   static int         iterno;
1767   static PetscScalar ustarl, ustarr, rarepr1, rarepr2;
1768 
1769   /* gascl1 = *gaml - 1.; */
1770   /* gascl2 = (*gaml + 1.) * .5; */
1771   /* gascl3 = gascl2 / *gaml; */
1772   gascl4 = 1. / (*gaml - 1.);
1773 
1774   /* gascr1 = *gamr - 1.; */
1775   /* gascr2 = (*gamr + 1.) * .5; */
1776   /* gascr3 = gascr2 / *gamr; */
1777   gascr4  = 1. / (*gamr - 1.);
1778   iterno  = 10;
1779   /*        find pstar: */
1780   cl      = PetscSqrtScalar(*gaml * *pl / *rl);
1781   cr      = PetscSqrtScalar(*gamr * *pr / *rr);
1782   wl      = *rl * cl;
1783   wr      = *rr * cr;
1784   /* csqrl = wl * wl; */
1785   /* csqrr = wr * wr; */
1786   *pstar  = (wl * *pr + wr * *pl) / (wl + wr);
1787   *pstar  = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
1788   pst     = *pl / *pr;
1789   skpr1   = cr * (pst - 1.) * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
1790   d__1    = (*gamr - 1.) / (*gamr * 2.);
1791   rarepr2 = gascr4 * 2. * cr * (1. - PetscPowScalar(pst, d__1));
1792   pst     = *pr / *pl;
1793   skpr2   = cl * (pst - 1.) * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
1794   d__1    = (*gaml - 1.) / (*gaml * 2.);
1795   rarepr1 = gascl4 * 2. * cl * (1. - PetscPowScalar(pst, d__1));
1796   durl    = *uxr - *uxl;
1797   if (PetscRealPart(*pr) < PetscRealPart(*pl)) {
1798     if (PetscRealPart(durl) >= PetscRealPart(rarepr1)) {
1799       iwave = 100;
1800     } else if (PetscRealPart(durl) <= PetscRealPart(-skpr1)) {
1801       iwave = 300;
1802     } else {
1803       iwave = 400;
1804     }
1805   } else {
1806     if (PetscRealPart(durl) >= PetscRealPart(rarepr2)) {
1807       iwave = 100;
1808     } else if (PetscRealPart(durl) <= PetscRealPart(-skpr2)) {
1809       iwave = 300;
1810     } else {
1811       iwave = 200;
1812     }
1813   }
1814   if (iwave == 100) {
1815     /*     1-wave: rarefaction wave, 3-wave: rarefaction wave */
1816     /*     case (100) */
1817     i__1 = iterno;
1818     for (i0 = 1; i0 <= i__1; ++i0) {
1819       d__1    = *pstar / *pl;
1820       d__2    = 1. / *gaml;
1821       *rstarl = *rl * PetscPowScalar(d__1, d__2);
1822       cstarl  = PetscSqrtScalar(*gaml * *pstar / *rstarl);
1823       ustarl  = *uxl - gascl4 * 2. * (cstarl - cl);
1824       zl      = *rstarl * cstarl;
1825       d__1    = *pstar / *pr;
1826       d__2    = 1. / *gamr;
1827       *rstarr = *rr * PetscPowScalar(d__1, d__2);
1828       cstarr  = PetscSqrtScalar(*gamr * *pstar / *rstarr);
1829       ustarr  = *uxr + gascr4 * 2. * (cstarr - cr);
1830       zr      = *rstarr * cstarr;
1831       dpstar  = zl * zr * (ustarr - ustarl) / (zl + zr);
1832       *pstar -= dpstar;
1833       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
1834       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1835 #if 0
1836         break;
1837 #endif
1838       }
1839     }
1840     /*     1-wave: shock wave, 3-wave: rarefaction wave */
1841   } else if (iwave == 200) {
1842     /*     case (200) */
1843     i__1 = iterno;
1844     for (i0 = 1; i0 <= i__1; ++i0) {
1845       pst     = *pstar / *pl;
1846       ustarl  = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
1847       zl      = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
1848       d__1    = *pstar / *pr;
1849       d__2    = 1. / *gamr;
1850       *rstarr = *rr * PetscPowScalar(d__1, d__2);
1851       cstarr  = PetscSqrtScalar(*gamr * *pstar / *rstarr);
1852       zr      = *rstarr * cstarr;
1853       ustarr  = *uxr + gascr4 * 2. * (cstarr - cr);
1854       dpstar  = zl * zr * (ustarr - ustarl) / (zl + zr);
1855       *pstar -= dpstar;
1856       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
1857       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1858 #if 0
1859         break;
1860 #endif
1861       }
1862     }
1863     /*     1-wave: shock wave, 3-wave: shock */
1864   } else if (iwave == 300) {
1865     /*     case (300) */
1866     i__1 = iterno;
1867     for (i0 = 1; i0 <= i__1; ++i0) {
1868       pst    = *pstar / *pl;
1869       ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
1870       zl     = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
1871       pst    = *pstar / *pr;
1872       ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
1873       zr     = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
1874       dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
1875       *pstar -= dpstar;
1876       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
1877       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1878 #if 0
1879         break;
1880 #endif
1881       }
1882     }
1883     /*     1-wave: rarefaction wave, 3-wave: shock */
1884   } else if (iwave == 400) {
1885     /*     case (400) */
1886     i__1 = iterno;
1887     for (i0 = 1; i0 <= i__1; ++i0) {
1888       d__1    = *pstar / *pl;
1889       d__2    = 1. / *gaml;
1890       *rstarl = *rl * PetscPowScalar(d__1, d__2);
1891       cstarl  = PetscSqrtScalar(*gaml * *pstar / *rstarl);
1892       ustarl  = *uxl - gascl4 * 2. * (cstarl - cl);
1893       zl      = *rstarl * cstarl;
1894       pst     = *pstar / *pr;
1895       ustarr  = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
1896       zr      = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
1897       dpstar  = zl * zr * (ustarr - ustarl) / (zl + zr);
1898       *pstar -= dpstar;
1899       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
1900       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1901 #if 0
1902               break;
1903 #endif
1904       }
1905     }
1906   }
1907 
1908   *ustar = (zl * ustarr + zr * ustarl) / (zl + zr);
1909   if (PetscRealPart(*pstar) > PetscRealPart(*pl)) {
1910     pst     = *pstar / *pl;
1911     *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + *gaml + 1.) * *rl;
1912   }
1913   if (PetscRealPart(*pstar) > PetscRealPart(*pr)) {
1914     pst     = *pstar / *pr;
1915     *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + *gamr + 1.) * *rr;
1916   }
1917   return iwave;
1918 }
1919 
1920 PetscScalar sign(PetscScalar x) {
1921   if (PetscRealPart(x) > 0) return 1.0;
1922   if (PetscRealPart(x) < 0) return -1.0;
1923   return 0.0;
1924 }
1925 /*        Riemann Solver */
1926 /* -------------------------------------------------------------------- */
1927 int riemannsolver(PetscScalar *xcen, PetscScalar *xp, PetscScalar *dtt, PetscScalar *rl, PetscScalar *uxl, PetscScalar *pl, PetscScalar *utl, PetscScalar *ubl, PetscScalar *gaml, PetscScalar *rho1l, PetscScalar *rr, PetscScalar *uxr, PetscScalar *pr, PetscScalar *utr, PetscScalar *ubr, PetscScalar *gamr, PetscScalar *rho1r, PetscScalar *rx, PetscScalar *uxm, PetscScalar *px, PetscScalar *utx, PetscScalar *ubx, PetscScalar *gam, PetscScalar *rho1) {
1928   /* System generated locals */
1929   PetscScalar d__1, d__2;
1930 
1931   /* Local variables */
1932   static PetscScalar s, c0, p0, r0, u0, w0, x0, x2, ri, cx, sgn0, wsp0, gasc1, gasc2, gasc3, gasc4;
1933   static PetscScalar cstar, pstar, rstar, ustar, xstar, wspst, ushock, streng, rstarl, rstarr, rstars;
1934   int                iwave;
1935 
1936   if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) {
1937     *rx  = *rl;
1938     *px  = *pl;
1939     *uxm = *uxl;
1940     *gam = *gaml;
1941     x2   = *xcen + *uxm * *dtt;
1942 
1943     if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
1944       *utx  = *utr;
1945       *ubx  = *ubr;
1946       *rho1 = *rho1r;
1947     } else {
1948       *utx  = *utl;
1949       *ubx  = *ubl;
1950       *rho1 = *rho1l;
1951     }
1952     return 0;
1953   }
1954   iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar);
1955 
1956   x2     = *xcen + ustar * *dtt;
1957   d__1   = *xp - x2;
1958   sgn0   = sign(d__1);
1959   /*            x is in 3-wave if sgn0 = 1 */
1960   /*            x is in 1-wave if sgn0 = -1 */
1961   r0     = cvmgm_(rl, rr, &sgn0);
1962   p0     = cvmgm_(pl, pr, &sgn0);
1963   u0     = cvmgm_(uxl, uxr, &sgn0);
1964   *gam   = cvmgm_(gaml, gamr, &sgn0);
1965   gasc1  = *gam - 1.;
1966   gasc2  = (*gam + 1.) * .5;
1967   gasc3  = gasc2 / *gam;
1968   gasc4  = 1. / (*gam - 1.);
1969   c0     = PetscSqrtScalar(*gam * p0 / r0);
1970   streng = pstar - p0;
1971   w0     = *gam * r0 * p0 * (gasc3 * streng / p0 + 1.);
1972   rstars = r0 / (1. - r0 * streng / w0);
1973   d__1   = p0 / pstar;
1974   d__2   = -1. / *gam;
1975   rstarr = r0 * PetscPowScalar(d__1, d__2);
1976   rstar  = cvmgm_(&rstarr, &rstars, &streng);
1977   w0     = PetscSqrtScalar(w0);
1978   cstar  = PetscSqrtScalar(*gam * pstar / rstar);
1979   wsp0   = u0 + sgn0 * c0;
1980   wspst  = ustar + sgn0 * cstar;
1981   ushock = ustar + sgn0 * w0 / rstar;
1982   wspst  = cvmgp_(&ushock, &wspst, &streng);
1983   wsp0   = cvmgp_(&ushock, &wsp0, &streng);
1984   x0     = *xcen + wsp0 * *dtt;
1985   xstar  = *xcen + wspst * *dtt;
1986   /*           using gas formula to evaluate rarefaction wave */
1987   /*            ri : reiman invariant */
1988   ri     = u0 - sgn0 * 2. * gasc4 * c0;
1989   cx     = sgn0 * .5 * gasc1 / gasc2 * ((*xp - *xcen) / *dtt - ri);
1990   *uxm   = ri + sgn0 * 2. * gasc4 * cx;
1991   s      = p0 / PetscPowScalar(r0, *gam);
1992   d__1   = cx * cx / (*gam * s);
1993   *rx    = PetscPowScalar(d__1, gasc4);
1994   *px    = cx * cx * *rx / *gam;
1995   d__1   = sgn0 * (x0 - *xp);
1996   *rx    = cvmgp_(rx, &r0, &d__1);
1997   d__1   = sgn0 * (x0 - *xp);
1998   *px    = cvmgp_(px, &p0, &d__1);
1999   d__1   = sgn0 * (x0 - *xp);
2000   *uxm   = cvmgp_(uxm, &u0, &d__1);
2001   d__1   = sgn0 * (xstar - *xp);
2002   *rx    = cvmgm_(rx, &rstar, &d__1);
2003   d__1   = sgn0 * (xstar - *xp);
2004   *px    = cvmgm_(px, &pstar, &d__1);
2005   d__1   = sgn0 * (xstar - *xp);
2006   *uxm   = cvmgm_(uxm, &ustar, &d__1);
2007   if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
2008     *utx  = *utr;
2009     *ubx  = *ubr;
2010     *rho1 = *rho1r;
2011   } else {
2012     *utx  = *utl;
2013     *ubx  = *ubl;
2014     *rho1 = *rho1l;
2015   }
2016   return iwave;
2017 }
2018 int godunovflux(const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, const int *ndim, const PetscReal *gamma) {
2019   /* System generated locals */
2020   int         i__1, iwave;
2021   PetscScalar d__1, d__2, d__3;
2022 
2023   /* Local variables */
2024   static int         k;
2025   static PetscScalar bn[3], fn, ft, tg[3], pl, rl, pm, pr, rr, xp, ubl, ubm, ubr, dtt, unm, tmp, utl, utm, uxl, utr, uxr, gaml, gamm, gamr, xcen, rhom, rho1l, rho1m, rho1r;
2026 
2027   /* Function Body */
2028   xcen = 0.;
2029   xp   = 0.;
2030   i__1 = *ndim;
2031   for (k = 1; k <= i__1; ++k) {
2032     tg[k - 1] = 0.;
2033     bn[k - 1] = 0.;
2034   }
2035   dtt = 1.;
2036   if (*ndim == 3) {
2037     if (nn[0] == 0. && nn[1] == 0.) {
2038       tg[0] = 1.;
2039     } else {
2040       tg[0] = -nn[1];
2041       tg[1] = nn[0];
2042     }
2043     /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2044     /*           tg=tg/tmp */
2045     bn[0] = -nn[2] * tg[1];
2046     bn[1] = nn[2] * tg[0];
2047     bn[2] = nn[0] * tg[1] - nn[1] * tg[0];
2048     /* Computing 2nd power */
2049     d__1  = bn[0];
2050     /* Computing 2nd power */
2051     d__2  = bn[1];
2052     /* Computing 2nd power */
2053     d__3  = bn[2];
2054     tmp   = PetscSqrtScalar(d__1 * d__1 + d__2 * d__2 + d__3 * d__3);
2055     i__1  = *ndim;
2056     for (k = 1; k <= i__1; ++k) bn[k - 1] /= tmp;
2057   } else if (*ndim == 2) {
2058     tg[0] = -nn[1];
2059     tg[1] = nn[0];
2060     /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2061     /*           tg=tg/tmp */
2062     bn[0] = 0.;
2063     bn[1] = 0.;
2064     bn[2] = 1.;
2065   }
2066   rl   = ul[0];
2067   rr   = ur[0];
2068   uxl  = 0.;
2069   uxr  = 0.;
2070   utl  = 0.;
2071   utr  = 0.;
2072   ubl  = 0.;
2073   ubr  = 0.;
2074   i__1 = *ndim;
2075   for (k = 1; k <= i__1; ++k) {
2076     uxl += ul[k] * nn[k - 1];
2077     uxr += ur[k] * nn[k - 1];
2078     utl += ul[k] * tg[k - 1];
2079     utr += ur[k] * tg[k - 1];
2080     ubl += ul[k] * bn[k - 1];
2081     ubr += ur[k] * bn[k - 1];
2082   }
2083   uxl /= rl;
2084   uxr /= rr;
2085   utl /= rl;
2086   utr /= rr;
2087   ubl /= rl;
2088   ubr /= rr;
2089 
2090   gaml  = *gamma;
2091   gamr  = *gamma;
2092   /* Computing 2nd power */
2093   d__1  = uxl;
2094   /* Computing 2nd power */
2095   d__2  = utl;
2096   /* Computing 2nd power */
2097   d__3  = ubl;
2098   pl    = (*gamma - 1.) * (ul[*ndim + 1] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2099   /* Computing 2nd power */
2100   d__1  = uxr;
2101   /* Computing 2nd power */
2102   d__2  = utr;
2103   /* Computing 2nd power */
2104   d__3  = ubr;
2105   pr    = (*gamma - 1.) * (ur[*ndim + 1] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2106   rho1l = rl;
2107   rho1r = rr;
2108 
2109   iwave = riemannsolver(&xcen, &xp, &dtt, &rl, &uxl, &pl, &utl, &ubl, &gaml, &rho1l, &rr, &uxr, &pr, &utr, &ubr, &gamr, &rho1r, &rhom, &unm, &pm, &utm, &ubm, &gamm, &rho1m);
2110 
2111   flux[0] = rhom * unm;
2112   fn      = rhom * unm * unm + pm;
2113   ft      = rhom * unm * utm;
2114   /*           flux(2)=fn*nn(1)+ft*nn(2) */
2115   /*           flux(3)=fn*tg(1)+ft*tg(2) */
2116   flux[1] = fn * nn[0] + ft * tg[0];
2117   flux[2] = fn * nn[1] + ft * tg[1];
2118   /*           flux(2)=rhom*unm*(unm)+pm */
2119   /*           flux(3)=rhom*(unm)*utm */
2120   if (*ndim == 3) flux[3] = rhom * unm * ubm;
2121   flux[*ndim + 1] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm;
2122   return iwave;
2123 } /* godunovflux_ */
2124 
2125 /* Subroutine to set up the initial conditions for the */
2126 /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */
2127 /* ----------------------------------------------------------------------- */
2128 int projecteqstate(PetscReal wc[], const PetscReal ueq[], PetscReal lv[][3]) {
2129   int j, k;
2130   /*      Wc=matmul(lv,Ueq) 3 vars */
2131   for (k = 0; k < 3; ++k) {
2132     wc[k] = 0.;
2133     for (j = 0; j < 3; ++j) wc[k] += lv[k][j] * ueq[j];
2134   }
2135   return 0;
2136 }
2137 /* ----------------------------------------------------------------------- */
2138 int projecttoprim(PetscReal v[], const PetscReal wc[], PetscReal rv[][3]) {
2139   int k, j;
2140   /*      V=matmul(rv,WC) 3 vars */
2141   for (k = 0; k < 3; ++k) {
2142     v[k] = 0.;
2143     for (j = 0; j < 3; ++j) v[k] += rv[k][j] * wc[j];
2144   }
2145   return 0;
2146 }
2147 /* ---------------------------------------------------------------------- */
2148 int eigenvectors(PetscReal rv[][3], PetscReal lv[][3], const PetscReal ueq[], PetscReal gamma) {
2149   int       j, k;
2150   PetscReal rho, csnd, p0;
2151   /* PetscScalar u; */
2152 
2153   for (k = 0; k < 3; ++k)
2154     for (j = 0; j < 3; ++j) {
2155       lv[k][j] = 0.;
2156       rv[k][j] = 0.;
2157     }
2158   rho      = ueq[0];
2159   /* u = ueq[1]; */
2160   p0       = ueq[2];
2161   csnd     = PetscSqrtReal(gamma * p0 / rho);
2162   lv[0][1] = rho * .5;
2163   lv[0][2] = -.5 / csnd;
2164   lv[1][0] = csnd;
2165   lv[1][2] = -1. / csnd;
2166   lv[2][1] = rho * .5;
2167   lv[2][2] = .5 / csnd;
2168   rv[0][0] = -1. / csnd;
2169   rv[1][0] = 1. / rho;
2170   rv[2][0] = -csnd;
2171   rv[0][1] = 1. / csnd;
2172   rv[0][2] = 1. / csnd;
2173   rv[1][2] = 1. / rho;
2174   rv[2][2] = csnd;
2175   return 0;
2176 }
2177 
2178 int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx) {
2179   PetscReal p0, u0, wcp[3], wc[3];
2180   PetscReal lv[3][3];
2181   PetscReal vp[3];
2182   PetscReal rv[3][3];
2183   PetscReal eps, ueq[3], rho0, twopi;
2184 
2185   /* Function Body */
2186   twopi  = 2. * PETSC_PI;
2187   eps    = 1e-4;    /* perturbation */
2188   rho0   = 1e3;     /* density of water */
2189   p0     = 101325.; /* init pressure of 1 atm (?) */
2190   u0     = 0.;
2191   ueq[0] = rho0;
2192   ueq[1] = u0;
2193   ueq[2] = p0;
2194   /* Project initial state to characteristic variables */
2195   eigenvectors(rv, lv, ueq, gamma);
2196   projecteqstate(wc, ueq, lv);
2197   wcp[0] = wc[0];
2198   wcp[1] = wc[1];
2199   wcp[2] = wc[2] + eps * PetscCosReal(coord[0] * 2. * twopi / Lx);
2200   projecttoprim(vp, wcp, rv);
2201   ux->r     = vp[0];         /* density */
2202   ux->ru[0] = vp[0] * vp[1]; /* x momentum */
2203   ux->ru[1] = 0.;
2204 #if defined DIM > 2
2205   if (dim > 2) ux->ru[2] = 0.;
2206 #endif
2207   /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */
2208   ux->E = vp[2] / (gamma - 1.) + 0.5 * vp[0] * vp[1] * vp[1];
2209   return 0;
2210 }
2211 
2212 /*TEST
2213 
2214   testset:
2215     args: -dm_plex_adj_cone -dm_plex_adj_closure 0
2216 
2217     test:
2218       suffix: adv_2d_tri_0
2219       requires: triangle
2220       TODO: how did this ever get in main when there is no support for this
2221       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
2222 
2223     test:
2224       suffix: adv_2d_tri_1
2225       requires: triangle
2226       TODO: how did this ever get in main when there is no support for this
2227       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
2228 
2229     test:
2230       suffix: tut_1
2231       requires: exodusii
2232       nsize: 1
2233       args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2234 
2235     test:
2236       suffix: tut_2
2237       requires: exodusii
2238       nsize: 1
2239       args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw
2240 
2241     test:
2242       suffix: tut_3
2243       requires: exodusii
2244       nsize: 4
2245       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
2246 
2247     test:
2248       suffix: tut_4
2249       requires: exodusii
2250       nsize: 4
2251       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
2252 
2253   testset:
2254     args: -dm_plex_adj_cone -dm_plex_adj_closure 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1
2255 
2256     # 2D Advection 0-10
2257     test:
2258       suffix: 0
2259       requires: exodusii
2260       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2261 
2262     test:
2263       suffix: 1
2264       requires: exodusii
2265       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
2266 
2267     test:
2268       suffix: 2
2269       requires: exodusii
2270       nsize: 2
2271       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2272 
2273     test:
2274       suffix: 3
2275       requires: exodusii
2276       nsize: 2
2277       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
2278 
2279     test:
2280       suffix: 4
2281       requires: exodusii
2282       nsize: 8
2283       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo
2284 
2285     test:
2286       suffix: 5
2287       requires: exodusii
2288       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1
2289 
2290     test:
2291       suffix: 7
2292       requires: exodusii
2293       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
2294 
2295     test:
2296       suffix: 8
2297       requires: exodusii
2298       nsize: 2
2299       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
2300 
2301     test:
2302       suffix: 9
2303       requires: exodusii
2304       nsize: 8
2305       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
2306 
2307     test:
2308       suffix: 10
2309       requires: exodusii
2310       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo
2311 
2312     # 2D Shallow water
2313     test:
2314       suffix: sw_0
2315       requires: exodusii
2316       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -bc_wall 100,101 -physics sw -ufv_cfl 5 -petscfv_type   leastsquares -petsclimiter_type sin -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 -monitor height,energy
2317 
2318     test:
2319       suffix: sw_hll
2320       args: -ufv_vtk_interval 0 -bc_wall 1,2,3,4 -physics sw -ufv_cfl 3 -petscfv_type leastsquares -petsclimiter_type sin -ts_max_steps 5 -ts_ssp_type rks2 -ts_ssp_nstages 10 -monitor height,energy -grid_bounds 0,5,0,5 -dm_plex_box_faces 25,25 -sw_riemann hll
2321 
2322     # 2D Advection: p4est
2323     test:
2324       suffix: p4est_advec_2d
2325       requires: p4est
2326       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
2327 
2328     # Advection in a box
2329     test:
2330       suffix: adv_2d_quad_0
2331       args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2332 
2333     test:
2334       suffix: adv_2d_quad_1
2335       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
2336       timeoutfactor: 3
2337 
2338     test:
2339       suffix: adv_2d_quad_p4est_0
2340       requires: p4est
2341       args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2342 
2343     test:
2344       suffix: adv_2d_quad_p4est_1
2345       requires: p4est
2346       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
2347       timeoutfactor: 3
2348 
2349     test:
2350       suffix: adv_2d_quad_p4est_adapt_0
2351       requires: p4est !__float128 #broken for quad precision
2352       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
2353       timeoutfactor: 3
2354 
2355     test:
2356       suffix: adv_0
2357       requires: exodusii
2358       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
2359 
2360     test:
2361       suffix: shock_0
2362       requires: p4est !single !complex
2363       args: -dm_plex_box_faces 2,1 -grid_bounds -1,1.,0.,1 -grid_skew_60 \
2364       -dm_type p4est -dm_forest_partition_overlap 1 -dm_forest_maximum_refinement 6 -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 \
2365       -ufv_use_amr -refine_vec_tagger_box 0.5,inf -coarsen_vec_tagger_box 0,1.e-2 -refine_tag_view -coarsen_tag_view \
2366       -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. \
2367       -petscfv_type leastsquares -petsclimiter_type minmod -petscfv_compute_gradients 0 \
2368       -ts_max_time 0.5 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
2369       -ufv_vtk_basename ${wPETSC_DIR}/ex11 -ufv_vtk_interval 0 -monitor density,energy
2370       timeoutfactor: 3
2371 
2372     # Test GLVis visualization of PetscFV fields
2373     test:
2374       suffix: glvis_adv_2d_tet
2375       args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 \
2376             -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0 \
2377             -ts_monitor_solution glvis: -ts_max_steps 0
2378 
2379     test:
2380       suffix: glvis_adv_2d_quad
2381       args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 -bc_inflow 1,2,4 -bc_outflow 3 \
2382             -dm_refine 5 -dm_plex_separate_marker \
2383             -ts_monitor_solution glvis: -ts_max_steps 0
2384 
2385 TEST*/
2386