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