xref: /petsc/src/ts/tutorials/ex11.c (revision 6c2b77d522d8aa5c8b27f04fddd7150d0d6755fb)
1c4762a1bSJed Brown static char help[] = "Second Order TVD Finite Volume Example.\n";
2c4762a1bSJed Brown /*F
3c4762a1bSJed Brown 
4c4762a1bSJed Brown We use a second order TVD finite volume method to evolve a system of PDEs. Our simple upwinded residual evaluation loops
5c4762a1bSJed Brown over all mesh faces and uses a Riemann solver to produce the flux given the face geometry and cell values,
6c4762a1bSJed Brown \begin{equation}
7c4762a1bSJed Brown   f_i = \mathrm{riemann}(\mathrm{phys}, p_\mathrm{centroid}, \hat n, x^L, x^R)
8c4762a1bSJed Brown \end{equation}
9c4762a1bSJed Brown and then update the cell values given the cell volume.
10c4762a1bSJed Brown \begin{eqnarray}
11c4762a1bSJed Brown     f^L_i &-=& \frac{f_i}{vol^L} \\
12c4762a1bSJed Brown     f^R_i &+=& \frac{f_i}{vol^R}
13c4762a1bSJed Brown \end{eqnarray}
14c4762a1bSJed Brown 
15c4762a1bSJed Brown As an example, we can consider the shallow water wave equation,
16c4762a1bSJed Brown \begin{eqnarray}
17c4762a1bSJed Brown      h_t + \nabla\cdot \left( uh                              \right) &=& 0 \\
18c4762a1bSJed Brown   (uh)_t + \nabla\cdot \left( u\otimes uh + \frac{g h^2}{2} I \right) &=& 0
19c4762a1bSJed Brown \end{eqnarray}
20c4762a1bSJed Brown where $h$ is wave height, $u$ is wave velocity, and $g$ is the acceleration due to gravity.
21c4762a1bSJed Brown 
22c4762a1bSJed Brown A representative Riemann solver for the shallow water equations is given in the PhysicsRiemann_SW() function,
23c4762a1bSJed Brown \begin{eqnarray}
24c4762a1bSJed Brown   f^{L,R}_h    &=& uh^{L,R} \cdot \hat n \\
25c4762a1bSJed Brown   f^{L,R}_{uh} &=& \frac{f^{L,R}_h}{h^{L,R}} uh^{L,R} + g (h^{L,R})^2 \hat n \\
26c4762a1bSJed Brown   c^{L,R}      &=& \sqrt{g h^{L,R}} \\
27c4762a1bSJed Brown   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) \\
28c4762a1bSJed Brown   f_i          &=& \frac{A_\mathrm{face}}{2} \left( f^L_i + f^R_i + s \left( x^L_i - x^R_i \right) \right)
29c4762a1bSJed Brown \end{eqnarray}
30c4762a1bSJed Brown where $c$ is the local gravity wave speed and $f_i$ is a Rusanov flux.
31c4762a1bSJed Brown 
32c4762a1bSJed Brown The more sophisticated residual evaluation in RHSFunctionLocal_LS() uses a least-squares fit to a quadratic polynomial
33c4762a1bSJed Brown over a neighborhood of the given element.
34c4762a1bSJed Brown 
35c4762a1bSJed Brown The mesh is read in from an ExodusII file, usually generated by Cubit.
36c4762a1bSJed Brown F*/
37c4762a1bSJed Brown #include <petscdmplex.h>
38c4762a1bSJed Brown #include <petscdmforest.h>
39c4762a1bSJed Brown #include <petscds.h>
40c4762a1bSJed Brown #include <petscts.h>
41c4762a1bSJed Brown 
42c4762a1bSJed Brown #define DIM 2 /* Geometric dimension */
43c4762a1bSJed Brown 
448d2c18caSMukkund Sunjii static PetscFunctionList PhysicsList, PhysicsRiemannList_SW;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /* Represents continuum physical equations. */
47c4762a1bSJed Brown typedef struct _n_Physics *Physics;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /* Physical model includes boundary conditions, initial conditions, and functionals of interest. It is
50c4762a1bSJed Brown  * discretization-independent, but its members depend on the scenario being solved. */
51c4762a1bSJed Brown typedef struct _n_Model *Model;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown /* 'User' implements a discretization of a continuous model. */
54c4762a1bSJed Brown typedef struct _n_User *User;
55c4762a1bSJed Brown typedef PetscErrorCode (*SolutionFunction)(Model, PetscReal, const PetscReal *, PetscScalar *, void *);
5645480ffeSMatthew G. Knepley typedef PetscErrorCode (*SetUpBCFunction)(DM, PetscDS, Physics);
57c4762a1bSJed Brown typedef PetscErrorCode (*FunctionalFunction)(Model, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *);
58c4762a1bSJed Brown typedef PetscErrorCode (*SetupFields)(Physics, PetscSection);
59c4762a1bSJed Brown static PetscErrorCode ModelSolutionSetDefault(Model, SolutionFunction, void *);
60c4762a1bSJed Brown static PetscErrorCode ModelFunctionalRegister(Model, const char *, PetscInt *, FunctionalFunction, void *);
61c4762a1bSJed Brown static PetscErrorCode OutputVTK(DM, const char *, PetscViewer *);
62c4762a1bSJed Brown 
63c4762a1bSJed Brown struct FieldDescription {
64c4762a1bSJed Brown   const char *name;
65c4762a1bSJed Brown   PetscInt    dof;
66c4762a1bSJed Brown };
67c4762a1bSJed Brown 
68c4762a1bSJed Brown typedef struct _n_FunctionalLink *FunctionalLink;
69c4762a1bSJed Brown struct _n_FunctionalLink {
70c4762a1bSJed Brown   char              *name;
71c4762a1bSJed Brown   FunctionalFunction func;
72c4762a1bSJed Brown   void              *ctx;
73c4762a1bSJed Brown   PetscInt           offset;
74c4762a1bSJed Brown   FunctionalLink     next;
75c4762a1bSJed Brown };
76c4762a1bSJed Brown 
77c4762a1bSJed Brown struct _n_Physics {
78c4762a1bSJed Brown   PetscRiemannFunc               riemann;
79c4762a1bSJed Brown   PetscInt                       dof;      /* number of degrees of freedom per cell */
80c4762a1bSJed Brown   PetscReal                      maxspeed; /* kludge to pick initial time step, need to add monitoring and step control */
81c4762a1bSJed Brown   void                          *data;
82c4762a1bSJed Brown   PetscInt                       nfields;
83c4762a1bSJed Brown   const struct FieldDescription *field_desc;
84c4762a1bSJed Brown };
85c4762a1bSJed Brown 
86c4762a1bSJed Brown struct _n_Model {
87c4762a1bSJed Brown   MPI_Comm         comm; /* Does not do collective communicaton, but some error conditions can be collective */
88c4762a1bSJed Brown   Physics          physics;
89c4762a1bSJed Brown   FunctionalLink   functionalRegistry;
90c4762a1bSJed Brown   PetscInt         maxComputed;
91c4762a1bSJed Brown   PetscInt         numMonitored;
92c4762a1bSJed Brown   FunctionalLink  *functionalMonitored;
93c4762a1bSJed Brown   PetscInt         numCall;
94c4762a1bSJed Brown   FunctionalLink  *functionalCall;
95c4762a1bSJed Brown   SolutionFunction solution;
96c4762a1bSJed Brown   SetUpBCFunction  setupbc;
97c4762a1bSJed Brown   void            *solutionctx;
98c4762a1bSJed Brown   PetscReal        maxspeed; /* estimate of global maximum speed (for CFL calculation) */
99c4762a1bSJed Brown   PetscReal        bounds[2 * DIM];
100c4762a1bSJed Brown   PetscErrorCode (*errorIndicator)(PetscInt, PetscReal, PetscInt, const PetscScalar[], const PetscScalar[], PetscReal *, void *);
101c4762a1bSJed Brown   void *errorCtx;
102c4762a1bSJed Brown };
103c4762a1bSJed Brown 
104c4762a1bSJed Brown struct _n_User {
105c4762a1bSJed Brown   PetscInt  vtkInterval;                        /* For monitor */
106c4762a1bSJed Brown   char      outputBasename[PETSC_MAX_PATH_LEN]; /* Basename for output files */
107c4762a1bSJed Brown   PetscInt  monitorStepOffset;
108c4762a1bSJed Brown   Model     model;
109c4762a1bSJed Brown   PetscBool vtkmon;
110c4762a1bSJed Brown };
111c4762a1bSJed Brown 
112d71ae5a4SJacob Faibussowitsch static inline PetscReal DotDIMReal(const PetscReal *x, const PetscReal *y)
113d71ae5a4SJacob Faibussowitsch {
114c4762a1bSJed Brown   PetscInt  i;
115c4762a1bSJed Brown   PetscReal prod = 0.0;
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   for (i = 0; i < DIM; i++) prod += x[i] * y[i];
118c4762a1bSJed Brown   return prod;
119c4762a1bSJed Brown }
120d71ae5a4SJacob Faibussowitsch static inline PetscReal NormDIM(const PetscReal *x)
121d71ae5a4SJacob Faibussowitsch {
1229371c9d4SSatish Balay   return PetscSqrtReal(PetscAbsReal(DotDIMReal(x, x)));
1239371c9d4SSatish Balay }
124c4762a1bSJed Brown 
125d71ae5a4SJacob Faibussowitsch static inline PetscReal Dot2Real(const PetscReal *x, const PetscReal *y)
126d71ae5a4SJacob Faibussowitsch {
1279371c9d4SSatish Balay   return x[0] * y[0] + x[1] * y[1];
1289371c9d4SSatish Balay }
129d71ae5a4SJacob Faibussowitsch static inline PetscReal Norm2Real(const PetscReal *x)
130d71ae5a4SJacob Faibussowitsch {
1319371c9d4SSatish Balay   return PetscSqrtReal(PetscAbsReal(Dot2Real(x, x)));
1329371c9d4SSatish Balay }
133d71ae5a4SJacob Faibussowitsch static inline void Normalize2Real(PetscReal *x)
134d71ae5a4SJacob Faibussowitsch {
1359371c9d4SSatish Balay   PetscReal a = 1. / Norm2Real(x);
1369371c9d4SSatish Balay   x[0] *= a;
1379371c9d4SSatish Balay   x[1] *= a;
1389371c9d4SSatish Balay }
139d71ae5a4SJacob Faibussowitsch static inline void Waxpy2Real(PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w)
140d71ae5a4SJacob Faibussowitsch {
1419371c9d4SSatish Balay   w[0] = a * x[0] + y[0];
1429371c9d4SSatish Balay   w[1] = a * x[1] + y[1];
1439371c9d4SSatish Balay }
144d71ae5a4SJacob Faibussowitsch static inline void Scale2Real(PetscReal a, const PetscReal *x, PetscReal *y)
145d71ae5a4SJacob Faibussowitsch {
1469371c9d4SSatish Balay   y[0] = a * x[0];
1479371c9d4SSatish Balay   y[1] = a * x[1];
1489371c9d4SSatish Balay }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown /******************* Advect ********************/
1519371c9d4SSatish Balay typedef enum {
1529371c9d4SSatish Balay   ADVECT_SOL_TILTED,
1539371c9d4SSatish Balay   ADVECT_SOL_BUMP,
1549371c9d4SSatish Balay   ADVECT_SOL_BUMP_CAVITY
1559371c9d4SSatish Balay } AdvectSolType;
156c4762a1bSJed Brown static const char *const AdvectSolTypes[] = {"TILTED", "BUMP", "BUMP_CAVITY", "AdvectSolType", "ADVECT_SOL_", 0};
1579371c9d4SSatish Balay typedef enum {
1589371c9d4SSatish Balay   ADVECT_SOL_BUMP_CONE,
1599371c9d4SSatish Balay   ADVECT_SOL_BUMP_COS
1609371c9d4SSatish Balay } AdvectSolBumpType;
161c4762a1bSJed Brown static const char *const AdvectSolBumpTypes[] = {"CONE", "COS", "AdvectSolBumpType", "ADVECT_SOL_BUMP_", 0};
162c4762a1bSJed Brown 
163c4762a1bSJed Brown typedef struct {
164c4762a1bSJed Brown   PetscReal wind[DIM];
165c4762a1bSJed Brown } Physics_Advect_Tilted;
166c4762a1bSJed Brown typedef struct {
167c4762a1bSJed Brown   PetscReal         center[DIM];
168c4762a1bSJed Brown   PetscReal         radius;
169c4762a1bSJed Brown   AdvectSolBumpType type;
170c4762a1bSJed Brown } Physics_Advect_Bump;
171c4762a1bSJed Brown 
172c4762a1bSJed Brown typedef struct {
173c4762a1bSJed Brown   PetscReal     inflowState;
174c4762a1bSJed Brown   AdvectSolType soltype;
1759371c9d4SSatish Balay   union
1769371c9d4SSatish Balay   {
177c4762a1bSJed Brown     Physics_Advect_Tilted tilted;
178c4762a1bSJed Brown     Physics_Advect_Bump   bump;
179c4762a1bSJed Brown   } sol;
180c4762a1bSJed Brown   struct {
181c4762a1bSJed Brown     PetscInt Solution;
182c4762a1bSJed Brown     PetscInt Error;
183c4762a1bSJed Brown   } functional;
184c4762a1bSJed Brown } Physics_Advect;
185c4762a1bSJed Brown 
1869371c9d4SSatish Balay static const struct FieldDescription PhysicsFields_Advect[] = {
1879371c9d4SSatish Balay   {"U",  1},
1889371c9d4SSatish Balay   {NULL, 0}
1899371c9d4SSatish Balay };
190c4762a1bSJed Brown 
191d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsBoundary_Advect_Inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
192d71ae5a4SJacob Faibussowitsch {
193c4762a1bSJed Brown   Physics         phys   = (Physics)ctx;
194c4762a1bSJed Brown   Physics_Advect *advect = (Physics_Advect *)phys->data;
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   PetscFunctionBeginUser;
197c4762a1bSJed Brown   xG[0] = advect->inflowState;
198c4762a1bSJed Brown   PetscFunctionReturn(0);
199c4762a1bSJed Brown }
200c4762a1bSJed Brown 
201d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
202d71ae5a4SJacob Faibussowitsch {
203c4762a1bSJed Brown   PetscFunctionBeginUser;
204c4762a1bSJed Brown   xG[0] = xI[0];
205c4762a1bSJed Brown   PetscFunctionReturn(0);
206c4762a1bSJed Brown }
207c4762a1bSJed Brown 
208d71ae5a4SJacob Faibussowitsch 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)
209d71ae5a4SJacob Faibussowitsch {
210c4762a1bSJed Brown   Physics_Advect *advect = (Physics_Advect *)phys->data;
211c4762a1bSJed Brown   PetscReal       wind[DIM], wn;
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   switch (advect->soltype) {
214c4762a1bSJed Brown   case ADVECT_SOL_TILTED: {
215c4762a1bSJed Brown     Physics_Advect_Tilted *tilted = &advect->sol.tilted;
216c4762a1bSJed Brown     wind[0]                       = tilted->wind[0];
217c4762a1bSJed Brown     wind[1]                       = tilted->wind[1];
218c4762a1bSJed Brown   } break;
219c4762a1bSJed Brown   case ADVECT_SOL_BUMP:
220c4762a1bSJed Brown     wind[0] = -qp[1];
221c4762a1bSJed Brown     wind[1] = qp[0];
222c4762a1bSJed Brown     break;
2239371c9d4SSatish Balay   case ADVECT_SOL_BUMP_CAVITY: {
224c4762a1bSJed Brown     PetscInt  i;
225c4762a1bSJed Brown     PetscReal comp2[3] = {0., 0., 0.}, rad2;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown     rad2 = 0.;
228c4762a1bSJed Brown     for (i = 0; i < dim; i++) {
229c4762a1bSJed Brown       comp2[i] = qp[i] * qp[i];
230c4762a1bSJed Brown       rad2 += comp2[i];
231c4762a1bSJed Brown     }
232c4762a1bSJed Brown 
233c4762a1bSJed Brown     wind[0] = -qp[1];
234c4762a1bSJed Brown     wind[1] = qp[0];
235c4762a1bSJed Brown     if (rad2 > 1.) {
236c4762a1bSJed Brown       PetscInt  maxI     = 0;
237c4762a1bSJed Brown       PetscReal maxComp2 = comp2[0];
238c4762a1bSJed Brown 
239c4762a1bSJed Brown       for (i = 1; i < dim; i++) {
240c4762a1bSJed Brown         if (comp2[i] > maxComp2) {
241c4762a1bSJed Brown           maxI     = i;
242c4762a1bSJed Brown           maxComp2 = comp2[i];
243c4762a1bSJed Brown         }
244c4762a1bSJed Brown       }
245c4762a1bSJed Brown       wind[maxI] = 0.;
246c4762a1bSJed Brown     }
2479371c9d4SSatish Balay   } break;
2489371c9d4SSatish Balay   default: {
249c4762a1bSJed Brown     PetscInt i;
250c4762a1bSJed Brown     for (i = 0; i < DIM; ++i) wind[i] = 0.0;
251c4762a1bSJed Brown   }
25298921bdaSJacob Faibussowitsch     /* default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for solution type %s",AdvectSolBumpTypes[advect->soltype]); */
253c4762a1bSJed Brown   }
254c4762a1bSJed Brown   wn      = Dot2Real(wind, n);
255c4762a1bSJed Brown   flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn;
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsSolution_Advect(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
259d71ae5a4SJacob Faibussowitsch {
260c4762a1bSJed Brown   Physics         phys   = (Physics)ctx;
261c4762a1bSJed Brown   Physics_Advect *advect = (Physics_Advect *)phys->data;
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   PetscFunctionBeginUser;
264c4762a1bSJed Brown   switch (advect->soltype) {
265c4762a1bSJed Brown   case ADVECT_SOL_TILTED: {
266c4762a1bSJed Brown     PetscReal              x0[DIM];
267c4762a1bSJed Brown     Physics_Advect_Tilted *tilted = &advect->sol.tilted;
268c4762a1bSJed Brown     Waxpy2Real(-time, tilted->wind, x, x0);
269c4762a1bSJed Brown     if (x0[1] > 0) u[0] = 1. * x[0] + 3. * x[1];
270c4762a1bSJed Brown     else u[0] = advect->inflowState;
271c4762a1bSJed Brown   } break;
272c4762a1bSJed Brown   case ADVECT_SOL_BUMP_CAVITY:
273c4762a1bSJed Brown   case ADVECT_SOL_BUMP: {
274c4762a1bSJed Brown     Physics_Advect_Bump *bump = &advect->sol.bump;
275c4762a1bSJed Brown     PetscReal            x0[DIM], v[DIM], r, cost, sint;
276c4762a1bSJed Brown     cost  = PetscCosReal(time);
277c4762a1bSJed Brown     sint  = PetscSinReal(time);
278c4762a1bSJed Brown     x0[0] = cost * x[0] + sint * x[1];
279c4762a1bSJed Brown     x0[1] = -sint * x[0] + cost * x[1];
280c4762a1bSJed Brown     Waxpy2Real(-1, bump->center, x0, v);
281c4762a1bSJed Brown     r = Norm2Real(v);
282c4762a1bSJed Brown     switch (bump->type) {
283d71ae5a4SJacob Faibussowitsch     case ADVECT_SOL_BUMP_CONE:
284d71ae5a4SJacob Faibussowitsch       u[0] = PetscMax(1 - r / bump->radius, 0);
285d71ae5a4SJacob Faibussowitsch       break;
286d71ae5a4SJacob Faibussowitsch     case ADVECT_SOL_BUMP_COS:
287d71ae5a4SJacob Faibussowitsch       u[0] = 0.5 + 0.5 * PetscCosReal(PetscMin(r / bump->radius, 1) * PETSC_PI);
288d71ae5a4SJacob Faibussowitsch       break;
289c4762a1bSJed Brown     }
290c4762a1bSJed Brown   } break;
291d71ae5a4SJacob Faibussowitsch   default:
292d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type");
293c4762a1bSJed Brown   }
294c4762a1bSJed Brown   PetscFunctionReturn(0);
295c4762a1bSJed Brown }
296c4762a1bSJed Brown 
297d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsFunctional_Advect(Model mod, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx)
298d71ae5a4SJacob Faibussowitsch {
299c4762a1bSJed Brown   Physics         phys      = (Physics)ctx;
300c4762a1bSJed Brown   Physics_Advect *advect    = (Physics_Advect *)phys->data;
301391faea4SSatish Balay   PetscScalar     yexact[1] = {0.0};
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   PetscFunctionBeginUser;
3049566063dSJacob Faibussowitsch   PetscCall(PhysicsSolution_Advect(mod, time, x, yexact, phys));
305c4762a1bSJed Brown   f[advect->functional.Solution] = PetscRealPart(y[0]);
306c4762a1bSJed Brown   f[advect->functional.Error]    = PetscAbsScalar(y[0] - yexact[0]);
307c4762a1bSJed Brown   PetscFunctionReturn(0);
308c4762a1bSJed Brown }
309c4762a1bSJed Brown 
310d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetUpBC_Advect(DM dm, PetscDS prob, Physics phys)
311d71ae5a4SJacob Faibussowitsch {
312c4762a1bSJed Brown   const PetscInt inflowids[] = {100, 200, 300}, outflowids[] = {101};
31345480ffeSMatthew G. Knepley   DMLabel        label;
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   PetscFunctionBeginUser;
316c4762a1bSJed Brown   /* Register "canned" boundary conditions and defaults for where to apply. */
3179566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "Face Sets", &label));
318dd39110bSPierre Jolivet   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));
319dd39110bSPierre Jolivet   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));
320c4762a1bSJed Brown   PetscFunctionReturn(0);
321c4762a1bSJed Brown }
322c4762a1bSJed Brown 
323d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCreate_Advect(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject)
324d71ae5a4SJacob Faibussowitsch {
325c4762a1bSJed Brown   Physics_Advect *advect;
326c4762a1bSJed Brown 
327c4762a1bSJed Brown   PetscFunctionBeginUser;
328c4762a1bSJed Brown   phys->field_desc = PhysicsFields_Advect;
329c4762a1bSJed Brown   phys->riemann    = (PetscRiemannFunc)PhysicsRiemann_Advect;
3309566063dSJacob Faibussowitsch   PetscCall(PetscNew(&advect));
331c4762a1bSJed Brown   phys->data   = advect;
332c4762a1bSJed Brown   mod->setupbc = SetUpBC_Advect;
333c4762a1bSJed Brown 
334d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Advect options");
335c4762a1bSJed Brown   {
336c4762a1bSJed Brown     PetscInt two = 2, dof = 1;
337c4762a1bSJed Brown     advect->soltype = ADVECT_SOL_TILTED;
3389566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-advect_sol_type", "solution type", "", AdvectSolTypes, (PetscEnum)advect->soltype, (PetscEnum *)&advect->soltype, NULL));
339c4762a1bSJed Brown     switch (advect->soltype) {
340c4762a1bSJed Brown     case ADVECT_SOL_TILTED: {
341c4762a1bSJed Brown       Physics_Advect_Tilted *tilted = &advect->sol.tilted;
342c4762a1bSJed Brown       two                           = 2;
343c4762a1bSJed Brown       tilted->wind[0]               = 0.0;
344c4762a1bSJed Brown       tilted->wind[1]               = 1.0;
3459566063dSJacob Faibussowitsch       PetscCall(PetscOptionsRealArray("-advect_tilted_wind", "background wind vx,vy", "", tilted->wind, &two, NULL));
346c4762a1bSJed Brown       advect->inflowState = -2.0;
3479566063dSJacob Faibussowitsch       PetscCall(PetscOptionsRealArray("-advect_tilted_inflow", "Inflow state", "", &advect->inflowState, &dof, NULL));
348c4762a1bSJed Brown       phys->maxspeed = Norm2Real(tilted->wind);
349c4762a1bSJed Brown     } break;
350c4762a1bSJed Brown     case ADVECT_SOL_BUMP_CAVITY:
351c4762a1bSJed Brown     case ADVECT_SOL_BUMP: {
352c4762a1bSJed Brown       Physics_Advect_Bump *bump = &advect->sol.bump;
353c4762a1bSJed Brown       two                       = 2;
354c4762a1bSJed Brown       bump->center[0]           = 2.;
355c4762a1bSJed Brown       bump->center[1]           = 0.;
3569566063dSJacob Faibussowitsch       PetscCall(PetscOptionsRealArray("-advect_bump_center", "location of center of bump x,y", "", bump->center, &two, NULL));
357c4762a1bSJed Brown       bump->radius = 0.9;
3589566063dSJacob Faibussowitsch       PetscCall(PetscOptionsReal("-advect_bump_radius", "radius of bump", "", bump->radius, &bump->radius, NULL));
359c4762a1bSJed Brown       bump->type = ADVECT_SOL_BUMP_CONE;
3609566063dSJacob Faibussowitsch       PetscCall(PetscOptionsEnum("-advect_bump_type", "type of bump", "", AdvectSolBumpTypes, (PetscEnum)bump->type, (PetscEnum *)&bump->type, NULL));
361c4762a1bSJed Brown       phys->maxspeed = 3.; /* radius of mesh, kludge */
362c4762a1bSJed Brown     } break;
363c4762a1bSJed Brown     }
364c4762a1bSJed Brown   }
365d0609cedSBarry Smith   PetscOptionsHeadEnd();
366c4762a1bSJed Brown   /* Initial/transient solution with default boundary conditions */
3679566063dSJacob Faibussowitsch   PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_Advect, phys));
368c4762a1bSJed Brown   /* Register "canned" functionals */
3699566063dSJacob Faibussowitsch   PetscCall(ModelFunctionalRegister(mod, "Solution", &advect->functional.Solution, PhysicsFunctional_Advect, phys));
3709566063dSJacob Faibussowitsch   PetscCall(ModelFunctionalRegister(mod, "Error", &advect->functional.Error, PhysicsFunctional_Advect, phys));
371c4762a1bSJed Brown   PetscFunctionReturn(0);
372c4762a1bSJed Brown }
373c4762a1bSJed Brown 
374c4762a1bSJed Brown /******************* Shallow Water ********************/
375c4762a1bSJed Brown typedef struct {
376c4762a1bSJed Brown   PetscReal gravity;
377c4762a1bSJed Brown   PetscReal boundaryHeight;
378c4762a1bSJed Brown   struct {
379c4762a1bSJed Brown     PetscInt Height;
380c4762a1bSJed Brown     PetscInt Speed;
381c4762a1bSJed Brown     PetscInt Energy;
382c4762a1bSJed Brown   } functional;
383c4762a1bSJed Brown } Physics_SW;
384c4762a1bSJed Brown typedef struct {
385c4762a1bSJed Brown   PetscReal h;
386c4762a1bSJed Brown   PetscReal uh[DIM];
387c4762a1bSJed Brown } SWNode;
3889371c9d4SSatish Balay typedef union
3899371c9d4SSatish Balay {
390c4762a1bSJed Brown   SWNode    swnode;
391c4762a1bSJed Brown   PetscReal vals[DIM + 1];
392c4762a1bSJed Brown } SWNodeUnion;
393c4762a1bSJed Brown 
3949371c9d4SSatish Balay static const struct FieldDescription PhysicsFields_SW[] = {
3959371c9d4SSatish Balay   {"Height",   1  },
3969371c9d4SSatish Balay   {"Momentum", DIM},
3979371c9d4SSatish Balay   {NULL,       0  }
3989371c9d4SSatish Balay };
399c4762a1bSJed Brown 
400c4762a1bSJed Brown /*
401c4762a1bSJed Brown  * h_t + div(uh) = 0
402c4762a1bSJed Brown  * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0
403c4762a1bSJed Brown  *
404c4762a1bSJed Brown  * */
405d71ae5a4SJacob Faibussowitsch static PetscErrorCode SWFlux(Physics phys, const PetscReal *n, const SWNode *x, SWNode *f)
406d71ae5a4SJacob Faibussowitsch {
407c4762a1bSJed Brown   Physics_SW *sw = (Physics_SW *)phys->data;
408c4762a1bSJed Brown   PetscReal   uhn, u[DIM];
409c4762a1bSJed Brown   PetscInt    i;
410c4762a1bSJed Brown 
411c4762a1bSJed Brown   PetscFunctionBeginUser;
412c4762a1bSJed Brown   Scale2Real(1. / x->h, x->uh, u);
413c4762a1bSJed Brown   uhn  = x->uh[0] * n[0] + x->uh[1] * n[1];
414c4762a1bSJed Brown   f->h = uhn;
415c4762a1bSJed Brown   for (i = 0; i < DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
416c4762a1bSJed Brown   PetscFunctionReturn(0);
417c4762a1bSJed Brown }
418c4762a1bSJed Brown 
419d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
420d71ae5a4SJacob Faibussowitsch {
421c4762a1bSJed Brown   PetscFunctionBeginUser;
422c4762a1bSJed Brown   xG[0] = xI[0];
423c4762a1bSJed Brown   xG[1] = -xI[1];
424c4762a1bSJed Brown   xG[2] = -xI[2];
425c4762a1bSJed Brown   PetscFunctionReturn(0);
426c4762a1bSJed Brown }
427c4762a1bSJed Brown 
428d71ae5a4SJacob Faibussowitsch 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)
429d71ae5a4SJacob Faibussowitsch {
4308d2c18caSMukkund Sunjii   Physics_SW *sw = (Physics_SW *)phys->data;
4318d2c18caSMukkund Sunjii   PetscReal   aL, aR;
4328d2c18caSMukkund Sunjii   PetscReal   nn[DIM];
4338d2c18caSMukkund Sunjii #if !defined(PETSC_USE_COMPLEX)
4348d2c18caSMukkund Sunjii   const SWNode *uL = (const SWNode *)xL, *uR = (const SWNode *)xR;
4358d2c18caSMukkund Sunjii #else
4368d2c18caSMukkund Sunjii   SWNodeUnion   uLreal, uRreal;
4378d2c18caSMukkund Sunjii   const SWNode *uL = &uLreal.swnode;
4388d2c18caSMukkund Sunjii   const SWNode *uR = &uRreal.swnode;
4398d2c18caSMukkund Sunjii #endif
4408d2c18caSMukkund Sunjii   SWNodeUnion fL, fR;
4418d2c18caSMukkund Sunjii   PetscInt    i;
4428d2c18caSMukkund Sunjii   PetscReal   zero = 0.;
4438d2c18caSMukkund Sunjii 
4448d2c18caSMukkund Sunjii #if defined(PETSC_USE_COMPLEX)
4459371c9d4SSatish Balay   uLreal.swnode.h = 0;
4469371c9d4SSatish Balay   uRreal.swnode.h = 0;
4478d2c18caSMukkund Sunjii   for (i = 0; i < 1 + dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
4488d2c18caSMukkund Sunjii   for (i = 0; i < 1 + dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
4498d2c18caSMukkund Sunjii #endif
4508d2c18caSMukkund Sunjii   if (uL->h <= 0 || uR->h <= 0) {
4518d2c18caSMukkund Sunjii     for (i = 0; i < 1 + dim; i++) flux[i] = zero;
4528d2c18caSMukkund Sunjii     return;
4538d2c18caSMukkund Sunjii   } /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */
4548d2c18caSMukkund Sunjii   nn[0] = n[0];
4558d2c18caSMukkund Sunjii   nn[1] = n[1];
4568d2c18caSMukkund Sunjii   Normalize2Real(nn);
4578d2c18caSMukkund Sunjii   SWFlux(phys, nn, uL, &(fL.swnode));
4588d2c18caSMukkund Sunjii   SWFlux(phys, nn, uR, &(fR.swnode));
4598d2c18caSMukkund Sunjii   /* gravity wave speed */
4608d2c18caSMukkund Sunjii   aL = PetscSqrtReal(sw->gravity * uL->h);
4618d2c18caSMukkund Sunjii   aR = PetscSqrtReal(sw->gravity * uR->h);
4628d2c18caSMukkund Sunjii   // Defining u_tilda and v_tilda as u and v
4638d2c18caSMukkund Sunjii   PetscReal u_L, u_R;
4648d2c18caSMukkund Sunjii   u_L = Dot2Real(uL->uh, nn) / uL->h;
4658d2c18caSMukkund Sunjii   u_R = Dot2Real(uR->uh, nn) / uR->h;
4668d2c18caSMukkund Sunjii   PetscReal sL, sR;
4678d2c18caSMukkund Sunjii   sL = PetscMin(u_L - aL, u_R - aR);
4688d2c18caSMukkund Sunjii   sR = PetscMax(u_L + aL, u_R + aR);
4698d2c18caSMukkund Sunjii   if (sL > zero) {
470ad540459SPierre Jolivet     for (i = 0; i < dim + 1; i++) flux[i] = fL.vals[i] * Norm2Real(n);
4718d2c18caSMukkund Sunjii   } else if (sR < zero) {
472ad540459SPierre Jolivet     for (i = 0; i < dim + 1; i++) flux[i] = fR.vals[i] * Norm2Real(n);
4738d2c18caSMukkund Sunjii   } else {
474ad540459SPierre Jolivet     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);
4758d2c18caSMukkund Sunjii   }
4768d2c18caSMukkund Sunjii }
4778d2c18caSMukkund Sunjii 
478d71ae5a4SJacob Faibussowitsch 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)
479d71ae5a4SJacob Faibussowitsch {
480c4762a1bSJed Brown   Physics_SW *sw = (Physics_SW *)phys->data;
481c4762a1bSJed Brown   PetscReal   cL, cR, speed;
482c4762a1bSJed Brown   PetscReal   nn[DIM];
483c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX)
484c4762a1bSJed Brown   const SWNode *uL = (const SWNode *)xL, *uR = (const SWNode *)xR;
485c4762a1bSJed Brown #else
486c4762a1bSJed Brown   SWNodeUnion   uLreal, uRreal;
487c4762a1bSJed Brown   const SWNode *uL = &uLreal.swnode;
488c4762a1bSJed Brown   const SWNode *uR = &uRreal.swnode;
489c4762a1bSJed Brown #endif
490c4762a1bSJed Brown   SWNodeUnion fL, fR;
491c4762a1bSJed Brown   PetscInt    i;
492c4762a1bSJed Brown   PetscReal   zero = 0.;
493c4762a1bSJed Brown 
494c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
4959371c9d4SSatish Balay   uLreal.swnode.h = 0;
4969371c9d4SSatish Balay   uRreal.swnode.h = 0;
497c4762a1bSJed Brown   for (i = 0; i < 1 + dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
498c4762a1bSJed Brown   for (i = 0; i < 1 + dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
499c4762a1bSJed Brown #endif
5009371c9d4SSatish Balay   if (uL->h < 0 || uR->h < 0) {
5019371c9d4SSatish Balay     for (i = 0; i < 1 + dim; i++) flux[i] = zero / zero;
5029371c9d4SSatish Balay     return;
5039371c9d4SSatish Balay   } /* reconstructed thickness is negative */
504c4762a1bSJed Brown   nn[0] = n[0];
505c4762a1bSJed Brown   nn[1] = n[1];
506c4762a1bSJed Brown   Normalize2Real(nn);
507c4762a1bSJed Brown   SWFlux(phys, nn, uL, &(fL.swnode));
508c4762a1bSJed Brown   SWFlux(phys, nn, uR, &(fR.swnode));
509c4762a1bSJed Brown   cL    = PetscSqrtReal(sw->gravity * uL->h);
510c4762a1bSJed Brown   cR    = PetscSqrtReal(sw->gravity * uR->h); /* gravity wave speed */
511c4762a1bSJed Brown   speed = PetscMax(PetscAbsReal(Dot2Real(uL->uh, nn) / uL->h) + cL, PetscAbsReal(Dot2Real(uR->uh, nn) / uR->h) + cR);
512c4762a1bSJed Brown   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);
513c4762a1bSJed Brown }
514c4762a1bSJed Brown 
515d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsSolution_SW(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
516d71ae5a4SJacob Faibussowitsch {
517c4762a1bSJed Brown   PetscReal dx[2], r, sigma;
518c4762a1bSJed Brown 
519c4762a1bSJed Brown   PetscFunctionBeginUser;
52054c59aa7SJacob Faibussowitsch   PetscCheck(time == 0.0, mod->comm, PETSC_ERR_SUP, "No solution known for time %g", (double)time);
521c4762a1bSJed Brown   dx[0] = x[0] - 1.5;
522c4762a1bSJed Brown   dx[1] = x[1] - 1.0;
523c4762a1bSJed Brown   r     = Norm2Real(dx);
524c4762a1bSJed Brown   sigma = 0.5;
525c4762a1bSJed Brown   u[0]  = 1 + 2 * PetscExpReal(-PetscSqr(r) / (2 * PetscSqr(sigma)));
526c4762a1bSJed Brown   u[1]  = 0.0;
527c4762a1bSJed Brown   u[2]  = 0.0;
528c4762a1bSJed Brown   PetscFunctionReturn(0);
529c4762a1bSJed Brown }
530c4762a1bSJed Brown 
531d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsFunctional_SW(Model mod, PetscReal time, const PetscReal *coord, const PetscScalar *xx, PetscReal *f, void *ctx)
532d71ae5a4SJacob Faibussowitsch {
533c4762a1bSJed Brown   Physics       phys = (Physics)ctx;
534c4762a1bSJed Brown   Physics_SW   *sw   = (Physics_SW *)phys->data;
535c4762a1bSJed Brown   const SWNode *x    = (const SWNode *)xx;
536c4762a1bSJed Brown   PetscReal     u[2];
537c4762a1bSJed Brown   PetscReal     h;
538c4762a1bSJed Brown 
539c4762a1bSJed Brown   PetscFunctionBeginUser;
540c4762a1bSJed Brown   h = x->h;
541c4762a1bSJed Brown   Scale2Real(1. / x->h, x->uh, u);
542c4762a1bSJed Brown   f[sw->functional.Height] = h;
543c4762a1bSJed Brown   f[sw->functional.Speed]  = Norm2Real(u) + PetscSqrtReal(sw->gravity * h);
544c4762a1bSJed Brown   f[sw->functional.Energy] = 0.5 * (Dot2Real(x->uh, u) + sw->gravity * PetscSqr(h));
545c4762a1bSJed Brown   PetscFunctionReturn(0);
546c4762a1bSJed Brown }
547c4762a1bSJed Brown 
548d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetUpBC_SW(DM dm, PetscDS prob, Physics phys)
549d71ae5a4SJacob Faibussowitsch {
550c4762a1bSJed Brown   const PetscInt wallids[] = {100, 101, 200, 300};
55145480ffeSMatthew G. Knepley   DMLabel        label;
55245480ffeSMatthew G. Knepley 
553c4762a1bSJed Brown   PetscFunctionBeginUser;
5549566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "Face Sets", &label));
555dd39110bSPierre Jolivet   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));
556c4762a1bSJed Brown   PetscFunctionReturn(0);
557c4762a1bSJed Brown }
558c4762a1bSJed Brown 
559d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCreate_SW(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject)
560d71ae5a4SJacob Faibussowitsch {
561c4762a1bSJed Brown   Physics_SW *sw;
5628d2c18caSMukkund Sunjii   char        sw_riemann[64] = "rusanov";
563c4762a1bSJed Brown 
564c4762a1bSJed Brown   PetscFunctionBeginUser;
565c4762a1bSJed Brown   phys->field_desc = PhysicsFields_SW;
5669566063dSJacob Faibussowitsch   PetscCall(PetscNew(&sw));
567c4762a1bSJed Brown   phys->data   = sw;
568c4762a1bSJed Brown   mod->setupbc = SetUpBC_SW;
569c4762a1bSJed Brown 
5708d2c18caSMukkund Sunjii   PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov", PhysicsRiemann_SW_Rusanov);
5718d2c18caSMukkund Sunjii   PetscFunctionListAdd(&PhysicsRiemannList_SW, "hll", PhysicsRiemann_SW_HLL);
5728d2c18caSMukkund Sunjii 
573d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SW options");
574c4762a1bSJed Brown   {
5758d2c18caSMukkund Sunjii     void (*PhysicsRiemann_SW)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics);
576c4762a1bSJed Brown     sw->gravity = 1.0;
5779566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-sw_gravity", "Gravitational constant", "", sw->gravity, &sw->gravity, NULL));
5789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-sw_riemann", "Riemann solver", "", PhysicsRiemannList_SW, sw_riemann, sw_riemann, sizeof sw_riemann, NULL));
5799566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(PhysicsRiemannList_SW, sw_riemann, &PhysicsRiemann_SW));
5808d2c18caSMukkund Sunjii     phys->riemann = (PetscRiemannFunc)PhysicsRiemann_SW;
581c4762a1bSJed Brown   }
582d0609cedSBarry Smith   PetscOptionsHeadEnd();
583c4762a1bSJed Brown   phys->maxspeed = PetscSqrtReal(2.0 * sw->gravity); /* Mach 1 for depth of 2 */
584c4762a1bSJed Brown 
5859566063dSJacob Faibussowitsch   PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_SW, phys));
5869566063dSJacob Faibussowitsch   PetscCall(ModelFunctionalRegister(mod, "Height", &sw->functional.Height, PhysicsFunctional_SW, phys));
5879566063dSJacob Faibussowitsch   PetscCall(ModelFunctionalRegister(mod, "Speed", &sw->functional.Speed, PhysicsFunctional_SW, phys));
5889566063dSJacob Faibussowitsch   PetscCall(ModelFunctionalRegister(mod, "Energy", &sw->functional.Energy, PhysicsFunctional_SW, phys));
589c4762a1bSJed Brown 
590c4762a1bSJed Brown   PetscFunctionReturn(0);
591c4762a1bSJed Brown }
592c4762a1bSJed Brown 
593c4762a1bSJed Brown /******************* Euler Density Shock (EULER_IV_SHOCK,EULER_SS_SHOCK) ********************/
594c4762a1bSJed Brown /* An initial-value and self-similar solutions of the compressible Euler equations */
595c4762a1bSJed Brown /* Ravi Samtaney and D. I. Pullin */
596c4762a1bSJed Brown /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
5979371c9d4SSatish Balay typedef enum {
5989371c9d4SSatish Balay   EULER_PAR_GAMMA,
5999371c9d4SSatish Balay   EULER_PAR_RHOR,
6009371c9d4SSatish Balay   EULER_PAR_AMACH,
6019371c9d4SSatish Balay   EULER_PAR_ITANA,
6029371c9d4SSatish Balay   EULER_PAR_SIZE
6039371c9d4SSatish Balay } EulerParamIdx;
6049371c9d4SSatish Balay typedef enum {
6059371c9d4SSatish Balay   EULER_IV_SHOCK,
6069371c9d4SSatish Balay   EULER_SS_SHOCK,
6079371c9d4SSatish Balay   EULER_SHOCK_TUBE,
6089371c9d4SSatish Balay   EULER_LINEAR_WAVE
6099371c9d4SSatish Balay } EulerType;
610c4762a1bSJed Brown typedef struct {
611c4762a1bSJed Brown   PetscReal r;
612c4762a1bSJed Brown   PetscReal ru[DIM];
613c4762a1bSJed Brown   PetscReal E;
614c4762a1bSJed Brown } EulerNode;
6159371c9d4SSatish Balay typedef union
6169371c9d4SSatish Balay {
617c4762a1bSJed Brown   EulerNode eulernode;
618c4762a1bSJed Brown   PetscReal vals[DIM + 2];
619c4762a1bSJed Brown } EulerNodeUnion;
620c4762a1bSJed Brown typedef PetscErrorCode (*EquationOfState)(const PetscReal *, const EulerNode *, PetscReal *);
621c4762a1bSJed Brown typedef struct {
622c4762a1bSJed Brown   EulerType       type;
623c4762a1bSJed Brown   PetscReal       pars[EULER_PAR_SIZE];
624c4762a1bSJed Brown   EquationOfState sound;
625c4762a1bSJed Brown   struct {
626c4762a1bSJed Brown     PetscInt Density;
627c4762a1bSJed Brown     PetscInt Momentum;
628c4762a1bSJed Brown     PetscInt Energy;
629c4762a1bSJed Brown     PetscInt Pressure;
630c4762a1bSJed Brown     PetscInt Speed;
631c4762a1bSJed Brown   } monitor;
632c4762a1bSJed Brown } Physics_Euler;
633c4762a1bSJed Brown 
6349371c9d4SSatish Balay static const struct FieldDescription PhysicsFields_Euler[] = {
6359371c9d4SSatish Balay   {"Density",  1  },
6369371c9d4SSatish Balay   {"Momentum", DIM},
6379371c9d4SSatish Balay   {"Energy",   1  },
6389371c9d4SSatish Balay   {NULL,       0  }
6399371c9d4SSatish Balay };
640c4762a1bSJed Brown 
641c4762a1bSJed Brown /* initial condition */
642c4762a1bSJed Brown int                   initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx);
643d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
644d71ae5a4SJacob Faibussowitsch {
645c4762a1bSJed Brown   PetscInt       i;
646c4762a1bSJed Brown   Physics        phys = (Physics)ctx;
647c4762a1bSJed Brown   Physics_Euler *eu   = (Physics_Euler *)phys->data;
648c4762a1bSJed Brown   EulerNode     *uu   = (EulerNode *)u;
649c4762a1bSJed Brown   PetscReal      p0, gamma, c;
650c4762a1bSJed Brown   PetscFunctionBeginUser;
65154c59aa7SJacob Faibussowitsch   PetscCheck(time == 0.0, mod->comm, PETSC_ERR_SUP, "No solution known for time %g", (double)time);
652c4762a1bSJed Brown 
653c4762a1bSJed Brown   for (i = 0; i < DIM; i++) uu->ru[i] = 0.0; /* zero out initial velocity */
654c4762a1bSJed Brown   /* set E and rho */
655c4762a1bSJed Brown   gamma = eu->pars[EULER_PAR_GAMMA];
656c4762a1bSJed Brown 
657c4762a1bSJed Brown   if (eu->type == EULER_IV_SHOCK || eu->type == EULER_SS_SHOCK) {
658c4762a1bSJed Brown     /******************* Euler Density Shock ********************/
659c4762a1bSJed Brown     /* On initial-value and self-similar solutions of the compressible Euler equations */
660c4762a1bSJed Brown     /* Ravi Samtaney and D. I. Pullin */
661c4762a1bSJed Brown     /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
662c4762a1bSJed Brown     /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity,  */
663c4762a1bSJed Brown     p0 = 1.;
664c4762a1bSJed Brown     if (x[0] < 0.0 + x[1] * eu->pars[EULER_PAR_ITANA]) {
665c4762a1bSJed Brown       if (x[0] < mod->bounds[0] * 0.5) { /* left of shock (1) */
666c4762a1bSJed Brown         PetscReal amach, rho, press, gas1, p1;
667c4762a1bSJed Brown         amach     = eu->pars[EULER_PAR_AMACH];
668c4762a1bSJed Brown         rho       = 1.;
669c4762a1bSJed Brown         press     = p0;
670c4762a1bSJed Brown         p1        = press * (1.0 + 2.0 * gamma / (gamma + 1.0) * (amach * amach - 1.0));
671c4762a1bSJed Brown         gas1      = (gamma - 1.0) / (gamma + 1.0);
672c4762a1bSJed Brown         uu->r     = rho * (p1 / press + gas1) / (gas1 * p1 / press + 1.0);
673c4762a1bSJed Brown         uu->ru[0] = ((uu->r - rho) * PetscSqrtReal(gamma * press / rho) * amach);
674c4762a1bSJed Brown         uu->E     = p1 / (gamma - 1.0) + .5 / uu->r * uu->ru[0] * uu->ru[0];
6759371c9d4SSatish Balay       } else {      /* left of discontinuity (0) */
676c4762a1bSJed Brown         uu->r = 1.; /* rho = 1 */
677c4762a1bSJed Brown         uu->E = p0 / (gamma - 1.0);
678c4762a1bSJed Brown       }
6799371c9d4SSatish Balay     } else { /* right of discontinuity (2) */
680c4762a1bSJed Brown       uu->r = eu->pars[EULER_PAR_RHOR];
681c4762a1bSJed Brown       uu->E = p0 / (gamma - 1.0);
682c4762a1bSJed Brown     }
6839371c9d4SSatish Balay   } else if (eu->type == EULER_SHOCK_TUBE) {
684c4762a1bSJed Brown     /* 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. */
685c4762a1bSJed Brown     if (x[0] < 0.0) {
686c4762a1bSJed Brown       uu->r = 8.;
687c4762a1bSJed Brown       uu->E = 10. / (gamma - 1.);
6889371c9d4SSatish Balay     } else {
689c4762a1bSJed Brown       uu->r = 1.;
690c4762a1bSJed Brown       uu->E = 1. / (gamma - 1.);
691c4762a1bSJed Brown     }
6929371c9d4SSatish Balay   } else if (eu->type == EULER_LINEAR_WAVE) {
693c4762a1bSJed Brown     initLinearWave(uu, gamma, x, mod->bounds[1] - mod->bounds[0]);
6949371c9d4SSatish Balay   } else SETERRQ(mod->comm, PETSC_ERR_SUP, "Unknown type %d", eu->type);
695c4762a1bSJed Brown 
696c4762a1bSJed Brown   /* set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main; */
697c4762a1bSJed Brown   eu->sound(&gamma, uu, &c);
698c4762a1bSJed Brown   c = (uu->ru[0] / uu->r) + c;
699c4762a1bSJed Brown   if (c > phys->maxspeed) phys->maxspeed = c;
700c4762a1bSJed Brown 
701c4762a1bSJed Brown   PetscFunctionReturn(0);
702c4762a1bSJed Brown }
703c4762a1bSJed Brown 
704d71ae5a4SJacob Faibussowitsch static PetscErrorCode Pressure_PG(const PetscReal gamma, const EulerNode *x, PetscReal *p)
705d71ae5a4SJacob Faibussowitsch {
706c4762a1bSJed Brown   PetscReal ru2;
707c4762a1bSJed Brown 
708c4762a1bSJed Brown   PetscFunctionBeginUser;
709c4762a1bSJed Brown   ru2  = DotDIMReal(x->ru, x->ru);
710c4762a1bSJed Brown   (*p) = (x->E - 0.5 * ru2 / x->r) * (gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */
711c4762a1bSJed Brown   PetscFunctionReturn(0);
712c4762a1bSJed Brown }
713c4762a1bSJed Brown 
714d71ae5a4SJacob Faibussowitsch static PetscErrorCode SpeedOfSound_PG(const PetscReal *gamma, const EulerNode *x, PetscReal *c)
715d71ae5a4SJacob Faibussowitsch {
716c4762a1bSJed Brown   PetscReal p;
717c4762a1bSJed Brown 
718c4762a1bSJed Brown   PetscFunctionBeginUser;
719c4762a1bSJed Brown   Pressure_PG(*gamma, x, &p);
72054c59aa7SJacob Faibussowitsch   PetscCheck(p >= 0., PETSC_COMM_WORLD, PETSC_ERR_SUP, "negative pressure time %g -- NEED TO FIX!!!!!!", (double)p);
721c4762a1bSJed Brown   /* pars[EULER_PAR_GAMMA] = heat capacity ratio */
722c4762a1bSJed Brown   (*c) = PetscSqrtReal(*gamma * p / x->r);
723c4762a1bSJed Brown   PetscFunctionReturn(0);
724c4762a1bSJed Brown }
725c4762a1bSJed Brown 
726c4762a1bSJed Brown /*
727c4762a1bSJed Brown  * x = (rho,rho*(u_1),...,rho*e)^T
728c4762a1bSJed Brown  * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
729c4762a1bSJed Brown  *
730c4762a1bSJed Brown  * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
731c4762a1bSJed Brown  *
732c4762a1bSJed Brown  */
733d71ae5a4SJacob Faibussowitsch static PetscErrorCode EulerFlux(Physics phys, const PetscReal *n, const EulerNode *x, EulerNode *f)
734d71ae5a4SJacob Faibussowitsch {
735c4762a1bSJed Brown   Physics_Euler *eu = (Physics_Euler *)phys->data;
736c4762a1bSJed Brown   PetscReal      nu, p;
737c4762a1bSJed Brown   PetscInt       i;
738c4762a1bSJed Brown 
739c4762a1bSJed Brown   PetscFunctionBeginUser;
740c4762a1bSJed Brown   Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p);
741c4762a1bSJed Brown   nu   = DotDIMReal(x->ru, n);
742c4762a1bSJed Brown   f->r = nu;                                                     /* A rho u */
743c4762a1bSJed Brown   nu /= x->r;                                                    /* A u */
744c4762a1bSJed Brown   for (i = 0; i < DIM; i++) f->ru[i] = nu * x->ru[i] + n[i] * p; /* r u^2 + p */
745c4762a1bSJed Brown   f->E = nu * (x->E + p);                                        /* u(e+p) */
746c4762a1bSJed Brown   PetscFunctionReturn(0);
747c4762a1bSJed Brown }
748c4762a1bSJed Brown 
749c4762a1bSJed Brown /* PetscReal* => EulerNode* conversion */
750d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, PetscScalar *a_xG, void *ctx)
751d71ae5a4SJacob Faibussowitsch {
752c4762a1bSJed Brown   PetscInt         i;
753c4762a1bSJed Brown   const EulerNode *xI   = (const EulerNode *)a_xI;
754c4762a1bSJed Brown   EulerNode       *xG   = (EulerNode *)a_xG;
755c4762a1bSJed Brown   Physics          phys = (Physics)ctx;
756c4762a1bSJed Brown   Physics_Euler   *eu   = (Physics_Euler *)phys->data;
757c4762a1bSJed Brown   PetscFunctionBeginUser;
758c4762a1bSJed Brown   xG->r = xI->r;                                     /* ghost cell density - same */
759c4762a1bSJed Brown   xG->E = xI->E;                                     /* ghost cell energy - same */
760c4762a1bSJed Brown   if (n[1] != 0.) {                                  /* top and bottom */
761c4762a1bSJed Brown     xG->ru[0] = xI->ru[0];                           /* copy tang to wall */
762c4762a1bSJed Brown     xG->ru[1] = -xI->ru[1];                          /* reflect perp to t/b wall */
7639371c9d4SSatish Balay   } else {                                           /* sides */
764c4762a1bSJed Brown     for (i = 0; i < DIM; i++) xG->ru[i] = xI->ru[i]; /* copy */
765c4762a1bSJed Brown   }
766c4762a1bSJed Brown   if (eu->type == EULER_LINEAR_WAVE) { /* debug */
767c4762a1bSJed Brown #if 0
76863a3b9bcSJacob Faibussowitsch     PetscPrintf(PETSC_COMM_WORLD,"%s coord=%g,%g\n",PETSC_FUNCTION_NAME,(double)c[0],(double)c[1]);
769c4762a1bSJed Brown #endif
770c4762a1bSJed Brown   }
771c4762a1bSJed Brown   PetscFunctionReturn(0);
772c4762a1bSJed Brown }
773c4762a1bSJed Brown int godunovflux(const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, const int *ndim, const PetscReal *gamma);
774c4762a1bSJed Brown /* PetscReal* => EulerNode* conversion */
775d71ae5a4SJacob Faibussowitsch 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)
776d71ae5a4SJacob Faibussowitsch {
777c4762a1bSJed Brown   Physics_Euler *eu = (Physics_Euler *)phys->data;
778c4762a1bSJed Brown   PetscReal      cL, cR, speed, velL, velR, nn[DIM], s2;
779c4762a1bSJed Brown   PetscInt       i;
780c4762a1bSJed Brown   PetscErrorCode ierr;
781c4762a1bSJed Brown 
782d0609cedSBarry Smith   PetscFunctionBeginUser;
783c4762a1bSJed Brown   for (i = 0, s2 = 0.; i < DIM; i++) {
784c4762a1bSJed Brown     nn[i] = n[i];
785c4762a1bSJed Brown     s2 += nn[i] * nn[i];
786c4762a1bSJed Brown   }
787c4762a1bSJed Brown   s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */
788c4762a1bSJed Brown   for (i = 0.; i < DIM; i++) nn[i] /= s2;
789c4762a1bSJed Brown   if (0) { /* Rusanov */
790c4762a1bSJed Brown     const EulerNode *uL = (const EulerNode *)xL, *uR = (const EulerNode *)xR;
791c4762a1bSJed Brown     EulerNodeUnion   fL, fR;
792c4762a1bSJed Brown     EulerFlux(phys, nn, uL, &(fL.eulernode));
793c4762a1bSJed Brown     EulerFlux(phys, nn, uR, &(fR.eulernode));
7949371c9d4SSatish Balay     ierr = eu->sound(&eu->pars[EULER_PAR_GAMMA], uL, &cL);
7959371c9d4SSatish Balay     if (ierr) exit(13);
7969371c9d4SSatish Balay     ierr = eu->sound(&eu->pars[EULER_PAR_GAMMA], uR, &cR);
7979371c9d4SSatish Balay     if (ierr) exit(14);
798c4762a1bSJed Brown     velL  = DotDIMReal(uL->ru, nn) / uL->r;
799c4762a1bSJed Brown     velR  = DotDIMReal(uR->ru, nn) / uR->r;
800c4762a1bSJed Brown     speed = PetscMax(velR + cR, velL + cL);
801c4762a1bSJed Brown     for (i = 0; i < 2 + dim; i++) flux[i] = 0.5 * ((fL.vals[i] + fR.vals[i]) + speed * (xL[i] - xR[i])) * s2;
8029371c9d4SSatish Balay   } else {
803c4762a1bSJed Brown     int dim = DIM;
804c4762a1bSJed Brown     /* int iwave =  */
805c4762a1bSJed Brown     godunovflux(xL, xR, flux, nn, &dim, &eu->pars[EULER_PAR_GAMMA]);
806c4762a1bSJed Brown     for (i = 0; i < 2 + dim; i++) flux[i] *= s2;
807c4762a1bSJed Brown   }
808c4762a1bSJed Brown   PetscFunctionReturnVoid();
809c4762a1bSJed Brown }
810c4762a1bSJed Brown 
811d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsFunctional_Euler(Model mod, PetscReal time, const PetscReal *coord, const PetscScalar *xx, PetscReal *f, void *ctx)
812d71ae5a4SJacob Faibussowitsch {
813c4762a1bSJed Brown   Physics          phys = (Physics)ctx;
814c4762a1bSJed Brown   Physics_Euler   *eu   = (Physics_Euler *)phys->data;
815c4762a1bSJed Brown   const EulerNode *x    = (const EulerNode *)xx;
816c4762a1bSJed Brown   PetscReal        p;
817c4762a1bSJed Brown 
818c4762a1bSJed Brown   PetscFunctionBeginUser;
819c4762a1bSJed Brown   f[eu->monitor.Density]  = x->r;
820c4762a1bSJed Brown   f[eu->monitor.Momentum] = NormDIM(x->ru);
821c4762a1bSJed Brown   f[eu->monitor.Energy]   = x->E;
822c4762a1bSJed Brown   f[eu->monitor.Speed]    = NormDIM(x->ru) / x->r;
823c4762a1bSJed Brown   Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p);
824c4762a1bSJed Brown   f[eu->monitor.Pressure] = p;
825c4762a1bSJed Brown   PetscFunctionReturn(0);
826c4762a1bSJed Brown }
827c4762a1bSJed Brown 
828d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetUpBC_Euler(DM dm, PetscDS prob, Physics phys)
829d71ae5a4SJacob Faibussowitsch {
830c4762a1bSJed Brown   Physics_Euler *eu = (Physics_Euler *)phys->data;
83145480ffeSMatthew G. Knepley   DMLabel        label;
83245480ffeSMatthew G. Knepley 
833362febeeSStefano Zampini   PetscFunctionBeginUser;
8349566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "Face Sets", &label));
835c4762a1bSJed Brown   if (eu->type == EULER_LINEAR_WAVE) {
836c4762a1bSJed Brown     const PetscInt wallids[] = {100, 101};
837dd39110bSPierre Jolivet     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));
8389371c9d4SSatish Balay   } else {
839c4762a1bSJed Brown     const PetscInt wallids[] = {100, 101, 200, 300};
840dd39110bSPierre Jolivet     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));
841c4762a1bSJed Brown   }
842c4762a1bSJed Brown   PetscFunctionReturn(0);
843c4762a1bSJed Brown }
844c4762a1bSJed Brown 
845d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCreate_Euler(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject)
846d71ae5a4SJacob Faibussowitsch {
847c4762a1bSJed Brown   Physics_Euler *eu;
848c4762a1bSJed Brown 
849c4762a1bSJed Brown   PetscFunctionBeginUser;
850c4762a1bSJed Brown   phys->field_desc = PhysicsFields_Euler;
851c4762a1bSJed Brown   phys->riemann    = (PetscRiemannFunc)PhysicsRiemann_Euler_Godunov;
8529566063dSJacob Faibussowitsch   PetscCall(PetscNew(&eu));
853c4762a1bSJed Brown   phys->data   = eu;
854c4762a1bSJed Brown   mod->setupbc = SetUpBC_Euler;
855d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Euler options");
856c4762a1bSJed Brown   {
857c4762a1bSJed Brown     PetscReal alpha;
858c4762a1bSJed Brown     char      type[64] = "linear_wave";
859c4762a1bSJed Brown     PetscBool is;
860c4762a1bSJed Brown     eu->pars[EULER_PAR_GAMMA] = 1.4;
861c4762a1bSJed Brown     eu->pars[EULER_PAR_AMACH] = 2.02;
862c4762a1bSJed Brown     eu->pars[EULER_PAR_RHOR]  = 3.0;
863c4762a1bSJed Brown     eu->pars[EULER_PAR_ITANA] = 0.57735026918963; /* angle of Euler self similar (SS) shock */
8649566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-eu_gamma", "Heat capacity ratio", "", eu->pars[EULER_PAR_GAMMA], &eu->pars[EULER_PAR_GAMMA], NULL));
8659566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-eu_amach", "Shock speed (Mach)", "", eu->pars[EULER_PAR_AMACH], &eu->pars[EULER_PAR_AMACH], NULL));
8669566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-eu_rho2", "Density right of discontinuity", "", eu->pars[EULER_PAR_RHOR], &eu->pars[EULER_PAR_RHOR], NULL));
867c4762a1bSJed Brown     alpha = 60.;
8689566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-eu_alpha", "Angle of discontinuity", "", alpha, &alpha, NULL));
86963a3b9bcSJacob Faibussowitsch     PetscCheck(alpha > 0. && alpha <= 90., PETSC_COMM_WORLD, PETSC_ERR_SUP, "Alpha bust be > 0 and <= 90 (%g)", (double)alpha);
870c4762a1bSJed Brown     eu->pars[EULER_PAR_ITANA] = 1. / PetscTanReal(alpha * PETSC_PI / 180.0);
8719566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-eu_type", "Type of Euler test", "", type, type, sizeof(type), NULL));
8729566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(type, "linear_wave", &is));
873c4762a1bSJed Brown     if (is) {
87430602db0SMatthew G. Knepley       /* Remember this should be periodic */
875c4762a1bSJed Brown       eu->type = EULER_LINEAR_WAVE;
8769566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "linear_wave"));
8779371c9d4SSatish Balay     } else {
87854c59aa7SJacob Faibussowitsch       PetscCheck(DIM == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "DIM must be 2 unless linear wave test %s", type);
8799566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(type, "iv_shock", &is));
880c4762a1bSJed Brown       if (is) {
881c4762a1bSJed Brown         eu->type = EULER_IV_SHOCK;
8829566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "iv_shock"));
8839371c9d4SSatish Balay       } else {
8849566063dSJacob Faibussowitsch         PetscCall(PetscStrcmp(type, "ss_shock", &is));
885c4762a1bSJed Brown         if (is) {
886c4762a1bSJed Brown           eu->type = EULER_SS_SHOCK;
8879566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "ss_shock"));
8889371c9d4SSatish Balay         } else {
8899566063dSJacob Faibussowitsch           PetscCall(PetscStrcmp(type, "shock_tube", &is));
890c4762a1bSJed Brown           if (is) eu->type = EULER_SHOCK_TUBE;
89198921bdaSJacob Faibussowitsch           else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unknown Euler type %s", type);
8929566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "shock_tube"));
893c4762a1bSJed Brown         }
894c4762a1bSJed Brown       }
895c4762a1bSJed Brown     }
896c4762a1bSJed Brown   }
897d0609cedSBarry Smith   PetscOptionsHeadEnd();
898c4762a1bSJed Brown   eu->sound      = SpeedOfSound_PG;
899c4762a1bSJed Brown   phys->maxspeed = 0.; /* will get set in solution */
9009566063dSJacob Faibussowitsch   PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_Euler, phys));
9019566063dSJacob Faibussowitsch   PetscCall(ModelFunctionalRegister(mod, "Speed", &eu->monitor.Speed, PhysicsFunctional_Euler, phys));
9029566063dSJacob Faibussowitsch   PetscCall(ModelFunctionalRegister(mod, "Energy", &eu->monitor.Energy, PhysicsFunctional_Euler, phys));
9039566063dSJacob Faibussowitsch   PetscCall(ModelFunctionalRegister(mod, "Density", &eu->monitor.Density, PhysicsFunctional_Euler, phys));
9049566063dSJacob Faibussowitsch   PetscCall(ModelFunctionalRegister(mod, "Momentum", &eu->monitor.Momentum, PhysicsFunctional_Euler, phys));
9059566063dSJacob Faibussowitsch   PetscCall(ModelFunctionalRegister(mod, "Pressure", &eu->monitor.Pressure, PhysicsFunctional_Euler, phys));
906c4762a1bSJed Brown 
907c4762a1bSJed Brown   PetscFunctionReturn(0);
908c4762a1bSJed Brown }
909c4762a1bSJed Brown 
910d71ae5a4SJacob Faibussowitsch static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx)
911d71ae5a4SJacob Faibussowitsch {
912c4762a1bSJed Brown   PetscReal err = 0.;
913c4762a1bSJed Brown   PetscInt  i, j;
914c4762a1bSJed Brown 
915c4762a1bSJed Brown   PetscFunctionBeginUser;
916c4762a1bSJed Brown   for (i = 0; i < numComps; i++) {
917ad540459SPierre Jolivet     for (j = 0; j < dim; j++) err += PetscSqr(PetscRealPart(grad[i * dim + j]));
918c4762a1bSJed Brown   }
919c4762a1bSJed Brown   *error = volume * err;
920c4762a1bSJed Brown   PetscFunctionReturn(0);
921c4762a1bSJed Brown }
922c4762a1bSJed Brown 
923d71ae5a4SJacob Faibussowitsch PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
924d71ae5a4SJacob Faibussowitsch {
925c4762a1bSJed Brown   PetscSF      sfPoint;
926c4762a1bSJed Brown   PetscSection coordSection;
927c4762a1bSJed Brown   Vec          coordinates;
928c4762a1bSJed Brown   PetscSection sectionCell;
929c4762a1bSJed Brown   PetscScalar *part;
930c4762a1bSJed Brown   PetscInt     cStart, cEnd, c;
931c4762a1bSJed Brown   PetscMPIInt  rank;
932c4762a1bSJed Brown 
933c4762a1bSJed Brown   PetscFunctionBeginUser;
9349566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordSection));
9359566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
9369566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, dmCell));
9379566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sfPoint));
9389566063dSJacob Faibussowitsch   PetscCall(DMSetPointSF(*dmCell, sfPoint));
9399566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection));
9409566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(*dmCell, coordinates));
9419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
9429566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionCell));
9439566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd));
9449566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(sectionCell, cStart, cEnd));
94548a46eb9SPierre Jolivet   for (c = cStart; c < cEnd; ++c) PetscCall(PetscSectionSetDof(sectionCell, c, 1));
9469566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(sectionCell));
9479566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(*dmCell, sectionCell));
9489566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&sectionCell));
9499566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(*dmCell, partition));
9509566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*partition, "partition"));
9519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(*partition, &part));
952c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
953c4762a1bSJed Brown     PetscScalar *p;
954c4762a1bSJed Brown 
9559566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRef(*dmCell, c, part, &p));
956c4762a1bSJed Brown     p[0] = rank;
957c4762a1bSJed Brown   }
9589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(*partition, &part));
959c4762a1bSJed Brown   PetscFunctionReturn(0);
960c4762a1bSJed Brown }
961c4762a1bSJed Brown 
962d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
963d71ae5a4SJacob Faibussowitsch {
9643e9753d6SMatthew G. Knepley   DM                 plex, dmMass, dmFace, dmCell, dmCoord;
965c4762a1bSJed Brown   PetscSection       coordSection;
966c4762a1bSJed Brown   Vec                coordinates, facegeom, cellgeom;
967c4762a1bSJed Brown   PetscSection       sectionMass;
968c4762a1bSJed Brown   PetscScalar       *m;
969c4762a1bSJed Brown   const PetscScalar *fgeom, *cgeom, *coords;
970c4762a1bSJed Brown   PetscInt           vStart, vEnd, v;
971c4762a1bSJed Brown 
972c4762a1bSJed Brown   PetscFunctionBeginUser;
9739566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
9749566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordSection));
9759566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
9769566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmMass));
9779566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection));
9789566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(dmMass, coordinates));
9799566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionMass));
9809566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
9819566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(sectionMass, vStart, vEnd));
982c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
983c4762a1bSJed Brown     PetscInt numFaces;
984c4762a1bSJed Brown 
9859566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(dmMass, v, &numFaces));
9869566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(sectionMass, v, numFaces * numFaces));
987c4762a1bSJed Brown   }
9889566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(sectionMass));
9899566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dmMass, sectionMass));
9909566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&sectionMass));
9919566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dmMass, massMatrix));
9929566063dSJacob Faibussowitsch   PetscCall(VecGetArray(*massMatrix, &m));
9939566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGeometryFVM(plex, &facegeom, &cellgeom, NULL));
9949566063dSJacob Faibussowitsch   PetscCall(VecGetDM(facegeom, &dmFace));
9959566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(facegeom, &fgeom));
9969566063dSJacob Faibussowitsch   PetscCall(VecGetDM(cellgeom, &dmCell));
9979566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(cellgeom, &cgeom));
9989566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
9999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coordinates, &coords));
1000c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
1001c4762a1bSJed Brown     const PetscInt  *faces;
1002c4762a1bSJed Brown     PetscFVFaceGeom *fgA, *fgB, *cg;
1003c4762a1bSJed Brown     PetscScalar     *vertex;
1004c4762a1bSJed Brown     PetscInt         numFaces, sides[2], f, g;
1005c4762a1bSJed Brown 
10069566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRead(dmCoord, v, coords, &vertex));
10079566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(dmMass, v, &numFaces));
10089566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(dmMass, v, &faces));
1009c4762a1bSJed Brown     for (f = 0; f < numFaces; ++f) {
1010c4762a1bSJed Brown       sides[0] = faces[f];
10119566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA));
1012c4762a1bSJed Brown       for (g = 0; g < numFaces; ++g) {
1013c4762a1bSJed Brown         const PetscInt *cells = NULL;
1014c4762a1bSJed Brown         PetscReal       area  = 0.0;
1015c4762a1bSJed Brown         PetscInt        numCells;
1016c4762a1bSJed Brown 
1017c4762a1bSJed Brown         sides[1] = faces[g];
10189566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB));
10199566063dSJacob Faibussowitsch         PetscCall(DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells));
102054c59aa7SJacob Faibussowitsch         PetscCheck(numCells == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces");
10219566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg));
1022c4762a1bSJed Brown         area += PetscAbsScalar((vertex[0] - cg->centroid[0]) * (fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1]) * (fgA->centroid[0] - cg->centroid[0]));
1023c4762a1bSJed Brown         area += PetscAbsScalar((vertex[0] - cg->centroid[0]) * (fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1]) * (fgB->centroid[0] - cg->centroid[0]));
1024c4762a1bSJed Brown         m[f * numFaces + g] = Dot2Real(fgA->normal, fgB->normal) * area * 0.5;
10259566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells));
1026c4762a1bSJed Brown       }
1027c4762a1bSJed Brown     }
1028c4762a1bSJed Brown   }
10299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(facegeom, &fgeom));
10309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
10319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coordinates, &coords));
10329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(*massMatrix, &m));
10339566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmMass));
10349566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1035c4762a1bSJed Brown   PetscFunctionReturn(0);
1036c4762a1bSJed Brown }
1037c4762a1bSJed Brown 
1038c4762a1bSJed Brown /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1039d71ae5a4SJacob Faibussowitsch static PetscErrorCode ModelSolutionSetDefault(Model mod, SolutionFunction func, void *ctx)
1040d71ae5a4SJacob Faibussowitsch {
1041c4762a1bSJed Brown   PetscFunctionBeginUser;
1042c4762a1bSJed Brown   mod->solution    = func;
1043c4762a1bSJed Brown   mod->solutionctx = ctx;
1044c4762a1bSJed Brown   PetscFunctionReturn(0);
1045c4762a1bSJed Brown }
1046c4762a1bSJed Brown 
1047d71ae5a4SJacob Faibussowitsch static PetscErrorCode ModelFunctionalRegister(Model mod, const char *name, PetscInt *offset, FunctionalFunction func, void *ctx)
1048d71ae5a4SJacob Faibussowitsch {
1049c4762a1bSJed Brown   FunctionalLink link, *ptr;
1050c4762a1bSJed Brown   PetscInt       lastoffset = -1;
1051c4762a1bSJed Brown 
1052c4762a1bSJed Brown   PetscFunctionBeginUser;
1053c4762a1bSJed Brown   for (ptr = &mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
10549566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
10559566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &link->name));
1056c4762a1bSJed Brown   link->offset = lastoffset + 1;
1057c4762a1bSJed Brown   link->func   = func;
1058c4762a1bSJed Brown   link->ctx    = ctx;
1059c4762a1bSJed Brown   link->next   = NULL;
1060c4762a1bSJed Brown   *ptr         = link;
1061c4762a1bSJed Brown   *offset      = link->offset;
1062c4762a1bSJed Brown   PetscFunctionReturn(0);
1063c4762a1bSJed Brown }
1064c4762a1bSJed Brown 
1065d71ae5a4SJacob Faibussowitsch static PetscErrorCode ModelFunctionalSetFromOptions(Model mod, PetscOptionItems *PetscOptionsObject)
1066d71ae5a4SJacob Faibussowitsch {
1067c4762a1bSJed Brown   PetscInt       i, j;
1068c4762a1bSJed Brown   FunctionalLink link;
1069c4762a1bSJed Brown   char          *names[256];
1070c4762a1bSJed Brown 
1071c4762a1bSJed Brown   PetscFunctionBeginUser;
1072dd39110bSPierre Jolivet   mod->numMonitored = PETSC_STATIC_ARRAY_LENGTH(names);
10739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsStringArray("-monitor", "list of functionals to monitor", "", names, &mod->numMonitored, NULL));
1074c4762a1bSJed Brown   /* Create list of functionals that will be computed somehow */
10759566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mod->numMonitored, &mod->functionalMonitored));
1076c4762a1bSJed Brown   /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
10779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mod->numMonitored, &mod->functionalCall));
1078c4762a1bSJed Brown   mod->numCall = 0;
1079c4762a1bSJed Brown   for (i = 0; i < mod->numMonitored; i++) {
1080c4762a1bSJed Brown     for (link = mod->functionalRegistry; link; link = link->next) {
1081c4762a1bSJed Brown       PetscBool match;
10829566063dSJacob Faibussowitsch       PetscCall(PetscStrcasecmp(names[i], link->name, &match));
1083c4762a1bSJed Brown       if (match) break;
1084c4762a1bSJed Brown     }
108554c59aa7SJacob Faibussowitsch     PetscCheck(link, mod->comm, PETSC_ERR_USER, "No known functional '%s'", names[i]);
1086c4762a1bSJed Brown     mod->functionalMonitored[i] = link;
1087c4762a1bSJed Brown     for (j = 0; j < i; j++) {
1088c4762a1bSJed Brown       if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1089c4762a1bSJed Brown     }
1090c4762a1bSJed Brown     mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1091c4762a1bSJed Brown   next_name:
10929566063dSJacob Faibussowitsch     PetscCall(PetscFree(names[i]));
1093c4762a1bSJed Brown   }
1094c4762a1bSJed Brown 
1095c4762a1bSJed Brown   /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1096c4762a1bSJed Brown   mod->maxComputed = -1;
1097c4762a1bSJed Brown   for (link = mod->functionalRegistry; link; link = link->next) {
1098c4762a1bSJed Brown     for (i = 0; i < mod->numCall; i++) {
1099c4762a1bSJed Brown       FunctionalLink call = mod->functionalCall[i];
1100ad540459SPierre Jolivet       if (link->func == call->func && link->ctx == call->ctx) mod->maxComputed = PetscMax(mod->maxComputed, link->offset);
1101c4762a1bSJed Brown     }
1102c4762a1bSJed Brown   }
1103c4762a1bSJed Brown   PetscFunctionReturn(0);
1104c4762a1bSJed Brown }
1105c4762a1bSJed Brown 
1106d71ae5a4SJacob Faibussowitsch static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1107d71ae5a4SJacob Faibussowitsch {
1108c4762a1bSJed Brown   FunctionalLink l, next;
1109c4762a1bSJed Brown 
1110c4762a1bSJed Brown   PetscFunctionBeginUser;
1111c4762a1bSJed Brown   if (!link) PetscFunctionReturn(0);
1112c4762a1bSJed Brown   l     = *link;
1113c4762a1bSJed Brown   *link = NULL;
1114c4762a1bSJed Brown   for (; l; l = next) {
1115c4762a1bSJed Brown     next = l->next;
11169566063dSJacob Faibussowitsch     PetscCall(PetscFree(l->name));
11179566063dSJacob Faibussowitsch     PetscCall(PetscFree(l));
1118c4762a1bSJed Brown   }
1119c4762a1bSJed Brown   PetscFunctionReturn(0);
1120c4762a1bSJed Brown }
1121c4762a1bSJed Brown 
1122c4762a1bSJed Brown /* put the solution callback into a functional callback */
1123d71ae5a4SJacob Faibussowitsch static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx)
1124d71ae5a4SJacob Faibussowitsch {
1125c4762a1bSJed Brown   Model mod;
11267510d9b0SBarry Smith   PetscFunctionBeginUser;
1127c4762a1bSJed Brown   mod = (Model)modctx;
11289566063dSJacob Faibussowitsch   PetscCall((*mod->solution)(mod, time, x, u, mod->solutionctx));
1129c4762a1bSJed Brown   PetscFunctionReturn(0);
1130c4762a1bSJed Brown }
1131c4762a1bSJed Brown 
1132d71ae5a4SJacob Faibussowitsch PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
1133d71ae5a4SJacob Faibussowitsch {
1134c4762a1bSJed Brown   PetscErrorCode (*func[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1135c4762a1bSJed Brown   void *ctx[1];
1136c4762a1bSJed Brown   Model mod = user->model;
1137c4762a1bSJed Brown 
1138c4762a1bSJed Brown   PetscFunctionBeginUser;
1139c4762a1bSJed Brown   func[0] = SolutionFunctional;
1140c4762a1bSJed Brown   ctx[0]  = (void *)mod;
11419566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, func, ctx, INSERT_ALL_VALUES, X));
1142c4762a1bSJed Brown   PetscFunctionReturn(0);
1143c4762a1bSJed Brown }
1144c4762a1bSJed Brown 
1145d71ae5a4SJacob Faibussowitsch static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer)
1146d71ae5a4SJacob Faibussowitsch {
1147c4762a1bSJed Brown   PetscFunctionBeginUser;
11489566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer));
11499566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERVTK));
11509566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*viewer, filename));
1151c4762a1bSJed Brown   PetscFunctionReturn(0);
1152c4762a1bSJed Brown }
1153c4762a1bSJed Brown 
1154d71ae5a4SJacob Faibussowitsch static PetscErrorCode MonitorVTK(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx)
1155d71ae5a4SJacob Faibussowitsch {
1156c4762a1bSJed Brown   User        user = (User)ctx;
11573e9753d6SMatthew G. Knepley   DM          dm, plex;
1158c4762a1bSJed Brown   PetscViewer viewer;
1159c4762a1bSJed Brown   char        filename[PETSC_MAX_PATH_LEN], *ftable = NULL;
1160c4762a1bSJed Brown   PetscReal   xnorm;
1161c4762a1bSJed Brown 
1162c4762a1bSJed Brown   PetscFunctionBeginUser;
11639566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)X, "u"));
11649566063dSJacob Faibussowitsch   PetscCall(VecGetDM(X, &dm));
11659566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_INFINITY, &xnorm));
1166c4762a1bSJed Brown 
1167ad540459SPierre Jolivet   if (stepnum >= 0) stepnum += user->monitorStepOffset;
1168c4762a1bSJed Brown   if (stepnum >= 0) { /* No summary for final time */
1169c4762a1bSJed Brown     Model              mod = user->model;
11703e9753d6SMatthew G. Knepley     Vec                cellgeom;
1171c4762a1bSJed Brown     PetscInt           c, cStart, cEnd, fcount, i;
1172c4762a1bSJed Brown     size_t             ftableused, ftablealloc;
1173c4762a1bSJed Brown     const PetscScalar *cgeom, *x;
1174c4762a1bSJed Brown     DM                 dmCell;
1175c4762a1bSJed Brown     DMLabel            vtkLabel;
1176c4762a1bSJed Brown     PetscReal         *fmin, *fmax, *fintegral, *ftmp;
11773e9753d6SMatthew G. Knepley 
11789566063dSJacob Faibussowitsch     PetscCall(DMConvert(dm, DMPLEX, &plex));
11799566063dSJacob Faibussowitsch     PetscCall(DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL));
1180c4762a1bSJed Brown     fcount = mod->maxComputed + 1;
11819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(fcount, &fmin, fcount, &fmax, fcount, &fintegral, fcount, &ftmp));
1182c4762a1bSJed Brown     for (i = 0; i < fcount; i++) {
1183c4762a1bSJed Brown       fmin[i]      = PETSC_MAX_REAL;
1184c4762a1bSJed Brown       fmax[i]      = PETSC_MIN_REAL;
1185c4762a1bSJed Brown       fintegral[i] = 0;
1186c4762a1bSJed Brown     }
11879566063dSJacob Faibussowitsch     PetscCall(VecGetDM(cellgeom, &dmCell));
11889566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(dmCell, 0, &cStart, &cEnd));
11899566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(cellgeom, &cgeom));
11909566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(X, &x));
11919566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "vtk", &vtkLabel));
1192c4762a1bSJed Brown     for (c = cStart; c < cEnd; ++c) {
1193c4762a1bSJed Brown       PetscFVCellGeom   *cg;
1194c4762a1bSJed Brown       const PetscScalar *cx     = NULL;
1195c4762a1bSJed Brown       PetscInt           vtkVal = 0;
1196c4762a1bSJed Brown 
1197c4762a1bSJed Brown       /* not that these two routines as currently implemented work for any dm with a
1198c4762a1bSJed Brown        * localSection/globalSection */
11999566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg));
12009566063dSJacob Faibussowitsch       PetscCall(DMPlexPointGlobalRead(dm, c, x, &cx));
12019566063dSJacob Faibussowitsch       if (vtkLabel) PetscCall(DMLabelGetValue(vtkLabel, c, &vtkVal));
1202c4762a1bSJed Brown       if (!vtkVal || !cx) continue; /* ghost, or not a global cell */
1203c4762a1bSJed Brown       for (i = 0; i < mod->numCall; i++) {
1204c4762a1bSJed Brown         FunctionalLink flink = mod->functionalCall[i];
12059566063dSJacob Faibussowitsch         PetscCall((*flink->func)(mod, time, cg->centroid, cx, ftmp, flink->ctx));
1206c4762a1bSJed Brown       }
1207c4762a1bSJed Brown       for (i = 0; i < fcount; i++) {
1208c4762a1bSJed Brown         fmin[i] = PetscMin(fmin[i], ftmp[i]);
1209c4762a1bSJed Brown         fmax[i] = PetscMax(fmax[i], ftmp[i]);
1210c4762a1bSJed Brown         fintegral[i] += cg->volume * ftmp[i];
1211c4762a1bSJed Brown       }
1212c4762a1bSJed Brown     }
12139566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
12149566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X, &x));
12159566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&plex));
12169566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts)));
12179566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts)));
12189566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fintegral, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts)));
1219c4762a1bSJed Brown 
1220c4762a1bSJed Brown     ftablealloc = fcount * 100;
1221c4762a1bSJed Brown     ftableused  = 0;
12229566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ftablealloc, &ftable));
1223c4762a1bSJed Brown     for (i = 0; i < mod->numMonitored; i++) {
1224c4762a1bSJed Brown       size_t         countused;
1225c4762a1bSJed Brown       char           buffer[256], *p;
1226c4762a1bSJed Brown       FunctionalLink flink = mod->functionalMonitored[i];
1227c4762a1bSJed Brown       PetscInt       id    = flink->offset;
1228c4762a1bSJed Brown       if (i % 3) {
12299566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(buffer, "  ", 2));
1230c4762a1bSJed Brown         p = buffer + 2;
1231c4762a1bSJed Brown       } else if (i) {
1232c4762a1bSJed Brown         char newline[] = "\n";
12339566063dSJacob Faibussowitsch         PetscCall(PetscMemcpy(buffer, newline, sizeof(newline) - 1));
1234c4762a1bSJed Brown         p = buffer + sizeof(newline) - 1;
1235c4762a1bSJed Brown       } else {
1236c4762a1bSJed Brown         p = buffer;
1237c4762a1bSJed Brown       }
12389566063dSJacob Faibussowitsch       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]));
1239c4762a1bSJed Brown       countused--;
1240c4762a1bSJed Brown       countused += p - buffer;
1241c4762a1bSJed Brown       if (countused > ftablealloc - ftableused - 1) { /* reallocate */
1242c4762a1bSJed Brown         char *ftablenew;
1243c4762a1bSJed Brown         ftablealloc = 2 * ftablealloc + countused;
12449566063dSJacob Faibussowitsch         PetscCall(PetscMalloc(ftablealloc, &ftablenew));
12459566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(ftablenew, ftable, ftableused));
12469566063dSJacob Faibussowitsch         PetscCall(PetscFree(ftable));
1247c4762a1bSJed Brown         ftable = ftablenew;
1248c4762a1bSJed Brown       }
12499566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(ftable + ftableused, buffer, countused));
1250c4762a1bSJed Brown       ftableused += countused;
1251c4762a1bSJed Brown       ftable[ftableused] = 0;
1252c4762a1bSJed Brown     }
12539566063dSJacob Faibussowitsch     PetscCall(PetscFree4(fmin, fmax, fintegral, ftmp));
1254c4762a1bSJed Brown 
125563a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "% 3" PetscInt_FMT "  time %8.4g  |x| %8.4g  %s\n", stepnum, (double)time, (double)xnorm, ftable ? ftable : ""));
12569566063dSJacob Faibussowitsch     PetscCall(PetscFree(ftable));
1257c4762a1bSJed Brown   }
1258c4762a1bSJed Brown   if (user->vtkInterval < 1) PetscFunctionReturn(0);
1259c4762a1bSJed Brown   if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1260c4762a1bSJed Brown     if (stepnum == -1) { /* Final time is not multiple of normal time interval, write it anyway */
12619566063dSJacob Faibussowitsch       PetscCall(TSGetStepNumber(ts, &stepnum));
1262c4762a1bSJed Brown     }
126363a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(filename, sizeof filename, "%s-%03" PetscInt_FMT ".vtu", user->outputBasename, stepnum));
12649566063dSJacob Faibussowitsch     PetscCall(OutputVTK(dm, filename, &viewer));
12659566063dSJacob Faibussowitsch     PetscCall(VecView(X, viewer));
12669566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1267c4762a1bSJed Brown   }
1268c4762a1bSJed Brown   PetscFunctionReturn(0);
1269c4762a1bSJed Brown }
1270c4762a1bSJed Brown 
1271d71ae5a4SJacob Faibussowitsch static PetscErrorCode initializeTS(DM dm, User user, TS *ts)
1272d71ae5a4SJacob Faibussowitsch {
12737510d9b0SBarry Smith   PetscFunctionBeginUser;
12749566063dSJacob Faibussowitsch   PetscCall(TSCreate(PetscObjectComm((PetscObject)dm), ts));
12759566063dSJacob Faibussowitsch   PetscCall(TSSetType(*ts, TSSSP));
12769566063dSJacob Faibussowitsch   PetscCall(TSSetDM(*ts, dm));
12771baa6e33SBarry Smith   if (user->vtkmon) PetscCall(TSMonitorSet(*ts, MonitorVTK, user, NULL));
12789566063dSJacob Faibussowitsch   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user));
12799566063dSJacob Faibussowitsch   PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user));
12809566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(*ts, 2.0));
12819566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
1282c4762a1bSJed Brown   PetscFunctionReturn(0);
1283c4762a1bSJed Brown }
1284c4762a1bSJed Brown 
1285d71ae5a4SJacob Faibussowitsch static PetscErrorCode adaptToleranceFVM(PetscFV fvm, TS ts, Vec sol, VecTagger refineTag, VecTagger coarsenTag, User user, TS *tsNew, Vec *solNew)
1286d71ae5a4SJacob Faibussowitsch {
1287c4762a1bSJed Brown   DM                 dm, gradDM, plex, cellDM, adaptedDM = NULL;
1288c4762a1bSJed Brown   Vec                cellGeom, faceGeom;
1289c4762a1bSJed Brown   PetscBool          isForest, computeGradient;
1290c4762a1bSJed Brown   Vec                grad, locGrad, locX, errVec;
1291c4762a1bSJed Brown   PetscInt           cStart, cEnd, c, dim, nRefine, nCoarsen;
1292c4762a1bSJed Brown   PetscReal          minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2], minInd, maxInd, time;
1293c4762a1bSJed Brown   PetscScalar       *errArray;
1294c4762a1bSJed Brown   const PetscScalar *pointVals;
1295c4762a1bSJed Brown   const PetscScalar *pointGrads;
1296c4762a1bSJed Brown   const PetscScalar *pointGeom;
1297c4762a1bSJed Brown   DMLabel            adaptLabel = NULL;
1298c4762a1bSJed Brown   IS                 refineIS, coarsenIS;
1299c4762a1bSJed Brown 
13007510d9b0SBarry Smith   PetscFunctionBeginUser;
13019566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &time));
13029566063dSJacob Faibussowitsch   PetscCall(VecGetDM(sol, &dm));
13039566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
13049566063dSJacob Faibussowitsch   PetscCall(PetscFVGetComputeGradients(fvm, &computeGradient));
13059566063dSJacob Faibussowitsch   PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE));
13069566063dSJacob Faibussowitsch   PetscCall(DMIsForest(dm, &isForest));
13079566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
13089566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM));
13099566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(plex, &locX));
13109566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL));
13119566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX));
13129566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(plex, sol, INSERT_VALUES, locX));
13139566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(gradDM, &grad));
13149566063dSJacob Faibussowitsch   PetscCall(DMPlexReconstructGradientsFVM(plex, locX, grad));
13159566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(gradDM, &locGrad));
13169566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad));
13179566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad));
13189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&grad));
13199566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &cStart, &cEnd));
13209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(locGrad, &pointGrads));
13219566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(cellGeom, &pointGeom));
13229566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(locX, &pointVals));
13239566063dSJacob Faibussowitsch   PetscCall(VecGetDM(cellGeom, &cellDM));
13249566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
13259566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)plex), cEnd - cStart, PETSC_DETERMINE, &errVec));
13269566063dSJacob Faibussowitsch   PetscCall(VecSetUp(errVec));
13279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(errVec, &errArray));
1328c4762a1bSJed Brown   for (c = cStart; c < cEnd; c++) {
1329c4762a1bSJed Brown     PetscReal        errInd = 0.;
1330c4762a1bSJed Brown     PetscScalar     *pointGrad;
1331c4762a1bSJed Brown     PetscScalar     *pointVal;
1332c4762a1bSJed Brown     PetscFVCellGeom *cg;
1333c4762a1bSJed Brown 
13349566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRead(gradDM, c, pointGrads, &pointGrad));
13359566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRead(cellDM, c, pointGeom, &cg));
13369566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRead(plex, c, pointVals, &pointVal));
1337c4762a1bSJed Brown 
13389566063dSJacob Faibussowitsch     PetscCall((user->model->errorIndicator)(dim, cg->volume, user->model->physics->dof, pointVal, pointGrad, &errInd, user->model->errorCtx));
1339c4762a1bSJed Brown     errArray[c - cStart] = errInd;
1340c4762a1bSJed Brown     minMaxInd[0]         = PetscMin(minMaxInd[0], errInd);
1341c4762a1bSJed Brown     minMaxInd[1]         = PetscMax(minMaxInd[1], errInd);
1342c4762a1bSJed Brown   }
13439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(errVec, &errArray));
13449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(locX, &pointVals));
13459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(cellGeom, &pointGeom));
13469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(locGrad, &pointGrads));
13479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&locGrad));
13489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&locX));
13499566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1350c4762a1bSJed Brown 
13519566063dSJacob Faibussowitsch   PetscCall(VecTaggerComputeIS(refineTag, errVec, &refineIS, NULL));
13529566063dSJacob Faibussowitsch   PetscCall(VecTaggerComputeIS(coarsenTag, errVec, &coarsenIS, NULL));
13539566063dSJacob Faibussowitsch   PetscCall(ISGetSize(refineIS, &nRefine));
13549566063dSJacob Faibussowitsch   PetscCall(ISGetSize(coarsenIS, &nCoarsen));
13559566063dSJacob Faibussowitsch   if (nRefine) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS));
13569566063dSJacob Faibussowitsch   if (nCoarsen) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_COARSEN, coarsenIS));
13579566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&coarsenIS));
13589566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&refineIS));
13599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&errVec));
1360c4762a1bSJed Brown 
13619566063dSJacob Faibussowitsch   PetscCall(PetscFVSetComputeGradients(fvm, computeGradient));
1362c4762a1bSJed Brown   minMaxInd[1] = -minMaxInd[1];
13639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(minMaxInd, minMaxIndGlobal, 2, MPIU_REAL, MPI_MIN, PetscObjectComm((PetscObject)dm)));
1364c4762a1bSJed Brown   minInd = minMaxIndGlobal[0];
1365c4762a1bSJed Brown   maxInd = -minMaxIndGlobal[1];
136663a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(ts, "error indicator range (%E, %E)\n", (double)minInd, (double)maxInd));
1367c4762a1bSJed Brown   if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */
13689566063dSJacob Faibussowitsch     PetscCall(DMAdaptLabel(dm, adaptLabel, &adaptedDM));
1369c4762a1bSJed Brown   }
13709566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&adaptLabel));
1371c4762a1bSJed Brown   if (adaptedDM) {
137263a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(ts, "Adapted mesh, marking %" PetscInt_FMT " cells for refinement, and %" PetscInt_FMT " cells for coarsening\n", nRefine, nCoarsen));
13739566063dSJacob Faibussowitsch     if (tsNew) PetscCall(initializeTS(adaptedDM, user, tsNew));
1374c4762a1bSJed Brown     if (solNew) {
13759566063dSJacob Faibussowitsch       PetscCall(DMCreateGlobalVector(adaptedDM, solNew));
13769566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)*solNew, "solution"));
13779566063dSJacob Faibussowitsch       PetscCall(DMForestTransferVec(dm, sol, adaptedDM, *solNew, PETSC_TRUE, time));
1378c4762a1bSJed Brown     }
13799566063dSJacob Faibussowitsch     if (isForest) PetscCall(DMForestSetAdaptivityForest(adaptedDM, NULL)); /* clear internal references to the previous dm */
13809566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&adaptedDM));
1381c4762a1bSJed Brown   } else {
1382c4762a1bSJed Brown     if (tsNew) *tsNew = NULL;
1383c4762a1bSJed Brown     if (solNew) *solNew = NULL;
1384c4762a1bSJed Brown   }
1385c4762a1bSJed Brown   PetscFunctionReturn(0);
1386c4762a1bSJed Brown }
1387c4762a1bSJed Brown 
1388d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
1389d71ae5a4SJacob Faibussowitsch {
1390c4762a1bSJed Brown   MPI_Comm          comm;
1391c4762a1bSJed Brown   PetscDS           prob;
1392c4762a1bSJed Brown   PetscFV           fvm;
1393c4762a1bSJed Brown   PetscLimiter      limiter = NULL, noneLimiter = NULL;
1394c4762a1bSJed Brown   User              user;
1395c4762a1bSJed Brown   Model             mod;
1396c4762a1bSJed Brown   Physics           phys;
13973e9753d6SMatthew G. Knepley   DM                dm, plex;
1398c4762a1bSJed Brown   PetscReal         ftime, cfl, dt, minRadius;
1399c4762a1bSJed Brown   PetscInt          dim, nsteps;
1400c4762a1bSJed Brown   TS                ts;
1401c4762a1bSJed Brown   TSConvergedReason reason;
1402c4762a1bSJed Brown   Vec               X;
1403c4762a1bSJed Brown   PetscViewer       viewer;
1404b5a892a1SMatthew G. Knepley   PetscBool         vtkCellGeom, useAMR;
140530602db0SMatthew G. Knepley   PetscInt          adaptInterval;
1406c4762a1bSJed Brown   char              physname[256] = "advect";
1407c4762a1bSJed Brown   VecTagger         refineTag = NULL, coarsenTag = NULL;
1408c4762a1bSJed Brown 
1409327415f7SBarry Smith   PetscFunctionBeginUser;
14109566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1411c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
1412c4762a1bSJed Brown 
14139566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user));
14149566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user->model));
14159566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user->model->physics));
1416c4762a1bSJed Brown   mod           = user->model;
1417c4762a1bSJed Brown   phys          = mod->physics;
1418c4762a1bSJed Brown   mod->comm     = comm;
1419c4762a1bSJed Brown   useAMR        = PETSC_FALSE;
1420c4762a1bSJed Brown   adaptInterval = 1;
1421c4762a1bSJed Brown 
1422c4762a1bSJed Brown   /* Register physical models to be available on the command line */
14239566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PhysicsList, "advect", PhysicsCreate_Advect));
14249566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PhysicsList, "sw", PhysicsCreate_SW));
14259566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PhysicsList, "euler", PhysicsCreate_Euler));
1426c4762a1bSJed Brown 
1427d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "Unstructured Finite Volume Mesh Options", "");
1428c4762a1bSJed Brown   {
1429c4762a1bSJed Brown     cfl = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
14309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ufv_cfl", "CFL number per step", "", cfl, &cfl, NULL));
1431c4762a1bSJed Brown     user->vtkInterval = 1;
14329566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ufv_vtk_interval", "VTK output interval (0 to disable)", "", user->vtkInterval, &user->vtkInterval, NULL));
1433c4762a1bSJed Brown     user->vtkmon = PETSC_TRUE;
14349566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ufv_vtk_monitor", "Use VTKMonitor routine", "", user->vtkmon, &user->vtkmon, NULL));
1435c4762a1bSJed Brown     vtkCellGeom = PETSC_FALSE;
14369566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(user->outputBasename, "ex11"));
14379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-ufv_vtk_basename", "VTK output basename", "", user->outputBasename, user->outputBasename, sizeof(user->outputBasename), NULL));
14389566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ufv_vtk_cellgeom", "Write cell geometry (for debugging)", "", vtkCellGeom, &vtkCellGeom, NULL));
14399566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ufv_use_amr", "use local adaptive mesh refinement", "", useAMR, &useAMR, NULL));
14409566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ufv_adapt_interval", "time steps between AMR", "", adaptInterval, &adaptInterval, NULL));
1441c4762a1bSJed Brown   }
1442d0609cedSBarry Smith   PetscOptionsEnd();
1443c4762a1bSJed Brown 
1444c4762a1bSJed Brown   if (useAMR) {
1445c4762a1bSJed Brown     VecTaggerBox refineBox, coarsenBox;
1446c4762a1bSJed Brown 
1447c4762a1bSJed Brown     refineBox.min = refineBox.max = PETSC_MAX_REAL;
1448c4762a1bSJed Brown     coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL;
1449c4762a1bSJed Brown 
14509566063dSJacob Faibussowitsch     PetscCall(VecTaggerCreate(comm, &refineTag));
14519566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)refineTag, "refine_"));
14529566063dSJacob Faibussowitsch     PetscCall(VecTaggerSetType(refineTag, VECTAGGERABSOLUTE));
14539566063dSJacob Faibussowitsch     PetscCall(VecTaggerAbsoluteSetBox(refineTag, &refineBox));
14549566063dSJacob Faibussowitsch     PetscCall(VecTaggerSetFromOptions(refineTag));
14559566063dSJacob Faibussowitsch     PetscCall(VecTaggerSetUp(refineTag));
14569566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)refineTag, NULL, "-tag_view"));
1457c4762a1bSJed Brown 
14589566063dSJacob Faibussowitsch     PetscCall(VecTaggerCreate(comm, &coarsenTag));
14599566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)coarsenTag, "coarsen_"));
14609566063dSJacob Faibussowitsch     PetscCall(VecTaggerSetType(coarsenTag, VECTAGGERABSOLUTE));
14619566063dSJacob Faibussowitsch     PetscCall(VecTaggerAbsoluteSetBox(coarsenTag, &coarsenBox));
14629566063dSJacob Faibussowitsch     PetscCall(VecTaggerSetFromOptions(coarsenTag));
14639566063dSJacob Faibussowitsch     PetscCall(VecTaggerSetUp(coarsenTag));
14649566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)coarsenTag, NULL, "-tag_view"));
1465c4762a1bSJed Brown   }
1466c4762a1bSJed Brown 
1467d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "Unstructured Finite Volume Physics Options", "");
1468c4762a1bSJed Brown   {
1469c4762a1bSJed Brown     PetscErrorCode (*physcreate)(Model, Physics, PetscOptionItems *);
14709566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-physics", "Physics module to solve", "", PhysicsList, physname, physname, sizeof physname, NULL));
14719566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(PhysicsList, physname, &physcreate));
14729566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(phys, sizeof(struct _n_Physics)));
14739566063dSJacob Faibussowitsch     PetscCall((*physcreate)(mod, phys, PetscOptionsObject));
1474c4762a1bSJed Brown     /* Count number of fields and dofs */
1475c4762a1bSJed Brown     for (phys->nfields = 0, phys->dof = 0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
147654c59aa7SJacob Faibussowitsch     PetscCheck(phys->dof > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set dof", physname);
14779566063dSJacob Faibussowitsch     PetscCall(ModelFunctionalSetFromOptions(mod, PetscOptionsObject));
1478c4762a1bSJed Brown   }
1479d0609cedSBarry Smith   PetscOptionsEnd();
1480c4762a1bSJed Brown 
1481c4762a1bSJed Brown   /* Create mesh */
1482c4762a1bSJed Brown   {
148330602db0SMatthew G. Knepley     PetscInt i;
148430602db0SMatthew G. Knepley 
14859566063dSJacob Faibussowitsch     PetscCall(DMCreate(comm, &dm));
14869566063dSJacob Faibussowitsch     PetscCall(DMSetType(dm, DMPLEX));
14879566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(dm));
14889371c9d4SSatish Balay     for (i = 0; i < DIM; i++) {
14899371c9d4SSatish Balay       mod->bounds[2 * i]     = 0.;
14909371c9d4SSatish Balay       mod->bounds[2 * i + 1] = 1.;
14919371c9d4SSatish Balay     };
1492c4762a1bSJed Brown     dim = DIM;
149330602db0SMatthew G. Knepley     { /* a null name means just do a hex box */
149430602db0SMatthew G. Knepley       PetscInt  cells[3] = {1, 1, 1}, n = 3;
149530602db0SMatthew G. Knepley       PetscBool flg2, skew              = PETSC_FALSE;
1496c4762a1bSJed Brown       PetscInt  nret2 = 2 * DIM;
1497d0609cedSBarry Smith       PetscOptionsBegin(comm, NULL, "Rectangular mesh options", "");
14989566063dSJacob Faibussowitsch       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));
14999566063dSJacob Faibussowitsch       PetscCall(PetscOptionsBool("-grid_skew_60", "Skew grid for 60 degree shock mesh", "", skew, &skew, NULL));
15009566063dSJacob Faibussowitsch       PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", cells, &n, NULL));
1501d0609cedSBarry Smith       PetscOptionsEnd();
150230602db0SMatthew G. Knepley       /* TODO Rewrite this with Mark, and remove grid_bounds at that time */
1503c4762a1bSJed Brown       if (flg2) {
1504c4762a1bSJed Brown         PetscInt     dimEmbed, i;
1505c4762a1bSJed Brown         PetscInt     nCoords;
1506c4762a1bSJed Brown         PetscScalar *coords;
1507c4762a1bSJed Brown         Vec          coordinates;
1508c4762a1bSJed Brown 
15099566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
15109566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
15119566063dSJacob Faibussowitsch         PetscCall(VecGetLocalSize(coordinates, &nCoords));
151254c59aa7SJacob Faibussowitsch         PetscCheck(!(nCoords % dimEmbed), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Coordinate vector the wrong size");
15139566063dSJacob Faibussowitsch         PetscCall(VecGetArray(coordinates, &coords));
1514c4762a1bSJed Brown         for (i = 0; i < nCoords; i += dimEmbed) {
1515c4762a1bSJed Brown           PetscInt j;
1516c4762a1bSJed Brown 
1517c4762a1bSJed Brown           PetscScalar *coord = &coords[i];
1518c4762a1bSJed Brown           for (j = 0; j < dimEmbed; j++) {
1519c4762a1bSJed Brown             coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]);
1520c4762a1bSJed Brown             if (dim == 2 && cells[1] == 1 && j == 0 && skew) {
1521c4762a1bSJed Brown               if (cells[0] == 2 && i == 8) {
1522c4762a1bSJed Brown                 coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */
152330602db0SMatthew G. Knepley               } else if (cells[0] == 3) {
1524c4762a1bSJed Brown                 if (i == 2 || i == 10) coord[j] = mod->bounds[1] / 4.;
1525c4762a1bSJed Brown                 else if (i == 4) coord[j] = mod->bounds[1] / 2.;
1526c4762a1bSJed Brown                 else if (i == 12) coord[j] = 1.57735026918963 * mod->bounds[1] / 2.;
1527c4762a1bSJed Brown               }
1528c4762a1bSJed Brown             }
1529c4762a1bSJed Brown           }
1530c4762a1bSJed Brown         }
15319566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(coordinates, &coords));
15329566063dSJacob Faibussowitsch         PetscCall(DMSetCoordinatesLocal(dm, coordinates));
1533c4762a1bSJed Brown       }
1534c4762a1bSJed Brown     }
1535c4762a1bSJed Brown   }
15369566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-orig_dm_view"));
15379566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1538c4762a1bSJed Brown 
1539c4762a1bSJed Brown   /* set up BCs, functions, tags */
15409566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "Face Sets"));
1541c4762a1bSJed Brown   mod->errorIndicator = ErrorIndicator_Simple;
1542c4762a1bSJed Brown 
1543c4762a1bSJed Brown   {
1544c4762a1bSJed Brown     DM gdm;
1545c4762a1bSJed Brown 
15469566063dSJacob Faibussowitsch     PetscCall(DMPlexConstructGhostCells(dm, NULL, NULL, &gdm));
15479566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm));
1548c4762a1bSJed Brown     dm = gdm;
15499566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1550c4762a1bSJed Brown   }
1551c4762a1bSJed Brown 
15529566063dSJacob Faibussowitsch   PetscCall(PetscFVCreate(comm, &fvm));
15539566063dSJacob Faibussowitsch   PetscCall(PetscFVSetFromOptions(fvm));
15549566063dSJacob Faibussowitsch   PetscCall(PetscFVSetNumComponents(fvm, phys->dof));
15559566063dSJacob Faibussowitsch   PetscCall(PetscFVSetSpatialDimension(fvm, dim));
15569566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fvm, ""));
1557c4762a1bSJed Brown   {
1558c4762a1bSJed Brown     PetscInt f, dof;
1559c4762a1bSJed Brown     for (f = 0, dof = 0; f < phys->nfields; f++) {
1560c4762a1bSJed Brown       PetscInt newDof = phys->field_desc[f].dof;
1561c4762a1bSJed Brown 
1562c4762a1bSJed Brown       if (newDof == 1) {
15639566063dSJacob Faibussowitsch         PetscCall(PetscFVSetComponentName(fvm, dof, phys->field_desc[f].name));
15649371c9d4SSatish Balay       } else {
1565c4762a1bSJed Brown         PetscInt j;
1566c4762a1bSJed Brown 
1567c4762a1bSJed Brown         for (j = 0; j < newDof; j++) {
1568c4762a1bSJed Brown           char compName[256] = "Unknown";
1569c4762a1bSJed Brown 
157063a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(compName, sizeof(compName), "%s_%" PetscInt_FMT, phys->field_desc[f].name, j));
15719566063dSJacob Faibussowitsch           PetscCall(PetscFVSetComponentName(fvm, dof + j, compName));
1572c4762a1bSJed Brown         }
1573c4762a1bSJed Brown       }
1574c4762a1bSJed Brown       dof += newDof;
1575c4762a1bSJed Brown     }
1576c4762a1bSJed Brown   }
1577c4762a1bSJed Brown   /* FV is now structured with one field having all physics as components */
15789566063dSJacob Faibussowitsch   PetscCall(DMAddField(dm, NULL, (PetscObject)fvm));
15799566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
15809566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
15819566063dSJacob Faibussowitsch   PetscCall(PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann));
15829566063dSJacob Faibussowitsch   PetscCall(PetscDSSetContext(prob, 0, user->model->physics));
15839566063dSJacob Faibussowitsch   PetscCall((*mod->setupbc)(dm, prob, phys));
15849566063dSJacob Faibussowitsch   PetscCall(PetscDSSetFromOptions(prob));
1585c4762a1bSJed Brown   {
1586c4762a1bSJed Brown     char      convType[256];
1587c4762a1bSJed Brown     PetscBool flg;
1588c4762a1bSJed Brown 
1589d0609cedSBarry Smith     PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
15909566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex12", DMList, DMPLEX, convType, 256, &flg));
1591d0609cedSBarry Smith     PetscOptionsEnd();
1592c4762a1bSJed Brown     if (flg) {
1593c4762a1bSJed Brown       DM dmConv;
1594c4762a1bSJed Brown 
15959566063dSJacob Faibussowitsch       PetscCall(DMConvert(dm, convType, &dmConv));
1596c4762a1bSJed Brown       if (dmConv) {
15979566063dSJacob Faibussowitsch         PetscCall(DMViewFromOptions(dmConv, NULL, "-dm_conv_view"));
15989566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&dm));
1599c4762a1bSJed Brown         dm = dmConv;
16009566063dSJacob Faibussowitsch         PetscCall(DMSetFromOptions(dm));
1601c4762a1bSJed Brown       }
1602c4762a1bSJed Brown     }
1603c4762a1bSJed Brown   }
1604c4762a1bSJed Brown 
16059566063dSJacob Faibussowitsch   PetscCall(initializeTS(dm, user, &ts));
1606c4762a1bSJed Brown 
16079566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &X));
16089566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)X, "solution"));
16099566063dSJacob Faibussowitsch   PetscCall(SetInitialCondition(dm, X, user));
1610c4762a1bSJed Brown   if (useAMR) {
1611c4762a1bSJed Brown     PetscInt adaptIter;
1612c4762a1bSJed Brown 
1613c4762a1bSJed Brown     /* use no limiting when reconstructing gradients for adaptivity */
16149566063dSJacob Faibussowitsch     PetscCall(PetscFVGetLimiter(fvm, &limiter));
16159566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)limiter));
16169566063dSJacob Faibussowitsch     PetscCall(PetscLimiterCreate(PetscObjectComm((PetscObject)fvm), &noneLimiter));
16179566063dSJacob Faibussowitsch     PetscCall(PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE));
1618c4762a1bSJed Brown 
16199566063dSJacob Faibussowitsch     PetscCall(PetscFVSetLimiter(fvm, noneLimiter));
1620c4762a1bSJed Brown     for (adaptIter = 0;; ++adaptIter) {
1621c4762a1bSJed Brown       PetscLogDouble bytes;
1622c4762a1bSJed Brown       TS             tsNew = NULL;
1623c4762a1bSJed Brown 
16249566063dSJacob Faibussowitsch       PetscCall(PetscMemoryGetCurrentUsage(&bytes));
162563a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "refinement loop %" PetscInt_FMT ": memory used %g\n", adaptIter, (double)bytes));
16269566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(dm, NULL, "-initial_dm_view"));
16279566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(X, NULL, "-initial_vec_view"));
1628c4762a1bSJed Brown #if 0
1629c4762a1bSJed Brown       if (viewInitial) {
1630c4762a1bSJed Brown         PetscViewer viewer;
1631c4762a1bSJed Brown         char        buf[256];
1632c4762a1bSJed Brown         PetscBool   isHDF5, isVTK;
1633c4762a1bSJed Brown 
16349566063dSJacob Faibussowitsch         PetscCall(PetscViewerCreate(comm,&viewer));
16359566063dSJacob Faibussowitsch         PetscCall(PetscViewerSetType(viewer,PETSCVIEWERVTK));
16369566063dSJacob Faibussowitsch         PetscCall(PetscViewerSetOptionsPrefix(viewer,"initial_"));
16379566063dSJacob Faibussowitsch         PetscCall(PetscViewerSetFromOptions(viewer));
16389566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5));
16399566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK));
1640c4762a1bSJed Brown         if (isHDF5) {
164163a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(buf, 256, "ex11-initial-%" PetscInt_FMT ".h5", adaptIter));
1642c4762a1bSJed Brown         } else if (isVTK) {
164363a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(buf, 256, "ex11-initial-%" PetscInt_FMT ".vtu", adaptIter));
16449566063dSJacob Faibussowitsch           PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU));
1645c4762a1bSJed Brown         }
16469566063dSJacob Faibussowitsch         PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE));
16479566063dSJacob Faibussowitsch         PetscCall(PetscViewerFileSetName(viewer,buf));
1648c4762a1bSJed Brown         if (isHDF5) {
16499566063dSJacob Faibussowitsch           PetscCall(DMView(dm,viewer));
16509566063dSJacob Faibussowitsch           PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_UPDATE));
1651c4762a1bSJed Brown         }
16529566063dSJacob Faibussowitsch         PetscCall(VecView(X,viewer));
16539566063dSJacob Faibussowitsch         PetscCall(PetscViewerDestroy(&viewer));
1654c4762a1bSJed Brown       }
1655c4762a1bSJed Brown #endif
1656c4762a1bSJed Brown 
16579566063dSJacob Faibussowitsch       PetscCall(adaptToleranceFVM(fvm, ts, X, refineTag, coarsenTag, user, &tsNew, NULL));
1658c4762a1bSJed Brown       if (!tsNew) {
1659c4762a1bSJed Brown         break;
1660c4762a1bSJed Brown       } else {
16619566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&dm));
16629566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&X));
16639566063dSJacob Faibussowitsch         PetscCall(TSDestroy(&ts));
1664c4762a1bSJed Brown         ts = tsNew;
16659566063dSJacob Faibussowitsch         PetscCall(TSGetDM(ts, &dm));
16669566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)dm));
16679566063dSJacob Faibussowitsch         PetscCall(DMCreateGlobalVector(dm, &X));
16689566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetName((PetscObject)X, "solution"));
16699566063dSJacob Faibussowitsch         PetscCall(SetInitialCondition(dm, X, user));
1670c4762a1bSJed Brown       }
1671c4762a1bSJed Brown     }
1672c4762a1bSJed Brown     /* restore original limiter */
16739566063dSJacob Faibussowitsch     PetscCall(PetscFVSetLimiter(fvm, limiter));
1674c4762a1bSJed Brown   }
1675c4762a1bSJed Brown 
16769566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
1677c4762a1bSJed Brown   if (vtkCellGeom) {
1678c4762a1bSJed Brown     DM  dmCell;
1679c4762a1bSJed Brown     Vec cellgeom, partition;
1680c4762a1bSJed Brown 
16819566063dSJacob Faibussowitsch     PetscCall(DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL));
16829566063dSJacob Faibussowitsch     PetscCall(OutputVTK(dm, "ex11-cellgeom.vtk", &viewer));
16839566063dSJacob Faibussowitsch     PetscCall(VecView(cellgeom, viewer));
16849566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
16859566063dSJacob Faibussowitsch     PetscCall(CreatePartitionVec(dm, &dmCell, &partition));
16869566063dSJacob Faibussowitsch     PetscCall(OutputVTK(dmCell, "ex11-partition.vtk", &viewer));
16879566063dSJacob Faibussowitsch     PetscCall(VecView(partition, viewer));
16889566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
16899566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&partition));
16909566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmCell));
1691c4762a1bSJed Brown   }
1692c4762a1bSJed Brown   /* collect max maxspeed from all processes -- todo */
16939566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGeometryFVM(plex, NULL, NULL, &minRadius));
16949566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
16959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&phys->maxspeed, &mod->maxspeed, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts)));
169654c59aa7SJacob Faibussowitsch   PetscCheck(mod->maxspeed > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set maxspeed", physname);
1697c4762a1bSJed Brown   dt = cfl * minRadius / mod->maxspeed;
16989566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
16999566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1700c4762a1bSJed Brown   if (!useAMR) {
17019566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, X));
17029566063dSJacob Faibussowitsch     PetscCall(TSGetSolveTime(ts, &ftime));
17039566063dSJacob Faibussowitsch     PetscCall(TSGetStepNumber(ts, &nsteps));
1704c4762a1bSJed Brown   } else {
1705c4762a1bSJed Brown     PetscReal finalTime;
1706c4762a1bSJed Brown     PetscInt  adaptIter;
1707c4762a1bSJed Brown     TS        tsNew  = NULL;
1708c4762a1bSJed Brown     Vec       solNew = NULL;
1709c4762a1bSJed Brown 
17109566063dSJacob Faibussowitsch     PetscCall(TSGetMaxTime(ts, &finalTime));
17119566063dSJacob Faibussowitsch     PetscCall(TSSetMaxSteps(ts, adaptInterval));
17129566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, X));
17139566063dSJacob Faibussowitsch     PetscCall(TSGetSolveTime(ts, &ftime));
17149566063dSJacob Faibussowitsch     PetscCall(TSGetStepNumber(ts, &nsteps));
1715c4762a1bSJed Brown     for (adaptIter = 0; ftime < finalTime; adaptIter++) {
1716c4762a1bSJed Brown       PetscLogDouble bytes;
1717c4762a1bSJed Brown 
17189566063dSJacob Faibussowitsch       PetscCall(PetscMemoryGetCurrentUsage(&bytes));
171963a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "AMR time step loop %" PetscInt_FMT ": memory used %g\n", adaptIter, bytes));
17209566063dSJacob Faibussowitsch       PetscCall(PetscFVSetLimiter(fvm, noneLimiter));
17219566063dSJacob Faibussowitsch       PetscCall(adaptToleranceFVM(fvm, ts, X, refineTag, coarsenTag, user, &tsNew, &solNew));
17229566063dSJacob Faibussowitsch       PetscCall(PetscFVSetLimiter(fvm, limiter));
1723c4762a1bSJed Brown       if (tsNew) {
17249566063dSJacob Faibussowitsch         PetscCall(PetscInfo(ts, "AMR used\n"));
17259566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&dm));
17269566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&X));
17279566063dSJacob Faibussowitsch         PetscCall(TSDestroy(&ts));
1728c4762a1bSJed Brown         ts = tsNew;
1729c4762a1bSJed Brown         X  = solNew;
17309566063dSJacob Faibussowitsch         PetscCall(TSSetFromOptions(ts));
17319566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dm));
17329566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)dm));
17339566063dSJacob Faibussowitsch         PetscCall(DMConvert(dm, DMPLEX, &plex));
17349566063dSJacob Faibussowitsch         PetscCall(DMPlexGetGeometryFVM(dm, NULL, NULL, &minRadius));
17359566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&plex));
17369566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Allreduce(&phys->maxspeed, &mod->maxspeed, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts)));
173754c59aa7SJacob Faibussowitsch         PetscCheck(mod->maxspeed > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set maxspeed", physname);
1738c4762a1bSJed Brown         dt = cfl * minRadius / mod->maxspeed;
17399566063dSJacob Faibussowitsch         PetscCall(TSSetStepNumber(ts, nsteps));
17409566063dSJacob Faibussowitsch         PetscCall(TSSetTime(ts, ftime));
17419566063dSJacob Faibussowitsch         PetscCall(TSSetTimeStep(ts, dt));
1742c4762a1bSJed Brown       } else {
17439566063dSJacob Faibussowitsch         PetscCall(PetscInfo(ts, "AMR not used\n"));
1744c4762a1bSJed Brown       }
1745c4762a1bSJed Brown       user->monitorStepOffset = nsteps;
17469566063dSJacob Faibussowitsch       PetscCall(TSSetMaxSteps(ts, nsteps + adaptInterval));
17479566063dSJacob Faibussowitsch       PetscCall(TSSolve(ts, X));
17489566063dSJacob Faibussowitsch       PetscCall(TSGetSolveTime(ts, &ftime));
17499566063dSJacob Faibussowitsch       PetscCall(TSGetStepNumber(ts, &nsteps));
1750c4762a1bSJed Brown     }
1751c4762a1bSJed Brown   }
17529566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts, &reason));
175363a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, nsteps));
17549566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1755c4762a1bSJed Brown 
17569566063dSJacob Faibussowitsch   PetscCall(VecTaggerDestroy(&refineTag));
17579566063dSJacob Faibussowitsch   PetscCall(VecTaggerDestroy(&coarsenTag));
17589566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&PhysicsList));
17599566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&PhysicsRiemannList_SW));
17609566063dSJacob Faibussowitsch   PetscCall(FunctionalLinkDestroy(&user->model->functionalRegistry));
17619566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->model->functionalMonitored));
17629566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->model->functionalCall));
17639566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->model->physics->data));
17649566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->model->physics));
17659566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->model));
17669566063dSJacob Faibussowitsch   PetscCall(PetscFree(user));
17679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
17689566063dSJacob Faibussowitsch   PetscCall(PetscLimiterDestroy(&limiter));
17699566063dSJacob Faibussowitsch   PetscCall(PetscLimiterDestroy(&noneLimiter));
17709566063dSJacob Faibussowitsch   PetscCall(PetscFVDestroy(&fvm));
17719566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
17729566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1773b122ec5aSJacob Faibussowitsch   return 0;
1774c4762a1bSJed Brown }
1775c4762a1bSJed Brown 
1776c4762a1bSJed Brown /* Godunov fluxs */
1777d71ae5a4SJacob Faibussowitsch PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
1778d71ae5a4SJacob Faibussowitsch {
1779c4762a1bSJed Brown   /* System generated locals */
1780c4762a1bSJed Brown   PetscScalar ret_val;
1781c4762a1bSJed Brown 
1782ad540459SPierre Jolivet   if (PetscRealPart(*test) > 0.) goto L10;
1783c4762a1bSJed Brown   ret_val = *b;
1784c4762a1bSJed Brown   return ret_val;
1785c4762a1bSJed Brown L10:
1786c4762a1bSJed Brown   ret_val = *a;
1787c4762a1bSJed Brown   return ret_val;
1788c4762a1bSJed Brown } /* cvmgp_ */
1789c4762a1bSJed Brown 
1790d71ae5a4SJacob Faibussowitsch PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
1791d71ae5a4SJacob Faibussowitsch {
1792c4762a1bSJed Brown   /* System generated locals */
1793c4762a1bSJed Brown   PetscScalar ret_val;
1794c4762a1bSJed Brown 
1795ad540459SPierre Jolivet   if (PetscRealPart(*test) < 0.) goto L10;
1796c4762a1bSJed Brown   ret_val = *b;
1797c4762a1bSJed Brown   return ret_val;
1798c4762a1bSJed Brown L10:
1799c4762a1bSJed Brown   ret_val = *a;
1800c4762a1bSJed Brown   return ret_val;
1801c4762a1bSJed Brown } /* cvmgm_ */
1802c4762a1bSJed Brown 
1803d71ae5a4SJacob Faibussowitsch 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)
1804d71ae5a4SJacob Faibussowitsch {
1805c4762a1bSJed Brown   /* Initialized data */
1806c4762a1bSJed Brown 
1807c4762a1bSJed Brown   static PetscScalar smallp = 1e-8;
1808c4762a1bSJed Brown 
1809c4762a1bSJed Brown   /* System generated locals */
1810c4762a1bSJed Brown   int         i__1;
1811c4762a1bSJed Brown   PetscScalar d__1, d__2;
1812c4762a1bSJed Brown 
1813c4762a1bSJed Brown   /* Local variables */
1814c4762a1bSJed Brown   static int         i0;
1815c4762a1bSJed Brown   static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2;
1816c4762a1bSJed Brown   static int         iwave;
1817c4762a1bSJed Brown   static PetscScalar gascl4, gascr4, cstarl, dpstar, cstarr;
1818c4762a1bSJed Brown   /* static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascr1, gascr2, gascr3; */
1819c4762a1bSJed Brown   static int         iterno;
1820c4762a1bSJed Brown   static PetscScalar ustarl, ustarr, rarepr1, rarepr2;
1821c4762a1bSJed Brown 
1822c4762a1bSJed Brown   /* gascl1 = *gaml - 1.; */
1823c4762a1bSJed Brown   /* gascl2 = (*gaml + 1.) * .5; */
1824c4762a1bSJed Brown   /* gascl3 = gascl2 / *gaml; */
1825c4762a1bSJed Brown   gascl4 = 1. / (*gaml - 1.);
1826c4762a1bSJed Brown 
1827c4762a1bSJed Brown   /* gascr1 = *gamr - 1.; */
1828c4762a1bSJed Brown   /* gascr2 = (*gamr + 1.) * .5; */
1829c4762a1bSJed Brown   /* gascr3 = gascr2 / *gamr; */
1830c4762a1bSJed Brown   gascr4 = 1. / (*gamr - 1.);
1831c4762a1bSJed Brown   iterno = 10;
1832c4762a1bSJed Brown   /*        find pstar: */
1833c4762a1bSJed Brown   cl = PetscSqrtScalar(*gaml * *pl / *rl);
1834c4762a1bSJed Brown   cr = PetscSqrtScalar(*gamr * *pr / *rr);
1835c4762a1bSJed Brown   wl = *rl * cl;
1836c4762a1bSJed Brown   wr = *rr * cr;
1837c4762a1bSJed Brown   /* csqrl = wl * wl; */
1838c4762a1bSJed Brown   /* csqrr = wr * wr; */
1839c4762a1bSJed Brown   *pstar  = (wl * *pr + wr * *pl) / (wl + wr);
1840c4762a1bSJed Brown   *pstar  = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
1841c4762a1bSJed Brown   pst     = *pl / *pr;
1842c4762a1bSJed Brown   skpr1   = cr * (pst - 1.) * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
1843c4762a1bSJed Brown   d__1    = (*gamr - 1.) / (*gamr * 2.);
1844c4762a1bSJed Brown   rarepr2 = gascr4 * 2. * cr * (1. - PetscPowScalar(pst, d__1));
1845c4762a1bSJed Brown   pst     = *pr / *pl;
1846c4762a1bSJed Brown   skpr2   = cl * (pst - 1.) * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
1847c4762a1bSJed Brown   d__1    = (*gaml - 1.) / (*gaml * 2.);
1848c4762a1bSJed Brown   rarepr1 = gascl4 * 2. * cl * (1. - PetscPowScalar(pst, d__1));
1849c4762a1bSJed Brown   durl    = *uxr - *uxl;
1850c4762a1bSJed Brown   if (PetscRealPart(*pr) < PetscRealPart(*pl)) {
1851c4762a1bSJed Brown     if (PetscRealPart(durl) >= PetscRealPart(rarepr1)) {
1852c4762a1bSJed Brown       iwave = 100;
1853c4762a1bSJed Brown     } else if (PetscRealPart(durl) <= PetscRealPart(-skpr1)) {
1854c4762a1bSJed Brown       iwave = 300;
1855c4762a1bSJed Brown     } else {
1856c4762a1bSJed Brown       iwave = 400;
1857c4762a1bSJed Brown     }
1858c4762a1bSJed Brown   } else {
1859c4762a1bSJed Brown     if (PetscRealPart(durl) >= PetscRealPart(rarepr2)) {
1860c4762a1bSJed Brown       iwave = 100;
1861c4762a1bSJed Brown     } else if (PetscRealPart(durl) <= PetscRealPart(-skpr2)) {
1862c4762a1bSJed Brown       iwave = 300;
1863c4762a1bSJed Brown     } else {
1864c4762a1bSJed Brown       iwave = 200;
1865c4762a1bSJed Brown     }
1866c4762a1bSJed Brown   }
1867c4762a1bSJed Brown   if (iwave == 100) {
1868c4762a1bSJed Brown     /*     1-wave: rarefaction wave, 3-wave: rarefaction wave */
1869c4762a1bSJed Brown     /*     case (100) */
1870c4762a1bSJed Brown     i__1 = iterno;
1871c4762a1bSJed Brown     for (i0 = 1; i0 <= i__1; ++i0) {
1872c4762a1bSJed Brown       d__1    = *pstar / *pl;
1873c4762a1bSJed Brown       d__2    = 1. / *gaml;
1874c4762a1bSJed Brown       *rstarl = *rl * PetscPowScalar(d__1, d__2);
1875c4762a1bSJed Brown       cstarl  = PetscSqrtScalar(*gaml * *pstar / *rstarl);
1876c4762a1bSJed Brown       ustarl  = *uxl - gascl4 * 2. * (cstarl - cl);
1877c4762a1bSJed Brown       zl      = *rstarl * cstarl;
1878c4762a1bSJed Brown       d__1    = *pstar / *pr;
1879c4762a1bSJed Brown       d__2    = 1. / *gamr;
1880c4762a1bSJed Brown       *rstarr = *rr * PetscPowScalar(d__1, d__2);
1881c4762a1bSJed Brown       cstarr  = PetscSqrtScalar(*gamr * *pstar / *rstarr);
1882c4762a1bSJed Brown       ustarr  = *uxr + gascr4 * 2. * (cstarr - cr);
1883c4762a1bSJed Brown       zr      = *rstarr * cstarr;
1884c4762a1bSJed Brown       dpstar  = zl * zr * (ustarr - ustarl) / (zl + zr);
1885c4762a1bSJed Brown       *pstar -= dpstar;
1886c4762a1bSJed Brown       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
1887c4762a1bSJed Brown       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1888c4762a1bSJed Brown #if 0
1889c4762a1bSJed Brown         break;
1890c4762a1bSJed Brown #endif
1891c4762a1bSJed Brown       }
1892c4762a1bSJed Brown     }
1893c4762a1bSJed Brown     /*     1-wave: shock wave, 3-wave: rarefaction wave */
1894c4762a1bSJed Brown   } else if (iwave == 200) {
1895c4762a1bSJed Brown     /*     case (200) */
1896c4762a1bSJed Brown     i__1 = iterno;
1897c4762a1bSJed Brown     for (i0 = 1; i0 <= i__1; ++i0) {
1898c4762a1bSJed Brown       pst     = *pstar / *pl;
1899c4762a1bSJed Brown       ustarl  = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
1900c4762a1bSJed Brown       zl      = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
1901c4762a1bSJed Brown       d__1    = *pstar / *pr;
1902c4762a1bSJed Brown       d__2    = 1. / *gamr;
1903c4762a1bSJed Brown       *rstarr = *rr * PetscPowScalar(d__1, d__2);
1904c4762a1bSJed Brown       cstarr  = PetscSqrtScalar(*gamr * *pstar / *rstarr);
1905c4762a1bSJed Brown       zr      = *rstarr * cstarr;
1906c4762a1bSJed Brown       ustarr  = *uxr + gascr4 * 2. * (cstarr - cr);
1907c4762a1bSJed Brown       dpstar  = zl * zr * (ustarr - ustarl) / (zl + zr);
1908c4762a1bSJed Brown       *pstar -= dpstar;
1909c4762a1bSJed Brown       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
1910c4762a1bSJed Brown       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1911c4762a1bSJed Brown #if 0
1912c4762a1bSJed Brown         break;
1913c4762a1bSJed Brown #endif
1914c4762a1bSJed Brown       }
1915c4762a1bSJed Brown     }
1916c4762a1bSJed Brown     /*     1-wave: shock wave, 3-wave: shock */
1917c4762a1bSJed Brown   } else if (iwave == 300) {
1918c4762a1bSJed Brown     /*     case (300) */
1919c4762a1bSJed Brown     i__1 = iterno;
1920c4762a1bSJed Brown     for (i0 = 1; i0 <= i__1; ++i0) {
1921c4762a1bSJed Brown       pst    = *pstar / *pl;
1922c4762a1bSJed Brown       ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
1923c4762a1bSJed Brown       zl     = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
1924c4762a1bSJed Brown       pst    = *pstar / *pr;
1925c4762a1bSJed Brown       ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
1926c4762a1bSJed Brown       zr     = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
1927c4762a1bSJed Brown       dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
1928c4762a1bSJed Brown       *pstar -= dpstar;
1929c4762a1bSJed Brown       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
1930c4762a1bSJed Brown       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1931c4762a1bSJed Brown #if 0
1932c4762a1bSJed Brown         break;
1933c4762a1bSJed Brown #endif
1934c4762a1bSJed Brown       }
1935c4762a1bSJed Brown     }
1936c4762a1bSJed Brown     /*     1-wave: rarefaction wave, 3-wave: shock */
1937c4762a1bSJed Brown   } else if (iwave == 400) {
1938c4762a1bSJed Brown     /*     case (400) */
1939c4762a1bSJed Brown     i__1 = iterno;
1940c4762a1bSJed Brown     for (i0 = 1; i0 <= i__1; ++i0) {
1941c4762a1bSJed Brown       d__1    = *pstar / *pl;
1942c4762a1bSJed Brown       d__2    = 1. / *gaml;
1943c4762a1bSJed Brown       *rstarl = *rl * PetscPowScalar(d__1, d__2);
1944c4762a1bSJed Brown       cstarl  = PetscSqrtScalar(*gaml * *pstar / *rstarl);
1945c4762a1bSJed Brown       ustarl  = *uxl - gascl4 * 2. * (cstarl - cl);
1946c4762a1bSJed Brown       zl      = *rstarl * cstarl;
1947c4762a1bSJed Brown       pst     = *pstar / *pr;
1948c4762a1bSJed Brown       ustarr  = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
1949c4762a1bSJed Brown       zr      = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
1950c4762a1bSJed Brown       dpstar  = zl * zr * (ustarr - ustarl) / (zl + zr);
1951c4762a1bSJed Brown       *pstar -= dpstar;
1952c4762a1bSJed Brown       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
1953c4762a1bSJed Brown       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1954c4762a1bSJed Brown #if 0
1955c4762a1bSJed Brown               break;
1956c4762a1bSJed Brown #endif
1957c4762a1bSJed Brown       }
1958c4762a1bSJed Brown     }
1959c4762a1bSJed Brown   }
1960c4762a1bSJed Brown 
1961c4762a1bSJed Brown   *ustar = (zl * ustarr + zr * ustarl) / (zl + zr);
1962c4762a1bSJed Brown   if (PetscRealPart(*pstar) > PetscRealPart(*pl)) {
1963c4762a1bSJed Brown     pst     = *pstar / *pl;
19649371c9d4SSatish Balay     *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + *gaml + 1.) * *rl;
1965c4762a1bSJed Brown   }
1966c4762a1bSJed Brown   if (PetscRealPart(*pstar) > PetscRealPart(*pr)) {
1967c4762a1bSJed Brown     pst     = *pstar / *pr;
19689371c9d4SSatish Balay     *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + *gamr + 1.) * *rr;
1969c4762a1bSJed Brown   }
1970c4762a1bSJed Brown   return iwave;
1971c4762a1bSJed Brown }
1972c4762a1bSJed Brown 
1973d71ae5a4SJacob Faibussowitsch PetscScalar sign(PetscScalar x)
1974d71ae5a4SJacob Faibussowitsch {
1975c4762a1bSJed Brown   if (PetscRealPart(x) > 0) return 1.0;
1976c4762a1bSJed Brown   if (PetscRealPart(x) < 0) return -1.0;
1977c4762a1bSJed Brown   return 0.0;
1978c4762a1bSJed Brown }
1979c4762a1bSJed Brown /*        Riemann Solver */
1980c4762a1bSJed Brown /* -------------------------------------------------------------------- */
1981d71ae5a4SJacob Faibussowitsch 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)
1982d71ae5a4SJacob Faibussowitsch {
1983c4762a1bSJed Brown   /* System generated locals */
1984c4762a1bSJed Brown   PetscScalar d__1, d__2;
1985c4762a1bSJed Brown 
1986c4762a1bSJed Brown   /* Local variables */
1987c4762a1bSJed Brown   static PetscScalar s, c0, p0, r0, u0, w0, x0, x2, ri, cx, sgn0, wsp0, gasc1, gasc2, gasc3, gasc4;
1988c4762a1bSJed Brown   static PetscScalar cstar, pstar, rstar, ustar, xstar, wspst, ushock, streng, rstarl, rstarr, rstars;
1989c4762a1bSJed Brown   int                iwave;
1990c4762a1bSJed Brown 
1991c4762a1bSJed Brown   if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) {
1992c4762a1bSJed Brown     *rx  = *rl;
1993c4762a1bSJed Brown     *px  = *pl;
1994c4762a1bSJed Brown     *uxm = *uxl;
1995c4762a1bSJed Brown     *gam = *gaml;
1996c4762a1bSJed Brown     x2   = *xcen + *uxm * *dtt;
1997c4762a1bSJed Brown 
1998c4762a1bSJed Brown     if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
1999c4762a1bSJed Brown       *utx  = *utr;
2000c4762a1bSJed Brown       *ubx  = *ubr;
2001c4762a1bSJed Brown       *rho1 = *rho1r;
2002c4762a1bSJed Brown     } else {
2003c4762a1bSJed Brown       *utx  = *utl;
2004c4762a1bSJed Brown       *ubx  = *ubl;
2005c4762a1bSJed Brown       *rho1 = *rho1l;
2006c4762a1bSJed Brown     }
2007c4762a1bSJed Brown     return 0;
2008c4762a1bSJed Brown   }
2009c4762a1bSJed Brown   iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar);
2010c4762a1bSJed Brown 
2011c4762a1bSJed Brown   x2   = *xcen + ustar * *dtt;
2012c4762a1bSJed Brown   d__1 = *xp - x2;
2013c4762a1bSJed Brown   sgn0 = sign(d__1);
2014c4762a1bSJed Brown   /*            x is in 3-wave if sgn0 = 1 */
2015c4762a1bSJed Brown   /*            x is in 1-wave if sgn0 = -1 */
2016c4762a1bSJed Brown   r0     = cvmgm_(rl, rr, &sgn0);
2017c4762a1bSJed Brown   p0     = cvmgm_(pl, pr, &sgn0);
2018c4762a1bSJed Brown   u0     = cvmgm_(uxl, uxr, &sgn0);
2019c4762a1bSJed Brown   *gam   = cvmgm_(gaml, gamr, &sgn0);
2020c4762a1bSJed Brown   gasc1  = *gam - 1.;
2021c4762a1bSJed Brown   gasc2  = (*gam + 1.) * .5;
2022c4762a1bSJed Brown   gasc3  = gasc2 / *gam;
2023c4762a1bSJed Brown   gasc4  = 1. / (*gam - 1.);
2024c4762a1bSJed Brown   c0     = PetscSqrtScalar(*gam * p0 / r0);
2025c4762a1bSJed Brown   streng = pstar - p0;
2026c4762a1bSJed Brown   w0     = *gam * r0 * p0 * (gasc3 * streng / p0 + 1.);
2027c4762a1bSJed Brown   rstars = r0 / (1. - r0 * streng / w0);
2028c4762a1bSJed Brown   d__1   = p0 / pstar;
2029c4762a1bSJed Brown   d__2   = -1. / *gam;
2030c4762a1bSJed Brown   rstarr = r0 * PetscPowScalar(d__1, d__2);
2031c4762a1bSJed Brown   rstar  = cvmgm_(&rstarr, &rstars, &streng);
2032c4762a1bSJed Brown   w0     = PetscSqrtScalar(w0);
2033c4762a1bSJed Brown   cstar  = PetscSqrtScalar(*gam * pstar / rstar);
2034c4762a1bSJed Brown   wsp0   = u0 + sgn0 * c0;
2035c4762a1bSJed Brown   wspst  = ustar + sgn0 * cstar;
2036c4762a1bSJed Brown   ushock = ustar + sgn0 * w0 / rstar;
2037c4762a1bSJed Brown   wspst  = cvmgp_(&ushock, &wspst, &streng);
2038c4762a1bSJed Brown   wsp0   = cvmgp_(&ushock, &wsp0, &streng);
2039c4762a1bSJed Brown   x0     = *xcen + wsp0 * *dtt;
2040c4762a1bSJed Brown   xstar  = *xcen + wspst * *dtt;
2041c4762a1bSJed Brown   /*           using gas formula to evaluate rarefaction wave */
2042c4762a1bSJed Brown   /*            ri : reiman invariant */
2043c4762a1bSJed Brown   ri   = u0 - sgn0 * 2. * gasc4 * c0;
2044c4762a1bSJed Brown   cx   = sgn0 * .5 * gasc1 / gasc2 * ((*xp - *xcen) / *dtt - ri);
2045c4762a1bSJed Brown   *uxm = ri + sgn0 * 2. * gasc4 * cx;
2046c4762a1bSJed Brown   s    = p0 / PetscPowScalar(r0, *gam);
2047c4762a1bSJed Brown   d__1 = cx * cx / (*gam * s);
2048c4762a1bSJed Brown   *rx  = PetscPowScalar(d__1, gasc4);
2049c4762a1bSJed Brown   *px  = cx * cx * *rx / *gam;
2050c4762a1bSJed Brown   d__1 = sgn0 * (x0 - *xp);
2051c4762a1bSJed Brown   *rx  = cvmgp_(rx, &r0, &d__1);
2052c4762a1bSJed Brown   d__1 = sgn0 * (x0 - *xp);
2053c4762a1bSJed Brown   *px  = cvmgp_(px, &p0, &d__1);
2054c4762a1bSJed Brown   d__1 = sgn0 * (x0 - *xp);
2055c4762a1bSJed Brown   *uxm = cvmgp_(uxm, &u0, &d__1);
2056c4762a1bSJed Brown   d__1 = sgn0 * (xstar - *xp);
2057c4762a1bSJed Brown   *rx  = cvmgm_(rx, &rstar, &d__1);
2058c4762a1bSJed Brown   d__1 = sgn0 * (xstar - *xp);
2059c4762a1bSJed Brown   *px  = cvmgm_(px, &pstar, &d__1);
2060c4762a1bSJed Brown   d__1 = sgn0 * (xstar - *xp);
2061c4762a1bSJed Brown   *uxm = cvmgm_(uxm, &ustar, &d__1);
2062c4762a1bSJed Brown   if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
2063c4762a1bSJed Brown     *utx  = *utr;
2064c4762a1bSJed Brown     *ubx  = *ubr;
2065c4762a1bSJed Brown     *rho1 = *rho1r;
2066c4762a1bSJed Brown   } else {
2067c4762a1bSJed Brown     *utx  = *utl;
2068c4762a1bSJed Brown     *ubx  = *ubl;
2069c4762a1bSJed Brown     *rho1 = *rho1l;
2070c4762a1bSJed Brown   }
2071c4762a1bSJed Brown   return iwave;
2072c4762a1bSJed Brown }
2073d71ae5a4SJacob Faibussowitsch int godunovflux(const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, const int *ndim, const PetscReal *gamma)
2074d71ae5a4SJacob Faibussowitsch {
2075c4762a1bSJed Brown   /* System generated locals */
2076c4762a1bSJed Brown   int         i__1, iwave;
2077c4762a1bSJed Brown   PetscScalar d__1, d__2, d__3;
2078c4762a1bSJed Brown 
2079c4762a1bSJed Brown   /* Local variables */
2080c4762a1bSJed Brown   static int         k;
20819371c9d4SSatish Balay   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;
2082c4762a1bSJed Brown 
2083c4762a1bSJed Brown   /* Function Body */
2084c4762a1bSJed Brown   xcen = 0.;
2085c4762a1bSJed Brown   xp   = 0.;
2086c4762a1bSJed Brown   i__1 = *ndim;
2087c4762a1bSJed Brown   for (k = 1; k <= i__1; ++k) {
2088c4762a1bSJed Brown     tg[k - 1] = 0.;
2089c4762a1bSJed Brown     bn[k - 1] = 0.;
2090c4762a1bSJed Brown   }
2091c4762a1bSJed Brown   dtt = 1.;
2092c4762a1bSJed Brown   if (*ndim == 3) {
2093a795f3e7SBarry Smith     if (nn[0] == 0. && nn[1] == 0.) {
2094c4762a1bSJed Brown       tg[0] = 1.;
2095c4762a1bSJed Brown     } else {
2096a795f3e7SBarry Smith       tg[0] = -nn[1];
2097a795f3e7SBarry Smith       tg[1] = nn[0];
2098c4762a1bSJed Brown     }
2099c4762a1bSJed Brown     /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2100c4762a1bSJed Brown     /*           tg=tg/tmp */
2101a795f3e7SBarry Smith     bn[0] = -nn[2] * tg[1];
2102a795f3e7SBarry Smith     bn[1] = nn[2] * tg[0];
2103a795f3e7SBarry Smith     bn[2] = nn[0] * tg[1] - nn[1] * tg[0];
2104c4762a1bSJed Brown     /* Computing 2nd power */
2105c4762a1bSJed Brown     d__1 = bn[0];
2106c4762a1bSJed Brown     /* Computing 2nd power */
2107c4762a1bSJed Brown     d__2 = bn[1];
2108c4762a1bSJed Brown     /* Computing 2nd power */
2109c4762a1bSJed Brown     d__3 = bn[2];
2110c4762a1bSJed Brown     tmp  = PetscSqrtScalar(d__1 * d__1 + d__2 * d__2 + d__3 * d__3);
2111c4762a1bSJed Brown     i__1 = *ndim;
2112ad540459SPierre Jolivet     for (k = 1; k <= i__1; ++k) bn[k - 1] /= tmp;
2113c4762a1bSJed Brown   } else if (*ndim == 2) {
2114a795f3e7SBarry Smith     tg[0] = -nn[1];
2115a795f3e7SBarry Smith     tg[1] = nn[0];
2116c4762a1bSJed Brown     /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2117c4762a1bSJed Brown     /*           tg=tg/tmp */
2118c4762a1bSJed Brown     bn[0] = 0.;
2119c4762a1bSJed Brown     bn[1] = 0.;
2120c4762a1bSJed Brown     bn[2] = 1.;
2121c4762a1bSJed Brown   }
2122a795f3e7SBarry Smith   rl   = ul[0];
2123a795f3e7SBarry Smith   rr   = ur[0];
2124c4762a1bSJed Brown   uxl  = 0.;
2125c4762a1bSJed Brown   uxr  = 0.;
2126c4762a1bSJed Brown   utl  = 0.;
2127c4762a1bSJed Brown   utr  = 0.;
2128c4762a1bSJed Brown   ubl  = 0.;
2129c4762a1bSJed Brown   ubr  = 0.;
2130c4762a1bSJed Brown   i__1 = *ndim;
2131c4762a1bSJed Brown   for (k = 1; k <= i__1; ++k) {
2132a795f3e7SBarry Smith     uxl += ul[k] * nn[k - 1];
2133a795f3e7SBarry Smith     uxr += ur[k] * nn[k - 1];
2134a795f3e7SBarry Smith     utl += ul[k] * tg[k - 1];
2135a795f3e7SBarry Smith     utr += ur[k] * tg[k - 1];
2136a795f3e7SBarry Smith     ubl += ul[k] * bn[k - 1];
2137a795f3e7SBarry Smith     ubr += ur[k] * bn[k - 1];
2138c4762a1bSJed Brown   }
2139c4762a1bSJed Brown   uxl /= rl;
2140c4762a1bSJed Brown   uxr /= rr;
2141c4762a1bSJed Brown   utl /= rl;
2142c4762a1bSJed Brown   utr /= rr;
2143c4762a1bSJed Brown   ubl /= rl;
2144c4762a1bSJed Brown   ubr /= rr;
2145c4762a1bSJed Brown 
2146c4762a1bSJed Brown   gaml = *gamma;
2147c4762a1bSJed Brown   gamr = *gamma;
2148c4762a1bSJed Brown   /* Computing 2nd power */
2149c4762a1bSJed Brown   d__1 = uxl;
2150c4762a1bSJed Brown   /* Computing 2nd power */
2151c4762a1bSJed Brown   d__2 = utl;
2152c4762a1bSJed Brown   /* Computing 2nd power */
2153c4762a1bSJed Brown   d__3 = ubl;
2154a795f3e7SBarry Smith   pl   = (*gamma - 1.) * (ul[*ndim + 1] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2155c4762a1bSJed Brown   /* Computing 2nd power */
2156c4762a1bSJed Brown   d__1 = uxr;
2157c4762a1bSJed Brown   /* Computing 2nd power */
2158c4762a1bSJed Brown   d__2 = utr;
2159c4762a1bSJed Brown   /* Computing 2nd power */
2160c4762a1bSJed Brown   d__3  = ubr;
2161a795f3e7SBarry Smith   pr    = (*gamma - 1.) * (ur[*ndim + 1] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2162c4762a1bSJed Brown   rho1l = rl;
2163c4762a1bSJed Brown   rho1r = rr;
2164c4762a1bSJed Brown 
21659371c9d4SSatish Balay   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);
2166c4762a1bSJed Brown 
2167a795f3e7SBarry Smith   flux[0] = rhom * unm;
2168c4762a1bSJed Brown   fn      = rhom * unm * unm + pm;
2169c4762a1bSJed Brown   ft      = rhom * unm * utm;
2170c4762a1bSJed Brown   /*           flux(2)=fn*nn(1)+ft*nn(2) */
2171c4762a1bSJed Brown   /*           flux(3)=fn*tg(1)+ft*tg(2) */
2172a795f3e7SBarry Smith   flux[1] = fn * nn[0] + ft * tg[0];
2173a795f3e7SBarry Smith   flux[2] = fn * nn[1] + ft * tg[1];
2174c4762a1bSJed Brown   /*           flux(2)=rhom*unm*(unm)+pm */
2175c4762a1bSJed Brown   /*           flux(3)=rhom*(unm)*utm */
2176ad540459SPierre Jolivet   if (*ndim == 3) flux[3] = rhom * unm * ubm;
2177a795f3e7SBarry Smith   flux[*ndim + 1] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm;
2178c4762a1bSJed Brown   return iwave;
2179c4762a1bSJed Brown } /* godunovflux_ */
2180c4762a1bSJed Brown 
2181c4762a1bSJed Brown /* Subroutine to set up the initial conditions for the */
2182c4762a1bSJed Brown /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */
2183c4762a1bSJed Brown /* ----------------------------------------------------------------------- */
2184d71ae5a4SJacob Faibussowitsch int projecteqstate(PetscReal wc[], const PetscReal ueq[], PetscReal lv[][3])
2185d71ae5a4SJacob Faibussowitsch {
2186c4762a1bSJed Brown   int j, k;
2187c4762a1bSJed Brown   /*      Wc=matmul(lv,Ueq) 3 vars */
2188c4762a1bSJed Brown   for (k = 0; k < 3; ++k) {
2189c4762a1bSJed Brown     wc[k] = 0.;
2190ad540459SPierre Jolivet     for (j = 0; j < 3; ++j) wc[k] += lv[k][j] * ueq[j];
2191c4762a1bSJed Brown   }
2192c4762a1bSJed Brown   return 0;
2193c4762a1bSJed Brown }
2194c4762a1bSJed Brown /* ----------------------------------------------------------------------- */
2195d71ae5a4SJacob Faibussowitsch int projecttoprim(PetscReal v[], const PetscReal wc[], PetscReal rv[][3])
2196d71ae5a4SJacob Faibussowitsch {
2197c4762a1bSJed Brown   int k, j;
2198c4762a1bSJed Brown   /*      V=matmul(rv,WC) 3 vars */
2199c4762a1bSJed Brown   for (k = 0; k < 3; ++k) {
2200c4762a1bSJed Brown     v[k] = 0.;
2201ad540459SPierre Jolivet     for (j = 0; j < 3; ++j) v[k] += rv[k][j] * wc[j];
2202c4762a1bSJed Brown   }
2203c4762a1bSJed Brown   return 0;
2204c4762a1bSJed Brown }
2205c4762a1bSJed Brown /* ---------------------------------------------------------------------- */
2206d71ae5a4SJacob Faibussowitsch int eigenvectors(PetscReal rv[][3], PetscReal lv[][3], const PetscReal ueq[], PetscReal gamma)
2207d71ae5a4SJacob Faibussowitsch {
2208c4762a1bSJed Brown   int       j, k;
2209c4762a1bSJed Brown   PetscReal rho, csnd, p0;
2210c4762a1bSJed Brown   /* PetscScalar u; */
2211c4762a1bSJed Brown 
22129371c9d4SSatish Balay   for (k = 0; k < 3; ++k)
22139371c9d4SSatish Balay     for (j = 0; j < 3; ++j) {
22149371c9d4SSatish Balay       lv[k][j] = 0.;
22159371c9d4SSatish Balay       rv[k][j] = 0.;
22169371c9d4SSatish Balay     }
2217c4762a1bSJed Brown   rho = ueq[0];
2218c4762a1bSJed Brown   /* u = ueq[1]; */
2219c4762a1bSJed Brown   p0       = ueq[2];
2220c4762a1bSJed Brown   csnd     = PetscSqrtReal(gamma * p0 / rho);
2221c4762a1bSJed Brown   lv[0][1] = rho * .5;
2222c4762a1bSJed Brown   lv[0][2] = -.5 / csnd;
2223c4762a1bSJed Brown   lv[1][0] = csnd;
2224c4762a1bSJed Brown   lv[1][2] = -1. / csnd;
2225c4762a1bSJed Brown   lv[2][1] = rho * .5;
2226c4762a1bSJed Brown   lv[2][2] = .5 / csnd;
2227c4762a1bSJed Brown   rv[0][0] = -1. / csnd;
2228c4762a1bSJed Brown   rv[1][0] = 1. / rho;
2229c4762a1bSJed Brown   rv[2][0] = -csnd;
2230c4762a1bSJed Brown   rv[0][1] = 1. / csnd;
2231c4762a1bSJed Brown   rv[0][2] = 1. / csnd;
2232c4762a1bSJed Brown   rv[1][2] = 1. / rho;
2233c4762a1bSJed Brown   rv[2][2] = csnd;
2234c4762a1bSJed Brown   return 0;
2235c4762a1bSJed Brown }
2236c4762a1bSJed Brown 
2237d71ae5a4SJacob Faibussowitsch int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx)
2238d71ae5a4SJacob Faibussowitsch {
2239c4762a1bSJed Brown   PetscReal p0, u0, wcp[3], wc[3];
2240c4762a1bSJed Brown   PetscReal lv[3][3];
2241c4762a1bSJed Brown   PetscReal vp[3];
2242c4762a1bSJed Brown   PetscReal rv[3][3];
2243c4762a1bSJed Brown   PetscReal eps, ueq[3], rho0, twopi;
2244c4762a1bSJed Brown 
2245c4762a1bSJed Brown   /* Function Body */
2246c4762a1bSJed Brown   twopi  = 2. * PETSC_PI;
2247c4762a1bSJed Brown   eps    = 1e-4;    /* perturbation */
2248c4762a1bSJed Brown   rho0   = 1e3;     /* density of water */
2249c4762a1bSJed Brown   p0     = 101325.; /* init pressure of 1 atm (?) */
2250c4762a1bSJed Brown   u0     = 0.;
2251c4762a1bSJed Brown   ueq[0] = rho0;
2252c4762a1bSJed Brown   ueq[1] = u0;
2253c4762a1bSJed Brown   ueq[2] = p0;
2254c4762a1bSJed Brown   /* Project initial state to characteristic variables */
2255c4762a1bSJed Brown   eigenvectors(rv, lv, ueq, gamma);
2256c4762a1bSJed Brown   projecteqstate(wc, ueq, lv);
2257c4762a1bSJed Brown   wcp[0] = wc[0];
2258c4762a1bSJed Brown   wcp[1] = wc[1];
2259c4762a1bSJed Brown   wcp[2] = wc[2] + eps * PetscCosReal(coord[0] * 2. * twopi / Lx);
2260c4762a1bSJed Brown   projecttoprim(vp, wcp, rv);
2261c4762a1bSJed Brown   ux->r     = vp[0];         /* density */
2262c4762a1bSJed Brown   ux->ru[0] = vp[0] * vp[1]; /* x momentum */
2263c4762a1bSJed Brown   ux->ru[1] = 0.;
2264c4762a1bSJed Brown #if defined DIM > 2
2265c4762a1bSJed Brown   if (dim > 2) ux->ru[2] = 0.;
2266c4762a1bSJed Brown #endif
2267c4762a1bSJed Brown   /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */
2268c4762a1bSJed Brown   ux->E = vp[2] / (gamma - 1.) + 0.5 * vp[0] * vp[1] * vp[1];
2269c4762a1bSJed Brown   return 0;
2270c4762a1bSJed Brown }
2271c4762a1bSJed Brown 
2272c4762a1bSJed Brown /*TEST
2273c4762a1bSJed Brown 
227430602db0SMatthew G. Knepley   testset:
227530602db0SMatthew G. Knepley     args: -dm_plex_adj_cone -dm_plex_adj_closure 0
227630602db0SMatthew G. Knepley 
227730602db0SMatthew G. Knepley     test:
227830602db0SMatthew G. Knepley       suffix: adv_2d_tri_0
227930602db0SMatthew G. Knepley       requires: triangle
228030602db0SMatthew G. Knepley       TODO: how did this ever get in main when there is no support for this
228130602db0SMatthew G. Knepley       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
228230602db0SMatthew G. Knepley 
228330602db0SMatthew G. Knepley     test:
228430602db0SMatthew G. Knepley       suffix: adv_2d_tri_1
228530602db0SMatthew G. Knepley       requires: triangle
228630602db0SMatthew G. Knepley       TODO: how did this ever get in main when there is no support for this
228730602db0SMatthew G. Knepley       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
228830602db0SMatthew G. Knepley 
228930602db0SMatthew G. Knepley     test:
229030602db0SMatthew G. Knepley       suffix: tut_1
229130602db0SMatthew G. Knepley       requires: exodusii
229230602db0SMatthew G. Knepley       nsize: 1
229330602db0SMatthew G. Knepley       args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
229430602db0SMatthew G. Knepley 
229530602db0SMatthew G. Knepley     test:
229630602db0SMatthew G. Knepley       suffix: tut_2
229730602db0SMatthew G. Knepley       requires: exodusii
229830602db0SMatthew G. Knepley       nsize: 1
229930602db0SMatthew G. Knepley       args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw
230030602db0SMatthew G. Knepley 
230130602db0SMatthew G. Knepley     test:
230230602db0SMatthew G. Knepley       suffix: tut_3
230330602db0SMatthew G. Knepley       requires: exodusii
230430602db0SMatthew G. Knepley       nsize: 4
2305e600fa54SMatthew G. Knepley       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
230630602db0SMatthew G. Knepley 
230730602db0SMatthew G. Knepley     test:
230830602db0SMatthew G. Knepley       suffix: tut_4
230930602db0SMatthew G. Knepley       requires: exodusii
231030602db0SMatthew G. Knepley       nsize: 4
2311e600fa54SMatthew G. Knepley       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
231230602db0SMatthew G. Knepley 
231330602db0SMatthew G. Knepley   testset:
231430602db0SMatthew G. Knepley     args: -dm_plex_adj_cone -dm_plex_adj_closure 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1
231530602db0SMatthew G. Knepley 
2316c4762a1bSJed Brown     # 2D Advection 0-10
2317c4762a1bSJed Brown     test:
2318c4762a1bSJed Brown       suffix: 0
2319c4762a1bSJed Brown       requires: exodusii
232030602db0SMatthew G. Knepley       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2321c4762a1bSJed Brown 
2322c4762a1bSJed Brown     test:
2323c4762a1bSJed Brown       suffix: 1
2324c4762a1bSJed Brown       requires: exodusii
232530602db0SMatthew G. Knepley       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
2326c4762a1bSJed Brown 
2327c4762a1bSJed Brown     test:
2328c4762a1bSJed Brown       suffix: 2
2329c4762a1bSJed Brown       requires: exodusii
2330c4762a1bSJed Brown       nsize: 2
2331e600fa54SMatthew G. Knepley       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2332c4762a1bSJed Brown 
2333c4762a1bSJed Brown     test:
2334c4762a1bSJed Brown       suffix: 3
2335c4762a1bSJed Brown       requires: exodusii
2336c4762a1bSJed Brown       nsize: 2
2337e600fa54SMatthew G. Knepley       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
2338c4762a1bSJed Brown 
2339c4762a1bSJed Brown     test:
2340c4762a1bSJed Brown       suffix: 4
2341c4762a1bSJed Brown       requires: exodusii
2342*6c2b77d5SStefano Zampini       nsize: 4
2343*6c2b77d5SStefano Zampini       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -petscpartitioner_type simple
2344c4762a1bSJed Brown 
2345c4762a1bSJed Brown     test:
2346c4762a1bSJed Brown       suffix: 5
2347c4762a1bSJed Brown       requires: exodusii
234830602db0SMatthew G. Knepley       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1
2349c4762a1bSJed Brown 
2350c4762a1bSJed Brown     test:
2351c4762a1bSJed Brown       suffix: 7
2352c4762a1bSJed Brown       requires: exodusii
235330602db0SMatthew G. Knepley       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
2354c4762a1bSJed Brown 
2355c4762a1bSJed Brown     test:
2356c4762a1bSJed Brown       suffix: 8
2357c4762a1bSJed Brown       requires: exodusii
2358c4762a1bSJed Brown       nsize: 2
2359e600fa54SMatthew G. Knepley       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
2360c4762a1bSJed Brown 
2361c4762a1bSJed Brown     test:
2362c4762a1bSJed Brown       suffix: 9
2363c4762a1bSJed Brown       requires: exodusii
2364c4762a1bSJed Brown       nsize: 8
2365e600fa54SMatthew G. Knepley       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
2366c4762a1bSJed Brown 
2367c4762a1bSJed Brown     test:
2368c4762a1bSJed Brown       suffix: 10
2369c4762a1bSJed Brown       requires: exodusii
237030602db0SMatthew G. Knepley       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo
2371c4762a1bSJed Brown 
2372c4762a1bSJed Brown   # 2D Shallow water
2373993a5519SMatthew G. Knepley   testset:
2374993a5519SMatthew G. Knepley     args: -physics sw -ufv_vtk_interval 0 -dm_plex_adj_cone -dm_plex_adj_closure 0
2375993a5519SMatthew G. Knepley 
2376c4762a1bSJed Brown     test:
2377c4762a1bSJed Brown       suffix: sw_0
2378c4762a1bSJed Brown       requires: exodusii
2379993a5519SMatthew G. Knepley       args: -bc_wall 100,101 -ufv_cfl 5 -petscfv_type leastsquares -petsclimiter_type sin \
2380993a5519SMatthew G. Knepley             -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo \
2381993a5519SMatthew G. Knepley             -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
2382993a5519SMatthew G. Knepley             -monitor height,energy
2383993a5519SMatthew G. Knepley 
2384993a5519SMatthew G. Knepley     test:
2385993a5519SMatthew G. Knepley       suffix: sw_1
2386993a5519SMatthew G. Knepley       nsize: 2
2387993a5519SMatthew G. Knepley       args: -bc_wall 1,3 -ufv_cfl 5 -petsclimiter_type sin \
2388993a5519SMatthew G. Knepley             -dm_plex_shape annulus -dm_plex_simplex 0 -dm_plex_box_faces 24,12 -dm_plex_box_lower 0,1 -dm_plex_box_upper 1,3 -dm_distribute_overlap 1 \
2389993a5519SMatthew G. Knepley             -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
2390993a5519SMatthew G. Knepley             -monitor height,energy
2391c4762a1bSJed Brown 
23928d2c18caSMukkund Sunjii     test:
23938d2c18caSMukkund Sunjii       suffix: sw_hll
2394993a5519SMatthew G. Knepley       args: -sw_riemann hll -bc_wall 1,2,3,4 -ufv_cfl 3 -petscfv_type leastsquares -petsclimiter_type sin \
2395993a5519SMatthew G. Knepley             -grid_bounds 0,5,0,5 -dm_plex_simplex 0 -dm_plex_box_faces 25,25 \
2396993a5519SMatthew G. Knepley             -ts_max_steps 5 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
2397993a5519SMatthew G. Knepley             -monitor height,energy
2398993a5519SMatthew G. Knepley 
2399993a5519SMatthew G. Knepley   testset:
2400993a5519SMatthew G. Knepley     args: -dm_plex_adj_cone -dm_plex_adj_closure 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1
24018d2c18caSMukkund Sunjii 
2402c4762a1bSJed Brown     # 2D Advection: p4est
2403c4762a1bSJed Brown     test:
2404c4762a1bSJed Brown       suffix: p4est_advec_2d
2405c4762a1bSJed Brown       requires: p4est
240630602db0SMatthew G. Knepley       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
2407c4762a1bSJed Brown 
2408c4762a1bSJed Brown     # Advection in a box
2409c4762a1bSJed Brown     test:
2410c4762a1bSJed Brown       suffix: adv_2d_quad_0
2411c4762a1bSJed Brown       args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2412c4762a1bSJed Brown 
2413c4762a1bSJed Brown     test:
2414c4762a1bSJed Brown       suffix: adv_2d_quad_1
2415c4762a1bSJed Brown       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
2416c4762a1bSJed Brown       timeoutfactor: 3
2417c4762a1bSJed Brown 
2418c4762a1bSJed Brown     test:
2419c4762a1bSJed Brown       suffix: adv_2d_quad_p4est_0
2420c4762a1bSJed Brown       requires: p4est
2421c4762a1bSJed Brown       args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2422c4762a1bSJed Brown 
2423c4762a1bSJed Brown     test:
2424c4762a1bSJed Brown       suffix: adv_2d_quad_p4est_1
2425c4762a1bSJed Brown       requires: p4est
2426c4762a1bSJed Brown       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
2427c4762a1bSJed Brown       timeoutfactor: 3
2428c4762a1bSJed Brown 
2429c4762a1bSJed Brown     test:
2430c4762a1bSJed Brown       suffix: adv_2d_quad_p4est_adapt_0
2431c4762a1bSJed Brown       requires: p4est !__float128 #broken for quad precision
2432c4762a1bSJed Brown       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
2433c4762a1bSJed Brown       timeoutfactor: 3
2434c4762a1bSJed Brown 
2435c4762a1bSJed Brown     test:
2436c4762a1bSJed Brown       suffix: adv_0
2437c4762a1bSJed Brown       requires: exodusii
243830602db0SMatthew G. Knepley       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
2439c4762a1bSJed Brown 
2440c4762a1bSJed Brown     test:
2441c4762a1bSJed Brown       suffix: shock_0
2442c4762a1bSJed Brown       requires: p4est !single !complex
244330602db0SMatthew G. Knepley       args: -dm_plex_box_faces 2,1 -grid_bounds -1,1.,0.,1 -grid_skew_60 \
244430602db0SMatthew G. Knepley       -dm_type p4est -dm_forest_partition_overlap 1 -dm_forest_maximum_refinement 6 -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 \
244530602db0SMatthew G. Knepley       -ufv_use_amr -refine_vec_tagger_box 0.5,inf -coarsen_vec_tagger_box 0,1.e-2 -refine_tag_view -coarsen_tag_view \
244630602db0SMatthew G. Knepley       -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. \
244730602db0SMatthew G. Knepley       -petscfv_type leastsquares -petsclimiter_type minmod -petscfv_compute_gradients 0 \
244830602db0SMatthew G. Knepley       -ts_max_time 0.5 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
244930602db0SMatthew G. Knepley       -ufv_vtk_basename ${wPETSC_DIR}/ex11 -ufv_vtk_interval 0 -monitor density,energy
2450c4762a1bSJed Brown       timeoutfactor: 3
2451c4762a1bSJed Brown 
2452c4762a1bSJed Brown     # Test GLVis visualization of PetscFV fields
2453c4762a1bSJed Brown     test:
2454c4762a1bSJed Brown       suffix: glvis_adv_2d_tet
245530602db0SMatthew G. Knepley       args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 \
245630602db0SMatthew G. Knepley             -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0 \
245730602db0SMatthew G. Knepley             -ts_monitor_solution glvis: -ts_max_steps 0
2458c4762a1bSJed Brown 
2459c4762a1bSJed Brown     test:
2460c4762a1bSJed Brown       suffix: glvis_adv_2d_quad
246130602db0SMatthew G. Knepley       args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 -bc_inflow 1,2,4 -bc_outflow 3 \
246230602db0SMatthew G. Knepley             -dm_refine 5 -dm_plex_separate_marker \
246330602db0SMatthew G. Knepley             -ts_monitor_solution glvis: -ts_max_steps 0
2464c4762a1bSJed Brown 
2465c4762a1bSJed Brown TEST*/
2466