xref: /petsc/src/ts/tutorials/ex11.h (revision a7d931c6983080384301a9270cd0ca2d035918fb)
1de8de29fSMatthew G. Knepley #include <petscdm.h>
2de8de29fSMatthew G. Knepley #include <petscdmceed.h>
3de8de29fSMatthew G. Knepley 
4*a7d931c6SMatthew G. Knepley #ifdef __CUDACC_RTC__
5*a7d931c6SMatthew G. Knepley   #define PETSC_HAVE_LIBCEED
6*a7d931c6SMatthew G. Knepley // Define Petsc types to be equal to Ceed types
7*a7d931c6SMatthew G. Knepley typedef CeedInt    PetscInt;
8*a7d931c6SMatthew G. Knepley typedef CeedScalar PetscReal;
9*a7d931c6SMatthew G. Knepley typedef CeedScalar PetscScalar;
10*a7d931c6SMatthew G. Knepley typedef CeedInt    PetscErrorCode;
11*a7d931c6SMatthew G. Knepley   // Define things we are missing from PETSc headers
12*a7d931c6SMatthew G. Knepley   #undef PETSC_SUCCESS
13*a7d931c6SMatthew G. Knepley   #define PETSC_SUCCESS   0
14*a7d931c6SMatthew G. Knepley   #define PETSC_COMM_SELF MPI_COMM_SELF
15*a7d931c6SMatthew G. Knepley   #undef PetscFunctionBeginUser
16*a7d931c6SMatthew G. Knepley   #define PetscFunctionBeginUser
17*a7d931c6SMatthew G. Knepley   #undef PetscFunctionReturn
18*a7d931c6SMatthew G. Knepley   #define PetscFunctionReturn(x) return x
19*a7d931c6SMatthew G. Knepley   #undef PetscCall
20*a7d931c6SMatthew G. Knepley   #define PetscCall(a)              a
21*a7d931c6SMatthew G. Knepley   #define PetscFunctionReturnVoid() return
22*a7d931c6SMatthew G. Knepley   //   Math definitions
23*a7d931c6SMatthew G. Knepley   #undef PetscSqrtReal
24*a7d931c6SMatthew G. Knepley   #define PetscSqrtReal(x) sqrt(x)
25*a7d931c6SMatthew G. Knepley   #undef PetscSqrtScalar
26*a7d931c6SMatthew G. Knepley   #define PetscSqrtScalar(x) sqrt(x)
27*a7d931c6SMatthew G. Knepley   #undef PetscSqr
28*a7d931c6SMatthew G. Knepley   #define PetscSqr(x)          (x * x)
29*a7d931c6SMatthew G. Knepley   #define PetscSqrReal(x)      (x * x)
30*a7d931c6SMatthew G. Knepley   #define PetscAbsReal(x)      abs(x)
31*a7d931c6SMatthew G. Knepley   #define PetscAbsScalar(x)    abs(x)
32*a7d931c6SMatthew G. Knepley   #define PetscMax(x, y)       x > y ? x : y
33*a7d931c6SMatthew G. Knepley   #define PetscMin(x, y)       x < y ? x : y
34*a7d931c6SMatthew G. Knepley   #define PetscRealPart(a)     a
35*a7d931c6SMatthew G. Knepley   #define PetscPowScalar(a, b) pow(a, b)
36*a7d931c6SMatthew G. Knepley #endif
37*a7d931c6SMatthew G. Knepley 
38de8de29fSMatthew G. Knepley #define DIM 2 /* Geometric dimension */
39de8de29fSMatthew G. Knepley 
40de8de29fSMatthew G. Knepley /* Represents continuum physical equations. */
41de8de29fSMatthew G. Knepley typedef struct _n_Physics *Physics;
42de8de29fSMatthew G. Knepley 
43de8de29fSMatthew G. Knepley /* Physical model includes boundary conditions, initial conditions, and functionals of interest. It is
44de8de29fSMatthew G. Knepley  * discretization-independent, but its members depend on the scenario being solved. */
45de8de29fSMatthew G. Knepley typedef struct _n_Model *Model;
46de8de29fSMatthew G. Knepley 
47de8de29fSMatthew G. Knepley struct FieldDescription {
48de8de29fSMatthew G. Knepley   const char *name;
49de8de29fSMatthew G. Knepley   PetscInt    dof;
50de8de29fSMatthew G. Knepley };
51de8de29fSMatthew G. Knepley 
52de8de29fSMatthew G. Knepley struct _n_Physics {
53de8de29fSMatthew G. Knepley   void (*riemann)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *);
54de8de29fSMatthew G. Knepley   PetscInt                       dof;      /* number of degrees of freedom per cell */
55de8de29fSMatthew G. Knepley   PetscReal                      maxspeed; /* kludge to pick initial time step, need to add monitoring and step control */
56de8de29fSMatthew G. Knepley   void                          *data;
57de8de29fSMatthew G. Knepley   PetscInt                       nfields;
58de8de29fSMatthew G. Knepley   const struct FieldDescription *field_desc;
59de8de29fSMatthew G. Knepley };
60de8de29fSMatthew G. Knepley 
61de8de29fSMatthew G. Knepley typedef struct {
62de8de29fSMatthew G. Knepley   PetscReal gravity;
63de8de29fSMatthew G. Knepley   struct {
64de8de29fSMatthew G. Knepley     PetscInt Height;
65de8de29fSMatthew G. Knepley     PetscInt Speed;
66de8de29fSMatthew G. Knepley     PetscInt Energy;
67de8de29fSMatthew G. Knepley   } functional;
68de8de29fSMatthew G. Knepley } Physics_SW;
69de8de29fSMatthew G. Knepley 
70de8de29fSMatthew G. Knepley typedef struct {
71de8de29fSMatthew G. Knepley   PetscReal h;
72de8de29fSMatthew G. Knepley   PetscReal uh[DIM];
73de8de29fSMatthew G. Knepley } SWNode;
74de8de29fSMatthew G. Knepley typedef union
75de8de29fSMatthew G. Knepley {
76de8de29fSMatthew G. Knepley   SWNode    swnode;
77de8de29fSMatthew G. Knepley   PetscReal vals[DIM + 1];
78de8de29fSMatthew G. Knepley } SWNodeUnion;
79de8de29fSMatthew G. Knepley 
80de8de29fSMatthew G. Knepley typedef enum {
81de8de29fSMatthew G. Knepley   EULER_IV_SHOCK,
82de8de29fSMatthew G. Knepley   EULER_SS_SHOCK,
83de8de29fSMatthew G. Knepley   EULER_SHOCK_TUBE,
84de8de29fSMatthew G. Knepley   EULER_LINEAR_WAVE
85de8de29fSMatthew G. Knepley } EulerType;
86de8de29fSMatthew G. Knepley 
87de8de29fSMatthew G. Knepley typedef struct {
88de8de29fSMatthew G. Knepley   PetscReal gamma;
89de8de29fSMatthew G. Knepley   PetscReal rhoR;
90de8de29fSMatthew G. Knepley   PetscReal amach;
91de8de29fSMatthew G. Knepley   PetscReal itana;
92de8de29fSMatthew G. Knepley   EulerType type;
93de8de29fSMatthew G. Knepley   struct {
94de8de29fSMatthew G. Knepley     PetscInt Density;
95de8de29fSMatthew G. Knepley     PetscInt Momentum;
96de8de29fSMatthew G. Knepley     PetscInt Energy;
97de8de29fSMatthew G. Knepley     PetscInt Pressure;
98de8de29fSMatthew G. Knepley     PetscInt Speed;
99de8de29fSMatthew G. Knepley   } monitor;
100de8de29fSMatthew G. Knepley } Physics_Euler;
101de8de29fSMatthew G. Knepley 
102de8de29fSMatthew G. Knepley typedef struct {
103de8de29fSMatthew G. Knepley   PetscReal r;
104de8de29fSMatthew G. Knepley   PetscReal ru[DIM];
105de8de29fSMatthew G. Knepley   PetscReal E;
106de8de29fSMatthew G. Knepley } EulerNode;
107de8de29fSMatthew G. Knepley typedef union
108de8de29fSMatthew G. Knepley {
109de8de29fSMatthew G. Knepley   EulerNode eulernode;
110de8de29fSMatthew G. Knepley   PetscReal vals[DIM + 2];
111de8de29fSMatthew G. Knepley } EulerNodeUnion;
112de8de29fSMatthew G. Knepley 
113de8de29fSMatthew G. Knepley static inline PetscReal Dot2Real(const PetscReal *x, const PetscReal *y)
114de8de29fSMatthew G. Knepley {
115de8de29fSMatthew G. Knepley   return x[0] * y[0] + x[1] * y[1];
116de8de29fSMatthew G. Knepley }
117de8de29fSMatthew G. Knepley static inline PetscReal Norm2Real(const PetscReal *x)
118de8de29fSMatthew G. Knepley {
119de8de29fSMatthew G. Knepley   return PetscSqrtReal(PetscAbsReal(Dot2Real(x, x)));
120de8de29fSMatthew G. Knepley }
121de8de29fSMatthew G. Knepley static inline void Normalize2Real(PetscReal *x)
122de8de29fSMatthew G. Knepley {
123de8de29fSMatthew G. Knepley   PetscReal a = 1. / Norm2Real(x);
124de8de29fSMatthew G. Knepley   x[0] *= a;
125de8de29fSMatthew G. Knepley   x[1] *= a;
126de8de29fSMatthew G. Knepley }
127de8de29fSMatthew G. Knepley static inline void Scale2Real(PetscReal a, const PetscReal *x, PetscReal *y)
128de8de29fSMatthew G. Knepley {
129de8de29fSMatthew G. Knepley   y[0] = a * x[0];
130de8de29fSMatthew G. Knepley   y[1] = a * x[1];
131de8de29fSMatthew G. Knepley }
132de8de29fSMatthew G. Knepley 
133de8de29fSMatthew G. Knepley static inline PetscReal DotDIMReal(const PetscReal *x, const PetscReal *y)
134de8de29fSMatthew G. Knepley {
135de8de29fSMatthew G. Knepley   PetscInt  i;
136de8de29fSMatthew G. Knepley   PetscReal prod = 0.0;
137de8de29fSMatthew G. Knepley 
138de8de29fSMatthew G. Knepley   for (i = 0; i < DIM; i++) prod += x[i] * y[i];
139de8de29fSMatthew G. Knepley   return prod;
140de8de29fSMatthew G. Knepley }
141de8de29fSMatthew G. Knepley static inline PetscReal NormDIM(const PetscReal *x)
142de8de29fSMatthew G. Knepley {
143de8de29fSMatthew G. Knepley   return PetscSqrtReal(PetscAbsReal(DotDIMReal(x, x)));
144de8de29fSMatthew G. Knepley }
145de8de29fSMatthew G. Knepley static inline void Waxpy2Real(PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w)
146de8de29fSMatthew G. Knepley {
147de8de29fSMatthew G. Knepley   w[0] = a * x[0] + y[0];
148de8de29fSMatthew G. Knepley   w[1] = a * x[1] + y[1];
149de8de29fSMatthew G. Knepley }
150de8de29fSMatthew G. Knepley 
151de8de29fSMatthew G. Knepley /*
152de8de29fSMatthew G. Knepley  * h_t + div(uh) = 0
153de8de29fSMatthew G. Knepley  * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0
154de8de29fSMatthew G. Knepley  *
155de8de29fSMatthew G. Knepley  * */
156de8de29fSMatthew G. Knepley static PetscErrorCode SWFlux(Physics phys, const PetscReal *n, const SWNode *x, SWNode *f)
157de8de29fSMatthew G. Knepley {
158de8de29fSMatthew G. Knepley   Physics_SW *sw = (Physics_SW *)phys->data;
159de8de29fSMatthew G. Knepley   PetscReal   uhn, u[DIM];
160de8de29fSMatthew G. Knepley   PetscInt    i;
161de8de29fSMatthew G. Knepley 
162de8de29fSMatthew G. Knepley   PetscFunctionBeginUser;
163de8de29fSMatthew G. Knepley   Scale2Real(1. / x->h, x->uh, u);
164de8de29fSMatthew G. Knepley   uhn  = x->uh[0] * n[0] + x->uh[1] * n[1];
165de8de29fSMatthew G. Knepley   f->h = uhn;
166de8de29fSMatthew G. Knepley   for (i = 0; i < DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
167de8de29fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
168de8de29fSMatthew G. Knepley }
169de8de29fSMatthew G. Knepley 
170de8de29fSMatthew G. Knepley 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)
171de8de29fSMatthew G. Knepley {
172de8de29fSMatthew G. Knepley   Physics_SW *sw = (Physics_SW *)phys->data;
173de8de29fSMatthew G. Knepley   PetscReal   cL, cR, speed;
174de8de29fSMatthew G. Knepley   PetscReal   nn[DIM];
175de8de29fSMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
176de8de29fSMatthew G. Knepley   const SWNode *uL = (const SWNode *)xL, *uR = (const SWNode *)xR;
177de8de29fSMatthew G. Knepley #else
178de8de29fSMatthew G. Knepley   SWNodeUnion   uLreal, uRreal;
179de8de29fSMatthew G. Knepley   const SWNode *uL = &uLreal.swnode;
180de8de29fSMatthew G. Knepley   const SWNode *uR = &uRreal.swnode;
181de8de29fSMatthew G. Knepley #endif
182de8de29fSMatthew G. Knepley   SWNodeUnion    fL, fR;
183de8de29fSMatthew G. Knepley   PetscInt       i;
184de8de29fSMatthew G. Knepley   PetscReal      zero = 0.;
185*a7d931c6SMatthew G. Knepley   PetscErrorCode ierr;
186de8de29fSMatthew G. Knepley 
187de8de29fSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX)
188de8de29fSMatthew G. Knepley   uLreal.swnode.h = 0;
189de8de29fSMatthew G. Knepley   uRreal.swnode.h = 0;
190de8de29fSMatthew G. Knepley   for (i = 0; i < 1 + dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
191de8de29fSMatthew G. Knepley   for (i = 0; i < 1 + dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
192de8de29fSMatthew G. Knepley #endif
193*a7d931c6SMatthew G. Knepley 
194de8de29fSMatthew G. Knepley   if (uL->h < 0 || uR->h < 0) {
195*a7d931c6SMatthew G. Knepley     // reconstructed thickness is negative
196*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
197*a7d931c6SMatthew G. Knepley     for (i = 0; i < 1 + dim; ++i) flux[i] = zero / zero;
198*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPop());
199de8de29fSMatthew G. Knepley     return;
200*a7d931c6SMatthew G. Knepley   }
201*a7d931c6SMatthew G. Knepley 
202de8de29fSMatthew G. Knepley   nn[0] = n[0];
203de8de29fSMatthew G. Knepley   nn[1] = n[1];
204de8de29fSMatthew G. Knepley   Normalize2Real(nn);
205*a7d931c6SMatthew G. Knepley   ierr = SWFlux(phys, nn, uL, &fL.swnode);
206*a7d931c6SMatthew G. Knepley   if (ierr) {
207*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
208*a7d931c6SMatthew G. Knepley     for (i = 0; i < 1 + dim; ++i) fL.vals[i] = zero / zero;
209*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPop());
210*a7d931c6SMatthew G. Knepley   }
211*a7d931c6SMatthew G. Knepley   ierr = SWFlux(phys, nn, uR, &fR.swnode);
212*a7d931c6SMatthew G. Knepley   if (ierr) {
213*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
214*a7d931c6SMatthew G. Knepley     for (i = 0; i < 1 + dim; ++i) fR.vals[i] = zero / zero;
215*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPop());
216*a7d931c6SMatthew G. Knepley   }
217de8de29fSMatthew G. Knepley   cL    = PetscSqrtReal(sw->gravity * uL->h);
218de8de29fSMatthew G. Knepley   cR    = PetscSqrtReal(sw->gravity * uR->h); /* gravity wave speed */
219de8de29fSMatthew G. Knepley   speed = PetscMax(PetscAbsReal(Dot2Real(uL->uh, nn) / uL->h) + cL, PetscAbsReal(Dot2Real(uR->uh, nn) / uR->h) + cR);
220de8de29fSMatthew G. Knepley   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);
221de8de29fSMatthew G. Knepley #if 0
222de8de29fSMatthew G. Knepley   PetscPrintf(PETSC_COMM_SELF, "Rusanov Flux (%g)\n", sw->gravity);
223de8de29fSMatthew G. Knepley   for (PetscInt j = 0; j < 3; ++j) {
224de8de29fSMatthew G. Knepley     PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", flux[j]);
225de8de29fSMatthew G. Knepley   }
226de8de29fSMatthew G. Knepley #endif
227de8de29fSMatthew G. Knepley }
228de8de29fSMatthew G. Knepley 
229de8de29fSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
230de8de29fSMatthew G. Knepley CEED_QFUNCTION(PhysicsRiemann_SW_Rusanov_CEED)(void *ctx, CeedInt Q, const CeedScalar *const in[], CeedScalar *const out[])
231de8de29fSMatthew G. Knepley {
232de8de29fSMatthew G. Knepley   const CeedScalar *xL = in[0], *xR = in[1], *geom = in[2];
233de8de29fSMatthew G. Knepley   CeedScalar       *cL = out[0], *cR = out[1];
234de8de29fSMatthew G. Knepley   const Physics_SW *sw = (Physics_SW *)ctx;
235de8de29fSMatthew G. Knepley   struct _n_Physics phys;
236*a7d931c6SMatthew G. Knepley   #if 0
237*a7d931c6SMatthew G. Knepley   const CeedScalar *info = in[3];
238*a7d931c6SMatthew G. Knepley   #endif
239de8de29fSMatthew G. Knepley 
240de8de29fSMatthew G. Knepley   phys.data = (void *)sw;
241de8de29fSMatthew G. Knepley   CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i)
242de8de29fSMatthew G. Knepley   {
243de8de29fSMatthew G. Knepley     const CeedScalar qL[3] = {xL[i + Q * 0], xL[i + Q * 1], xL[i + Q * 2]};
244de8de29fSMatthew G. Knepley     const CeedScalar qR[3] = {xR[i + Q * 0], xR[i + Q * 1], xR[i + Q * 2]};
245de8de29fSMatthew G. Knepley     const CeedScalar n[2]  = {geom[i + Q * 0], geom[i + Q * 1]};
246de8de29fSMatthew G. Knepley     CeedScalar       flux[3];
247de8de29fSMatthew G. Knepley 
248de8de29fSMatthew G. Knepley   #if 0
249*a7d931c6SMatthew G. Knepley     PetscPrintf(PETSC_COMM_SELF, "Face %d Normal\n", (int)info[i + Q * 0]);
250de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < DIM; ++j) {
251de8de29fSMatthew G. Knepley       PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", n[j]);
252de8de29fSMatthew G. Knepley     }
253*a7d931c6SMatthew G. Knepley     PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: left state\n", (int)info[i + Q * 1]);
254de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < DIM + 1; ++j) {
255de8de29fSMatthew G. Knepley       PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", qL[j]);
256de8de29fSMatthew G. Knepley     }
257*a7d931c6SMatthew G. Knepley     PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: right state\n", (int)info[i + Q * 2]);
258de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < DIM + 1; ++j) {
259de8de29fSMatthew G. Knepley       PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", qR[j]);
260de8de29fSMatthew G. Knepley     }
261de8de29fSMatthew G. Knepley   #endif
262de8de29fSMatthew G. Knepley     PhysicsRiemann_SW_Rusanov(DIM, DIM + 1, NULL, n, qL, qR, 0, NULL, flux, &phys);
263de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < 3; ++j) {
264de8de29fSMatthew G. Knepley       cL[i + Q * j] = -flux[j] / geom[i + Q * 2];
265de8de29fSMatthew G. Knepley       cR[i + Q * j] = flux[j] / geom[i + Q * 3];
266de8de29fSMatthew G. Knepley     }
267de8de29fSMatthew G. Knepley   #if 0
268*a7d931c6SMatthew G. Knepley     PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: left flux\n", (int)info[i + Q * 1]);
269de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < DIM + 1; ++j) {
270de8de29fSMatthew G. Knepley       PetscPrintf(PETSC_COMM_SELF, "  | %g | (%g)\n", cL[i + Q * j], geom[i + Q * 2]);
271de8de29fSMatthew G. Knepley     }
272*a7d931c6SMatthew G. Knepley     PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: right flux\n", (int)info[i + Q * 2]);
273de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < DIM + 1; ++j) {
274de8de29fSMatthew G. Knepley       PetscPrintf(PETSC_COMM_SELF, "  | %g | (%g)\n", cR[i + Q * j], geom[i + Q * 3]);
275de8de29fSMatthew G. Knepley     }
276de8de29fSMatthew G. Knepley   #endif
277de8de29fSMatthew G. Knepley   }
278de8de29fSMatthew G. Knepley   return CEED_ERROR_SUCCESS;
279de8de29fSMatthew G. Knepley }
280de8de29fSMatthew G. Knepley #endif
281de8de29fSMatthew G. Knepley 
282de8de29fSMatthew G. Knepley 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)
283de8de29fSMatthew G. Knepley {
284de8de29fSMatthew G. Knepley   Physics_SW *sw = (Physics_SW *)phys->data;
285de8de29fSMatthew G. Knepley   PetscReal   aL, aR;
286de8de29fSMatthew G. Knepley   PetscReal   nn[DIM];
287de8de29fSMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
288de8de29fSMatthew G. Knepley   const SWNode *uL = (const SWNode *)xL, *uR = (const SWNode *)xR;
289de8de29fSMatthew G. Knepley #else
290de8de29fSMatthew G. Knepley   SWNodeUnion   uLreal, uRreal;
291de8de29fSMatthew G. Knepley   const SWNode *uL = &uLreal.swnode;
292de8de29fSMatthew G. Knepley   const SWNode *uR = &uRreal.swnode;
293de8de29fSMatthew G. Knepley #endif
294de8de29fSMatthew G. Knepley   SWNodeUnion    fL, fR;
295de8de29fSMatthew G. Knepley   PetscInt       i;
296de8de29fSMatthew G. Knepley   PetscReal      zero = 0.;
297*a7d931c6SMatthew G. Knepley   PetscErrorCode ierr;
298de8de29fSMatthew G. Knepley 
299de8de29fSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX)
300de8de29fSMatthew G. Knepley   uLreal.swnode.h = 0;
301de8de29fSMatthew G. Knepley   uRreal.swnode.h = 0;
302de8de29fSMatthew G. Knepley   for (i = 0; i < 1 + dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
303de8de29fSMatthew G. Knepley   for (i = 0; i < 1 + dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
304de8de29fSMatthew G. Knepley #endif
305de8de29fSMatthew G. Knepley   if (uL->h <= 0 || uR->h <= 0) {
306*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
307de8de29fSMatthew G. Knepley     for (i = 0; i < 1 + dim; i++) flux[i] = zero;
308*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPop());
309de8de29fSMatthew G. Knepley     return;
310de8de29fSMatthew G. Knepley   } /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */
311de8de29fSMatthew G. Knepley   nn[0] = n[0];
312de8de29fSMatthew G. Knepley   nn[1] = n[1];
313de8de29fSMatthew G. Knepley   Normalize2Real(nn);
314*a7d931c6SMatthew G. Knepley   ierr = SWFlux(phys, nn, uL, &fL.swnode);
315*a7d931c6SMatthew G. Knepley   if (ierr) {
316*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
317*a7d931c6SMatthew G. Knepley     for (i = 0; i < 1 + dim; ++i) fL.vals[i] = zero / zero;
318*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPop());
319*a7d931c6SMatthew G. Knepley   }
320*a7d931c6SMatthew G. Knepley   ierr = SWFlux(phys, nn, uR, &fR.swnode);
321*a7d931c6SMatthew G. Knepley   if (ierr) {
322*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
323*a7d931c6SMatthew G. Knepley     for (i = 0; i < 1 + dim; ++i) fR.vals[i] = zero / zero;
324*a7d931c6SMatthew G. Knepley     PetscCallVoid(PetscFPTrapPop());
325*a7d931c6SMatthew G. Knepley   }
326de8de29fSMatthew G. Knepley   /* gravity wave speed */
327de8de29fSMatthew G. Knepley   aL = PetscSqrtReal(sw->gravity * uL->h);
328de8de29fSMatthew G. Knepley   aR = PetscSqrtReal(sw->gravity * uR->h);
329de8de29fSMatthew G. Knepley   // Defining u_tilda and v_tilda as u and v
330de8de29fSMatthew G. Knepley   PetscReal u_L, u_R;
331de8de29fSMatthew G. Knepley   u_L = Dot2Real(uL->uh, nn) / uL->h;
332de8de29fSMatthew G. Knepley   u_R = Dot2Real(uR->uh, nn) / uR->h;
333de8de29fSMatthew G. Knepley   PetscReal sL, sR;
334de8de29fSMatthew G. Knepley   sL = PetscMin(u_L - aL, u_R - aR);
335de8de29fSMatthew G. Knepley   sR = PetscMax(u_L + aL, u_R + aR);
336de8de29fSMatthew G. Knepley   if (sL > zero) {
337de8de29fSMatthew G. Knepley     for (i = 0; i < dim + 1; i++) flux[i] = fL.vals[i] * Norm2Real(n);
338de8de29fSMatthew G. Knepley   } else if (sR < zero) {
339de8de29fSMatthew G. Knepley     for (i = 0; i < dim + 1; i++) flux[i] = fR.vals[i] * Norm2Real(n);
340de8de29fSMatthew G. Knepley   } else {
341de8de29fSMatthew G. Knepley     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);
342de8de29fSMatthew G. Knepley   }
343de8de29fSMatthew G. Knepley }
344de8de29fSMatthew G. Knepley 
345de8de29fSMatthew G. Knepley static PetscErrorCode Pressure_PG(const PetscReal gamma, const EulerNode *x, PetscReal *p)
346de8de29fSMatthew G. Knepley {
347de8de29fSMatthew G. Knepley   PetscReal ru2;
348de8de29fSMatthew G. Knepley 
349de8de29fSMatthew G. Knepley   PetscFunctionBeginUser;
350de8de29fSMatthew G. Knepley   ru2  = DotDIMReal(x->ru, x->ru);
351de8de29fSMatthew G. Knepley   (*p) = (x->E - 0.5 * ru2 / x->r) * (gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */
352de8de29fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
353de8de29fSMatthew G. Knepley }
354de8de29fSMatthew G. Knepley 
355de8de29fSMatthew G. Knepley static PetscErrorCode SpeedOfSound_PG(const PetscReal gamma, const EulerNode *x, PetscReal *c)
356de8de29fSMatthew G. Knepley {
357de8de29fSMatthew G. Knepley   PetscReal p;
358de8de29fSMatthew G. Knepley 
359de8de29fSMatthew G. Knepley   PetscFunctionBeginUser;
360de8de29fSMatthew G. Knepley   PetscCall(Pressure_PG(gamma, x, &p));
361de8de29fSMatthew G. Knepley   /* gamma = heat capacity ratio */
362de8de29fSMatthew G. Knepley   (*c) = PetscSqrtReal(gamma * p / x->r);
363de8de29fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
364de8de29fSMatthew G. Knepley }
365de8de29fSMatthew G. Knepley 
366de8de29fSMatthew G. Knepley /*
367de8de29fSMatthew G. Knepley  * x = (rho,rho*(u_1),...,rho*e)^T
368de8de29fSMatthew G. Knepley  * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
369de8de29fSMatthew G. Knepley  *
370de8de29fSMatthew G. Knepley  * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
371de8de29fSMatthew G. Knepley  *
372de8de29fSMatthew G. Knepley  */
373de8de29fSMatthew G. Knepley static PetscErrorCode EulerFlux(Physics phys, const PetscReal *n, const EulerNode *x, EulerNode *f)
374de8de29fSMatthew G. Knepley {
375de8de29fSMatthew G. Knepley   Physics_Euler *eu = (Physics_Euler *)phys->data;
376de8de29fSMatthew G. Knepley   PetscReal      nu, p;
377de8de29fSMatthew G. Knepley   PetscInt       i;
378de8de29fSMatthew G. Knepley 
379de8de29fSMatthew G. Knepley   PetscFunctionBeginUser;
380de8de29fSMatthew G. Knepley   PetscCall(Pressure_PG(eu->gamma, x, &p));
381de8de29fSMatthew G. Knepley   nu   = DotDIMReal(x->ru, n);
382de8de29fSMatthew G. Knepley   f->r = nu;                                                     /* A rho u */
383de8de29fSMatthew G. Knepley   nu /= x->r;                                                    /* A u */
384de8de29fSMatthew G. Knepley   for (i = 0; i < DIM; i++) f->ru[i] = nu * x->ru[i] + n[i] * p; /* r u^2 + p */
385de8de29fSMatthew G. Knepley   f->E = nu * (x->E + p);                                        /* u(e+p) */
386de8de29fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
387de8de29fSMatthew G. Knepley }
388de8de29fSMatthew G. Knepley 
389de8de29fSMatthew G. Knepley /* Godunov fluxs */
390de8de29fSMatthew G. Knepley static PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
391de8de29fSMatthew G. Knepley {
392de8de29fSMatthew G. Knepley   /* System generated locals */
393de8de29fSMatthew G. Knepley   PetscScalar ret_val;
394de8de29fSMatthew G. Knepley 
395de8de29fSMatthew G. Knepley   if (PetscRealPart(*test) > 0.) goto L10;
396de8de29fSMatthew G. Knepley   ret_val = *b;
397de8de29fSMatthew G. Knepley   return ret_val;
398de8de29fSMatthew G. Knepley L10:
399de8de29fSMatthew G. Knepley   ret_val = *a;
400de8de29fSMatthew G. Knepley   return ret_val;
401de8de29fSMatthew G. Knepley } /* cvmgp_ */
402de8de29fSMatthew G. Knepley 
403de8de29fSMatthew G. Knepley static PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
404de8de29fSMatthew G. Knepley {
405de8de29fSMatthew G. Knepley   /* System generated locals */
406de8de29fSMatthew G. Knepley   PetscScalar ret_val;
407de8de29fSMatthew G. Knepley 
408de8de29fSMatthew G. Knepley   if (PetscRealPart(*test) < 0.) goto L10;
409de8de29fSMatthew G. Knepley   ret_val = *b;
410de8de29fSMatthew G. Knepley   return ret_val;
411de8de29fSMatthew G. Knepley L10:
412de8de29fSMatthew G. Knepley   ret_val = *a;
413de8de29fSMatthew G. Knepley   return ret_val;
414de8de29fSMatthew G. Knepley } /* cvmgm_ */
415de8de29fSMatthew G. Knepley 
416de8de29fSMatthew G. Knepley static 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)
417de8de29fSMatthew G. Knepley {
418de8de29fSMatthew G. Knepley   /* Initialized data */
419de8de29fSMatthew G. Knepley 
420de8de29fSMatthew G. Knepley   static PetscScalar smallp = 1e-8;
421de8de29fSMatthew G. Knepley 
422de8de29fSMatthew G. Knepley   /* System generated locals */
423de8de29fSMatthew G. Knepley   int         i__1;
424de8de29fSMatthew G. Knepley   PetscScalar d__1, d__2;
425de8de29fSMatthew G. Knepley 
426de8de29fSMatthew G. Knepley   /* Local variables */
427de8de29fSMatthew G. Knepley   static int         i0;
428de8de29fSMatthew G. Knepley   static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2;
429de8de29fSMatthew G. Knepley   static int         iwave;
430de8de29fSMatthew G. Knepley   static PetscScalar gascl4, gascr4, cstarl, dpstar, cstarr;
431de8de29fSMatthew G. Knepley   /* static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascr1, gascr2, gascr3; */
432de8de29fSMatthew G. Knepley   static int         iterno;
433de8de29fSMatthew G. Knepley   static PetscScalar ustarl, ustarr, rarepr1, rarepr2;
434de8de29fSMatthew G. Knepley 
435de8de29fSMatthew G. Knepley   /* gascl1 = *gaml - 1.; */
436de8de29fSMatthew G. Knepley   /* gascl2 = (*gaml + 1.) * .5; */
437de8de29fSMatthew G. Knepley   /* gascl3 = gascl2 / *gaml; */
438de8de29fSMatthew G. Knepley   gascl4 = 1. / (*gaml - 1.);
439de8de29fSMatthew G. Knepley 
440de8de29fSMatthew G. Knepley   /* gascr1 = *gamr - 1.; */
441de8de29fSMatthew G. Knepley   /* gascr2 = (*gamr + 1.) * .5; */
442de8de29fSMatthew G. Knepley   /* gascr3 = gascr2 / *gamr; */
443de8de29fSMatthew G. Knepley   gascr4 = 1. / (*gamr - 1.);
444de8de29fSMatthew G. Knepley   iterno = 10;
445de8de29fSMatthew G. Knepley   /*        find pstar: */
446de8de29fSMatthew G. Knepley   cl = PetscSqrtScalar(*gaml * *pl / *rl);
447de8de29fSMatthew G. Knepley   cr = PetscSqrtScalar(*gamr * *pr / *rr);
448de8de29fSMatthew G. Knepley   wl = *rl * cl;
449de8de29fSMatthew G. Knepley   wr = *rr * cr;
450de8de29fSMatthew G. Knepley   /* csqrl = wl * wl; */
451de8de29fSMatthew G. Knepley   /* csqrr = wr * wr; */
452de8de29fSMatthew G. Knepley   *pstar  = (wl * *pr + wr * *pl) / (wl + wr);
453de8de29fSMatthew G. Knepley   *pstar  = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
454de8de29fSMatthew G. Knepley   pst     = *pl / *pr;
455de8de29fSMatthew G. Knepley   skpr1   = cr * (pst - 1.) * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
456de8de29fSMatthew G. Knepley   d__1    = (*gamr - 1.) / (*gamr * 2.);
457de8de29fSMatthew G. Knepley   rarepr2 = gascr4 * 2. * cr * (1. - PetscPowScalar(pst, d__1));
458de8de29fSMatthew G. Knepley   pst     = *pr / *pl;
459de8de29fSMatthew G. Knepley   skpr2   = cl * (pst - 1.) * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
460de8de29fSMatthew G. Knepley   d__1    = (*gaml - 1.) / (*gaml * 2.);
461de8de29fSMatthew G. Knepley   rarepr1 = gascl4 * 2. * cl * (1. - PetscPowScalar(pst, d__1));
462de8de29fSMatthew G. Knepley   durl    = *uxr - *uxl;
463de8de29fSMatthew G. Knepley   if (PetscRealPart(*pr) < PetscRealPart(*pl)) {
464de8de29fSMatthew G. Knepley     if (PetscRealPart(durl) >= PetscRealPart(rarepr1)) {
465de8de29fSMatthew G. Knepley       iwave = 100;
466de8de29fSMatthew G. Knepley     } else if (PetscRealPart(durl) <= PetscRealPart(-skpr1)) {
467de8de29fSMatthew G. Knepley       iwave = 300;
468de8de29fSMatthew G. Knepley     } else {
469de8de29fSMatthew G. Knepley       iwave = 400;
470de8de29fSMatthew G. Knepley     }
471de8de29fSMatthew G. Knepley   } else {
472de8de29fSMatthew G. Knepley     if (PetscRealPart(durl) >= PetscRealPart(rarepr2)) {
473de8de29fSMatthew G. Knepley       iwave = 100;
474de8de29fSMatthew G. Knepley     } else if (PetscRealPart(durl) <= PetscRealPart(-skpr2)) {
475de8de29fSMatthew G. Knepley       iwave = 300;
476de8de29fSMatthew G. Knepley     } else {
477de8de29fSMatthew G. Knepley       iwave = 200;
478de8de29fSMatthew G. Knepley     }
479de8de29fSMatthew G. Knepley   }
480de8de29fSMatthew G. Knepley   if (iwave == 100) {
481de8de29fSMatthew G. Knepley     /*     1-wave: rarefaction wave, 3-wave: rarefaction wave */
482de8de29fSMatthew G. Knepley     /*     case (100) */
483de8de29fSMatthew G. Knepley     i__1 = iterno;
484de8de29fSMatthew G. Knepley     for (i0 = 1; i0 <= i__1; ++i0) {
485de8de29fSMatthew G. Knepley       d__1    = *pstar / *pl;
486de8de29fSMatthew G. Knepley       d__2    = 1. / *gaml;
487de8de29fSMatthew G. Knepley       *rstarl = *rl * PetscPowScalar(d__1, d__2);
488de8de29fSMatthew G. Knepley       cstarl  = PetscSqrtScalar(*gaml * *pstar / *rstarl);
489de8de29fSMatthew G. Knepley       ustarl  = *uxl - gascl4 * 2. * (cstarl - cl);
490de8de29fSMatthew G. Knepley       zl      = *rstarl * cstarl;
491de8de29fSMatthew G. Knepley       d__1    = *pstar / *pr;
492de8de29fSMatthew G. Knepley       d__2    = 1. / *gamr;
493de8de29fSMatthew G. Knepley       *rstarr = *rr * PetscPowScalar(d__1, d__2);
494de8de29fSMatthew G. Knepley       cstarr  = PetscSqrtScalar(*gamr * *pstar / *rstarr);
495de8de29fSMatthew G. Knepley       ustarr  = *uxr + gascr4 * 2. * (cstarr - cr);
496de8de29fSMatthew G. Knepley       zr      = *rstarr * cstarr;
497de8de29fSMatthew G. Knepley       dpstar  = zl * zr * (ustarr - ustarl) / (zl + zr);
498de8de29fSMatthew G. Knepley       *pstar -= dpstar;
499de8de29fSMatthew G. Knepley       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
500de8de29fSMatthew G. Knepley       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
501de8de29fSMatthew G. Knepley #if 0
502de8de29fSMatthew G. Knepley         break;
503de8de29fSMatthew G. Knepley #endif
504de8de29fSMatthew G. Knepley       }
505de8de29fSMatthew G. Knepley     }
506de8de29fSMatthew G. Knepley     /*     1-wave: shock wave, 3-wave: rarefaction wave */
507de8de29fSMatthew G. Knepley   } else if (iwave == 200) {
508de8de29fSMatthew G. Knepley     /*     case (200) */
509de8de29fSMatthew G. Knepley     i__1 = iterno;
510de8de29fSMatthew G. Knepley     for (i0 = 1; i0 <= i__1; ++i0) {
511de8de29fSMatthew G. Knepley       pst     = *pstar / *pl;
512de8de29fSMatthew G. Knepley       ustarl  = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
513de8de29fSMatthew G. Knepley       zl      = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
514de8de29fSMatthew G. Knepley       d__1    = *pstar / *pr;
515de8de29fSMatthew G. Knepley       d__2    = 1. / *gamr;
516de8de29fSMatthew G. Knepley       *rstarr = *rr * PetscPowScalar(d__1, d__2);
517de8de29fSMatthew G. Knepley       cstarr  = PetscSqrtScalar(*gamr * *pstar / *rstarr);
518de8de29fSMatthew G. Knepley       zr      = *rstarr * cstarr;
519de8de29fSMatthew G. Knepley       ustarr  = *uxr + gascr4 * 2. * (cstarr - cr);
520de8de29fSMatthew G. Knepley       dpstar  = zl * zr * (ustarr - ustarl) / (zl + zr);
521de8de29fSMatthew G. Knepley       *pstar -= dpstar;
522de8de29fSMatthew G. Knepley       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
523de8de29fSMatthew G. Knepley       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
524de8de29fSMatthew G. Knepley #if 0
525de8de29fSMatthew G. Knepley         break;
526de8de29fSMatthew G. Knepley #endif
527de8de29fSMatthew G. Knepley       }
528de8de29fSMatthew G. Knepley     }
529de8de29fSMatthew G. Knepley     /*     1-wave: shock wave, 3-wave: shock */
530de8de29fSMatthew G. Knepley   } else if (iwave == 300) {
531de8de29fSMatthew G. Knepley     /*     case (300) */
532de8de29fSMatthew G. Knepley     i__1 = iterno;
533de8de29fSMatthew G. Knepley     for (i0 = 1; i0 <= i__1; ++i0) {
534de8de29fSMatthew G. Knepley       pst    = *pstar / *pl;
535de8de29fSMatthew G. Knepley       ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
536de8de29fSMatthew G. Knepley       zl     = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
537de8de29fSMatthew G. Knepley       pst    = *pstar / *pr;
538de8de29fSMatthew G. Knepley       ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
539de8de29fSMatthew G. Knepley       zr     = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
540de8de29fSMatthew G. Knepley       dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
541de8de29fSMatthew G. Knepley       *pstar -= dpstar;
542de8de29fSMatthew G. Knepley       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
543de8de29fSMatthew G. Knepley       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
544de8de29fSMatthew G. Knepley #if 0
545de8de29fSMatthew G. Knepley         break;
546de8de29fSMatthew G. Knepley #endif
547de8de29fSMatthew G. Knepley       }
548de8de29fSMatthew G. Knepley     }
549de8de29fSMatthew G. Knepley     /*     1-wave: rarefaction wave, 3-wave: shock */
550de8de29fSMatthew G. Knepley   } else if (iwave == 400) {
551de8de29fSMatthew G. Knepley     /*     case (400) */
552de8de29fSMatthew G. Knepley     i__1 = iterno;
553de8de29fSMatthew G. Knepley     for (i0 = 1; i0 <= i__1; ++i0) {
554de8de29fSMatthew G. Knepley       d__1    = *pstar / *pl;
555de8de29fSMatthew G. Knepley       d__2    = 1. / *gaml;
556de8de29fSMatthew G. Knepley       *rstarl = *rl * PetscPowScalar(d__1, d__2);
557de8de29fSMatthew G. Knepley       cstarl  = PetscSqrtScalar(*gaml * *pstar / *rstarl);
558de8de29fSMatthew G. Knepley       ustarl  = *uxl - gascl4 * 2. * (cstarl - cl);
559de8de29fSMatthew G. Knepley       zl      = *rstarl * cstarl;
560de8de29fSMatthew G. Knepley       pst     = *pstar / *pr;
561de8de29fSMatthew G. Knepley       ustarr  = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
562de8de29fSMatthew G. Knepley       zr      = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
563de8de29fSMatthew G. Knepley       dpstar  = zl * zr * (ustarr - ustarl) / (zl + zr);
564de8de29fSMatthew G. Knepley       *pstar -= dpstar;
565de8de29fSMatthew G. Knepley       *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp));
566de8de29fSMatthew G. Knepley       if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
567de8de29fSMatthew G. Knepley #if 0
568de8de29fSMatthew G. Knepley               break;
569de8de29fSMatthew G. Knepley #endif
570de8de29fSMatthew G. Knepley       }
571de8de29fSMatthew G. Knepley     }
572de8de29fSMatthew G. Knepley   }
573de8de29fSMatthew G. Knepley 
574de8de29fSMatthew G. Knepley   *ustar = (zl * ustarr + zr * ustarl) / (zl + zr);
575de8de29fSMatthew G. Knepley   if (PetscRealPart(*pstar) > PetscRealPart(*pl)) {
576de8de29fSMatthew G. Knepley     pst     = *pstar / *pl;
577de8de29fSMatthew G. Knepley     *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + *gaml + 1.) * *rl;
578de8de29fSMatthew G. Knepley   }
579de8de29fSMatthew G. Knepley   if (PetscRealPart(*pstar) > PetscRealPart(*pr)) {
580de8de29fSMatthew G. Knepley     pst     = *pstar / *pr;
581de8de29fSMatthew G. Knepley     *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + *gamr + 1.) * *rr;
582de8de29fSMatthew G. Knepley   }
583de8de29fSMatthew G. Knepley   return iwave;
584de8de29fSMatthew G. Knepley }
585de8de29fSMatthew G. Knepley 
586de8de29fSMatthew G. Knepley static PetscScalar sign(PetscScalar x)
587de8de29fSMatthew G. Knepley {
588de8de29fSMatthew G. Knepley   if (PetscRealPart(x) > 0) return 1.0;
589de8de29fSMatthew G. Knepley   if (PetscRealPart(x) < 0) return -1.0;
590de8de29fSMatthew G. Knepley   return 0.0;
591de8de29fSMatthew G. Knepley }
592de8de29fSMatthew G. Knepley /*        Riemann Solver */
593de8de29fSMatthew G. Knepley /* -------------------------------------------------------------------- */
594de8de29fSMatthew G. Knepley static 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)
595de8de29fSMatthew G. Knepley {
596de8de29fSMatthew G. Knepley   /* System generated locals */
597de8de29fSMatthew G. Knepley   PetscScalar d__1, d__2;
598de8de29fSMatthew G. Knepley 
599de8de29fSMatthew G. Knepley   /* Local variables */
600de8de29fSMatthew G. Knepley   static PetscScalar s, c0, p0, r0, u0, w0, x0, x2, ri, cx, sgn0, wsp0, gasc1, gasc2, gasc3, gasc4;
601de8de29fSMatthew G. Knepley   static PetscScalar cstar, pstar, rstar, ustar, xstar, wspst, ushock, streng, rstarl, rstarr, rstars;
602de8de29fSMatthew G. Knepley   int                iwave;
603de8de29fSMatthew G. Knepley 
604de8de29fSMatthew G. Knepley   if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) {
605de8de29fSMatthew G. Knepley     *rx  = *rl;
606de8de29fSMatthew G. Knepley     *px  = *pl;
607de8de29fSMatthew G. Knepley     *uxm = *uxl;
608de8de29fSMatthew G. Knepley     *gam = *gaml;
609de8de29fSMatthew G. Knepley     x2   = *xcen + *uxm * *dtt;
610de8de29fSMatthew G. Knepley 
611de8de29fSMatthew G. Knepley     if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
612de8de29fSMatthew G. Knepley       *utx  = *utr;
613de8de29fSMatthew G. Knepley       *ubx  = *ubr;
614de8de29fSMatthew G. Knepley       *rho1 = *rho1r;
615de8de29fSMatthew G. Knepley     } else {
616de8de29fSMatthew G. Knepley       *utx  = *utl;
617de8de29fSMatthew G. Knepley       *ubx  = *ubl;
618de8de29fSMatthew G. Knepley       *rho1 = *rho1l;
619de8de29fSMatthew G. Knepley     }
620de8de29fSMatthew G. Knepley     return 0;
621de8de29fSMatthew G. Knepley   }
622de8de29fSMatthew G. Knepley   iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar);
623de8de29fSMatthew G. Knepley 
624de8de29fSMatthew G. Knepley   x2   = *xcen + ustar * *dtt;
625de8de29fSMatthew G. Knepley   d__1 = *xp - x2;
626de8de29fSMatthew G. Knepley   sgn0 = sign(d__1);
627de8de29fSMatthew G. Knepley   /*            x is in 3-wave if sgn0 = 1 */
628de8de29fSMatthew G. Knepley   /*            x is in 1-wave if sgn0 = -1 */
629de8de29fSMatthew G. Knepley   r0     = cvmgm_(rl, rr, &sgn0);
630de8de29fSMatthew G. Knepley   p0     = cvmgm_(pl, pr, &sgn0);
631de8de29fSMatthew G. Knepley   u0     = cvmgm_(uxl, uxr, &sgn0);
632de8de29fSMatthew G. Knepley   *gam   = cvmgm_(gaml, gamr, &sgn0);
633de8de29fSMatthew G. Knepley   gasc1  = *gam - 1.;
634de8de29fSMatthew G. Knepley   gasc2  = (*gam + 1.) * .5;
635de8de29fSMatthew G. Knepley   gasc3  = gasc2 / *gam;
636de8de29fSMatthew G. Knepley   gasc4  = 1. / (*gam - 1.);
637de8de29fSMatthew G. Knepley   c0     = PetscSqrtScalar(*gam * p0 / r0);
638de8de29fSMatthew G. Knepley   streng = pstar - p0;
639de8de29fSMatthew G. Knepley   w0     = *gam * r0 * p0 * (gasc3 * streng / p0 + 1.);
640de8de29fSMatthew G. Knepley   rstars = r0 / (1. - r0 * streng / w0);
641de8de29fSMatthew G. Knepley   d__1   = p0 / pstar;
642de8de29fSMatthew G. Knepley   d__2   = -1. / *gam;
643de8de29fSMatthew G. Knepley   rstarr = r0 * PetscPowScalar(d__1, d__2);
644de8de29fSMatthew G. Knepley   rstar  = cvmgm_(&rstarr, &rstars, &streng);
645de8de29fSMatthew G. Knepley   w0     = PetscSqrtScalar(w0);
646de8de29fSMatthew G. Knepley   cstar  = PetscSqrtScalar(*gam * pstar / rstar);
647de8de29fSMatthew G. Knepley   wsp0   = u0 + sgn0 * c0;
648de8de29fSMatthew G. Knepley   wspst  = ustar + sgn0 * cstar;
649de8de29fSMatthew G. Knepley   ushock = ustar + sgn0 * w0 / rstar;
650de8de29fSMatthew G. Knepley   wspst  = cvmgp_(&ushock, &wspst, &streng);
651de8de29fSMatthew G. Knepley   wsp0   = cvmgp_(&ushock, &wsp0, &streng);
652de8de29fSMatthew G. Knepley   x0     = *xcen + wsp0 * *dtt;
653de8de29fSMatthew G. Knepley   xstar  = *xcen + wspst * *dtt;
654de8de29fSMatthew G. Knepley   /*           using gas formula to evaluate rarefaction wave */
655de8de29fSMatthew G. Knepley   /*            ri : reiman invariant */
656de8de29fSMatthew G. Knepley   ri   = u0 - sgn0 * 2. * gasc4 * c0;
657de8de29fSMatthew G. Knepley   cx   = sgn0 * .5 * gasc1 / gasc2 * ((*xp - *xcen) / *dtt - ri);
658de8de29fSMatthew G. Knepley   *uxm = ri + sgn0 * 2. * gasc4 * cx;
659de8de29fSMatthew G. Knepley   s    = p0 / PetscPowScalar(r0, *gam);
660de8de29fSMatthew G. Knepley   d__1 = cx * cx / (*gam * s);
661de8de29fSMatthew G. Knepley   *rx  = PetscPowScalar(d__1, gasc4);
662de8de29fSMatthew G. Knepley   *px  = cx * cx * *rx / *gam;
663de8de29fSMatthew G. Knepley   d__1 = sgn0 * (x0 - *xp);
664de8de29fSMatthew G. Knepley   *rx  = cvmgp_(rx, &r0, &d__1);
665de8de29fSMatthew G. Knepley   d__1 = sgn0 * (x0 - *xp);
666de8de29fSMatthew G. Knepley   *px  = cvmgp_(px, &p0, &d__1);
667de8de29fSMatthew G. Knepley   d__1 = sgn0 * (x0 - *xp);
668de8de29fSMatthew G. Knepley   *uxm = cvmgp_(uxm, &u0, &d__1);
669de8de29fSMatthew G. Knepley   d__1 = sgn0 * (xstar - *xp);
670de8de29fSMatthew G. Knepley   *rx  = cvmgm_(rx, &rstar, &d__1);
671de8de29fSMatthew G. Knepley   d__1 = sgn0 * (xstar - *xp);
672de8de29fSMatthew G. Knepley   *px  = cvmgm_(px, &pstar, &d__1);
673de8de29fSMatthew G. Knepley   d__1 = sgn0 * (xstar - *xp);
674de8de29fSMatthew G. Knepley   *uxm = cvmgm_(uxm, &ustar, &d__1);
675de8de29fSMatthew G. Knepley   if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
676de8de29fSMatthew G. Knepley     *utx  = *utr;
677de8de29fSMatthew G. Knepley     *ubx  = *ubr;
678de8de29fSMatthew G. Knepley     *rho1 = *rho1r;
679de8de29fSMatthew G. Knepley   } else {
680de8de29fSMatthew G. Knepley     *utx  = *utl;
681de8de29fSMatthew G. Knepley     *ubx  = *ubl;
682de8de29fSMatthew G. Knepley     *rho1 = *rho1l;
683de8de29fSMatthew G. Knepley   }
684de8de29fSMatthew G. Knepley   return iwave;
685de8de29fSMatthew G. Knepley }
686de8de29fSMatthew G. Knepley 
687de8de29fSMatthew G. Knepley static int godunovflux(const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, int ndim, PetscReal gamma)
688de8de29fSMatthew G. Knepley {
689de8de29fSMatthew G. Knepley   /* System generated locals */
690de8de29fSMatthew G. Knepley   int         i__1, iwave;
691de8de29fSMatthew G. Knepley   PetscScalar d__1, d__2, d__3;
692de8de29fSMatthew G. Knepley 
693de8de29fSMatthew G. Knepley   /* Local variables */
694de8de29fSMatthew G. Knepley   static int         k;
695de8de29fSMatthew G. Knepley   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;
696de8de29fSMatthew G. Knepley 
697de8de29fSMatthew G. Knepley   /* Function Body */
698de8de29fSMatthew G. Knepley   xcen = 0.;
699de8de29fSMatthew G. Knepley   xp   = 0.;
700de8de29fSMatthew G. Knepley   i__1 = ndim;
701de8de29fSMatthew G. Knepley   for (k = 1; k <= i__1; ++k) {
702de8de29fSMatthew G. Knepley     tg[k - 1] = 0.;
703de8de29fSMatthew G. Knepley     bn[k - 1] = 0.;
704de8de29fSMatthew G. Knepley   }
705de8de29fSMatthew G. Knepley   dtt = 1.;
706de8de29fSMatthew G. Knepley   if (ndim == 3) {
707de8de29fSMatthew G. Knepley     if (nn[0] == 0. && nn[1] == 0.) {
708de8de29fSMatthew G. Knepley       tg[0] = 1.;
709de8de29fSMatthew G. Knepley     } else {
710de8de29fSMatthew G. Knepley       tg[0] = -nn[1];
711de8de29fSMatthew G. Knepley       tg[1] = nn[0];
712de8de29fSMatthew G. Knepley     }
713de8de29fSMatthew G. Knepley     /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
714de8de29fSMatthew G. Knepley     /*           tg=tg/tmp */
715de8de29fSMatthew G. Knepley     bn[0] = -nn[2] * tg[1];
716de8de29fSMatthew G. Knepley     bn[1] = nn[2] * tg[0];
717de8de29fSMatthew G. Knepley     bn[2] = nn[0] * tg[1] - nn[1] * tg[0];
718de8de29fSMatthew G. Knepley     /* Computing 2nd power */
719de8de29fSMatthew G. Knepley     d__1 = bn[0];
720de8de29fSMatthew G. Knepley     /* Computing 2nd power */
721de8de29fSMatthew G. Knepley     d__2 = bn[1];
722de8de29fSMatthew G. Knepley     /* Computing 2nd power */
723de8de29fSMatthew G. Knepley     d__3 = bn[2];
724de8de29fSMatthew G. Knepley     tmp  = PetscSqrtScalar(d__1 * d__1 + d__2 * d__2 + d__3 * d__3);
725de8de29fSMatthew G. Knepley     i__1 = ndim;
726de8de29fSMatthew G. Knepley     for (k = 1; k <= i__1; ++k) bn[k - 1] /= tmp;
727de8de29fSMatthew G. Knepley   } else if (ndim == 2) {
728de8de29fSMatthew G. Knepley     tg[0] = -nn[1];
729de8de29fSMatthew G. Knepley     tg[1] = nn[0];
730de8de29fSMatthew G. Knepley     /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
731de8de29fSMatthew G. Knepley     /*           tg=tg/tmp */
732de8de29fSMatthew G. Knepley     bn[0] = 0.;
733de8de29fSMatthew G. Knepley     bn[1] = 0.;
734de8de29fSMatthew G. Knepley     bn[2] = 1.;
735de8de29fSMatthew G. Knepley   }
736de8de29fSMatthew G. Knepley   rl   = ul[0];
737de8de29fSMatthew G. Knepley   rr   = ur[0];
738de8de29fSMatthew G. Knepley   uxl  = 0.;
739de8de29fSMatthew G. Knepley   uxr  = 0.;
740de8de29fSMatthew G. Knepley   utl  = 0.;
741de8de29fSMatthew G. Knepley   utr  = 0.;
742de8de29fSMatthew G. Knepley   ubl  = 0.;
743de8de29fSMatthew G. Knepley   ubr  = 0.;
744de8de29fSMatthew G. Knepley   i__1 = ndim;
745de8de29fSMatthew G. Knepley   for (k = 1; k <= i__1; ++k) {
746de8de29fSMatthew G. Knepley     uxl += ul[k] * nn[k - 1];
747de8de29fSMatthew G. Knepley     uxr += ur[k] * nn[k - 1];
748de8de29fSMatthew G. Knepley     utl += ul[k] * tg[k - 1];
749de8de29fSMatthew G. Knepley     utr += ur[k] * tg[k - 1];
750de8de29fSMatthew G. Knepley     ubl += ul[k] * bn[k - 1];
751de8de29fSMatthew G. Knepley     ubr += ur[k] * bn[k - 1];
752de8de29fSMatthew G. Knepley   }
753de8de29fSMatthew G. Knepley   uxl /= rl;
754de8de29fSMatthew G. Knepley   uxr /= rr;
755de8de29fSMatthew G. Knepley   utl /= rl;
756de8de29fSMatthew G. Knepley   utr /= rr;
757de8de29fSMatthew G. Knepley   ubl /= rl;
758de8de29fSMatthew G. Knepley   ubr /= rr;
759de8de29fSMatthew G. Knepley 
760de8de29fSMatthew G. Knepley   gaml = gamma;
761de8de29fSMatthew G. Knepley   gamr = gamma;
762de8de29fSMatthew G. Knepley   /* Computing 2nd power */
763de8de29fSMatthew G. Knepley   d__1 = uxl;
764de8de29fSMatthew G. Knepley   /* Computing 2nd power */
765de8de29fSMatthew G. Knepley   d__2 = utl;
766de8de29fSMatthew G. Knepley   /* Computing 2nd power */
767de8de29fSMatthew G. Knepley   d__3 = ubl;
768de8de29fSMatthew G. Knepley   pl   = (gamma - 1.) * (ul[ndim + 1] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
769de8de29fSMatthew G. Knepley   /* Computing 2nd power */
770de8de29fSMatthew G. Knepley   d__1 = uxr;
771de8de29fSMatthew G. Knepley   /* Computing 2nd power */
772de8de29fSMatthew G. Knepley   d__2 = utr;
773de8de29fSMatthew G. Knepley   /* Computing 2nd power */
774de8de29fSMatthew G. Knepley   d__3  = ubr;
775de8de29fSMatthew G. Knepley   pr    = (gamma - 1.) * (ur[ndim + 1] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
776de8de29fSMatthew G. Knepley   rho1l = rl;
777de8de29fSMatthew G. Knepley   rho1r = rr;
778de8de29fSMatthew G. Knepley 
779de8de29fSMatthew G. Knepley   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);
780de8de29fSMatthew G. Knepley 
781de8de29fSMatthew G. Knepley   flux[0] = rhom * unm;
782de8de29fSMatthew G. Knepley   fn      = rhom * unm * unm + pm;
783de8de29fSMatthew G. Knepley   ft      = rhom * unm * utm;
784de8de29fSMatthew G. Knepley   /*           flux(2)=fn*nn(1)+ft*nn(2) */
785de8de29fSMatthew G. Knepley   /*           flux(3)=fn*tg(1)+ft*tg(2) */
786de8de29fSMatthew G. Knepley   flux[1] = fn * nn[0] + ft * tg[0];
787de8de29fSMatthew G. Knepley   flux[2] = fn * nn[1] + ft * tg[1];
788de8de29fSMatthew G. Knepley   /*           flux(2)=rhom*unm*(unm)+pm */
789de8de29fSMatthew G. Knepley   /*           flux(3)=rhom*(unm)*utm */
790de8de29fSMatthew G. Knepley   if (ndim == 3) flux[3] = rhom * unm * ubm;
791de8de29fSMatthew G. Knepley   flux[ndim + 1] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm;
792de8de29fSMatthew G. Knepley   return iwave;
793de8de29fSMatthew G. Knepley } /* godunovflux_ */
794de8de29fSMatthew G. Knepley 
795de8de29fSMatthew G. Knepley /* PetscReal* => EulerNode* conversion */
796de8de29fSMatthew G. Knepley 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)
797de8de29fSMatthew G. Knepley {
798de8de29fSMatthew G. Knepley   Physics_Euler  *eu    = (Physics_Euler *)phys->data;
799de8de29fSMatthew G. Knepley   const PetscReal gamma = eu->gamma;
800*a7d931c6SMatthew G. Knepley   PetscReal       zero  = 0.;
801de8de29fSMatthew G. Knepley   PetscReal       cL, cR, speed, velL, velR, nn[DIM], s2;
802de8de29fSMatthew G. Knepley   PetscInt        i;
803de8de29fSMatthew G. Knepley   PetscErrorCode  ierr;
804de8de29fSMatthew G. Knepley 
805de8de29fSMatthew G. Knepley   PetscFunctionBeginUser;
806de8de29fSMatthew G. Knepley   for (i = 0, s2 = 0.; i < DIM; i++) {
807de8de29fSMatthew G. Knepley     nn[i] = n[i];
808de8de29fSMatthew G. Knepley     s2 += nn[i] * nn[i];
809de8de29fSMatthew G. Knepley   }
810de8de29fSMatthew G. Knepley   s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */
811de8de29fSMatthew G. Knepley   for (i = 0.; i < DIM; i++) nn[i] /= s2;
812de8de29fSMatthew G. Knepley   if (0) { /* Rusanov */
813de8de29fSMatthew G. Knepley     const EulerNode *uL = (const EulerNode *)xL, *uR = (const EulerNode *)xR;
814de8de29fSMatthew G. Knepley     EulerNodeUnion   fL, fR;
815*a7d931c6SMatthew G. Knepley     ierr = EulerFlux(phys, nn, uL, &fL.eulernode);
816*a7d931c6SMatthew G. Knepley     if (ierr) {
817*a7d931c6SMatthew G. Knepley       PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
818*a7d931c6SMatthew G. Knepley       for (i = 0; i < 2 + dim; i++) fL.vals[i] = zero / zero;
819*a7d931c6SMatthew G. Knepley       PetscCallVoid(PetscFPTrapPop());
820*a7d931c6SMatthew G. Knepley     }
821*a7d931c6SMatthew G. Knepley     ierr = EulerFlux(phys, nn, uR, &fR.eulernode);
822*a7d931c6SMatthew G. Knepley     if (ierr) {
823*a7d931c6SMatthew G. Knepley       PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
824*a7d931c6SMatthew G. Knepley       for (i = 0; i < 2 + dim; i++) fR.vals[i] = zero / zero;
825*a7d931c6SMatthew G. Knepley       PetscCallVoid(PetscFPTrapPop());
826*a7d931c6SMatthew G. Knepley     }
827de8de29fSMatthew G. Knepley     ierr = SpeedOfSound_PG(gamma, uL, &cL);
828*a7d931c6SMatthew G. Knepley     if (ierr) {
829*a7d931c6SMatthew G. Knepley       PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
830*a7d931c6SMatthew G. Knepley       cL = zero / zero;
831*a7d931c6SMatthew G. Knepley       PetscCallVoid(PetscFPTrapPop());
832*a7d931c6SMatthew G. Knepley     }
833de8de29fSMatthew G. Knepley     ierr = SpeedOfSound_PG(gamma, uR, &cR);
834*a7d931c6SMatthew G. Knepley     if (ierr) {
835*a7d931c6SMatthew G. Knepley       PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
836*a7d931c6SMatthew G. Knepley       cR = zero / zero;
837*a7d931c6SMatthew G. Knepley       PetscCallVoid(PetscFPTrapPop());
838*a7d931c6SMatthew G. Knepley     }
839de8de29fSMatthew G. Knepley     velL  = DotDIMReal(uL->ru, nn) / uL->r;
840de8de29fSMatthew G. Knepley     velR  = DotDIMReal(uR->ru, nn) / uR->r;
841de8de29fSMatthew G. Knepley     speed = PetscMax(velR + cR, velL + cL);
842de8de29fSMatthew G. Knepley     for (i = 0; i < 2 + dim; i++) flux[i] = 0.5 * ((fL.vals[i] + fR.vals[i]) + speed * (xL[i] - xR[i])) * s2;
843de8de29fSMatthew G. Knepley   } else {
844de8de29fSMatthew G. Knepley     /* int iwave =  */
845de8de29fSMatthew G. Knepley     godunovflux(xL, xR, flux, nn, DIM, gamma);
846de8de29fSMatthew G. Knepley     for (i = 0; i < 2 + dim; i++) flux[i] *= s2;
847de8de29fSMatthew G. Knepley   }
848de8de29fSMatthew G. Knepley   PetscFunctionReturnVoid();
849de8de29fSMatthew G. Knepley }
850de8de29fSMatthew G. Knepley 
851de8de29fSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
852de8de29fSMatthew G. Knepley CEED_QFUNCTION(PhysicsRiemann_Euler_Godunov_CEED)(void *ctx, CeedInt Q, const CeedScalar *const in[], CeedScalar *const out[])
853de8de29fSMatthew G. Knepley {
854de8de29fSMatthew G. Knepley   const CeedScalar    *xL = in[0], *xR = in[1], *geom = in[2];
855de8de29fSMatthew G. Knepley   CeedScalar          *cL = out[0], *cR = out[1];
856de8de29fSMatthew G. Knepley   const Physics_Euler *eu = (Physics_Euler *)ctx;
857de8de29fSMatthew G. Knepley   struct _n_Physics    phys;
858de8de29fSMatthew G. Knepley 
859de8de29fSMatthew G. Knepley   phys.data = (void *)eu;
860de8de29fSMatthew G. Knepley   CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i)
861de8de29fSMatthew G. Knepley   {
862de8de29fSMatthew G. Knepley     const CeedScalar qL[DIM + 2] = {xL[i + Q * 0], xL[i + Q * 1], xL[i + Q * 2], xL[i + Q * 3]};
863de8de29fSMatthew G. Knepley     const CeedScalar qR[DIM + 2] = {xR[i + Q * 0], xR[i + Q * 1], xR[i + Q * 2], xR[i + Q * 3]};
864de8de29fSMatthew G. Knepley     const CeedScalar n[DIM]      = {geom[i + Q * 0], geom[i + Q * 1]};
865de8de29fSMatthew G. Knepley     CeedScalar       flux[DIM + 2];
866de8de29fSMatthew G. Knepley 
867de8de29fSMatthew G. Knepley   #if 0
868de8de29fSMatthew G. Knepley     PetscPrintf(PETSC_COMM_SELF, "Cell %d Normal\n", 0);
869de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < DIM; ++j) {
870de8de29fSMatthew G. Knepley       PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", n[j]);
871de8de29fSMatthew G. Knepley     }
872de8de29fSMatthew G. Knepley     PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: left state\n", 0);
873de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < DIM + 2; ++j) {
874de8de29fSMatthew G. Knepley       PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", qL[j]);
875de8de29fSMatthew G. Knepley     }
876de8de29fSMatthew G. Knepley     PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: right state\n", 0);
877de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < DIM + 2; ++j) {
878de8de29fSMatthew G. Knepley       PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", qR[j]);
879de8de29fSMatthew G. Knepley     }
880de8de29fSMatthew G. Knepley   #endif
881de8de29fSMatthew G. Knepley     PhysicsRiemann_Euler_Godunov(DIM, DIM + 2, NULL, n, qL, qR, 0, NULL, flux, &phys);
882de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < DIM + 2; ++j) {
883de8de29fSMatthew G. Knepley       cL[i + Q * j] = -flux[j] / geom[i + Q * 2];
884de8de29fSMatthew G. Knepley       cR[i + Q * j] = flux[j] / geom[i + Q * 3];
885de8de29fSMatthew G. Knepley     }
886de8de29fSMatthew G. Knepley   #if 0
887de8de29fSMatthew G. Knepley     PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: left flux\n", 0);
888de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < DIM + 2; ++j) {
889de8de29fSMatthew G. Knepley       PetscPrintf(PETSC_COMM_SELF, "  | %g | (%g)\n", cL[i + Q * j], geom[i + Q * 2]);
890de8de29fSMatthew G. Knepley     }
891de8de29fSMatthew G. Knepley     PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: right flux\n", 0);
892de8de29fSMatthew G. Knepley     for (CeedInt j = 0; j < DIM + 2; ++j) {
893de8de29fSMatthew G. Knepley       PetscPrintf(PETSC_COMM_SELF, "  | %g | (%g)\n", cR[i + Q * j], geom[i + Q * 3]);
894de8de29fSMatthew G. Knepley     }
895de8de29fSMatthew G. Knepley   #endif
896de8de29fSMatthew G. Knepley   }
897de8de29fSMatthew G. Knepley   return CEED_ERROR_SUCCESS;
898de8de29fSMatthew G. Knepley }
899de8de29fSMatthew G. Knepley #endif
900