xref: /honee/qfunctions/eulervortex.h (revision d8a22b9effbf5cac417f760aeb29833c41e3abc6)
1a515125bSLeila Ghaffari // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2a515125bSLeila Ghaffari // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3a515125bSLeila Ghaffari // reserved. See files LICENSE and NOTICE for details.
4a515125bSLeila Ghaffari //
5a515125bSLeila Ghaffari // This file is part of CEED, a collection of benchmarks, miniapps, software
6a515125bSLeila Ghaffari // libraries and APIs for efficient high-order finite element and spectral
7a515125bSLeila Ghaffari // element discretizations for exascale applications. For more information and
8a515125bSLeila Ghaffari // source code availability see http://github.com/ceed.
9a515125bSLeila Ghaffari //
10a515125bSLeila Ghaffari // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11a515125bSLeila Ghaffari // a collaborative effort of two U.S. Department of Energy organizations (Office
12a515125bSLeila Ghaffari // of Science and the National Nuclear Security Administration) responsible for
13a515125bSLeila Ghaffari // the planning and preparation of a capable exascale ecosystem, including
14a515125bSLeila Ghaffari // software, applications, hardware, advanced system engineering and early
15a515125bSLeila Ghaffari // testbed platforms, in support of the nation's exascale computing imperative.
16a515125bSLeila Ghaffari 
17a515125bSLeila Ghaffari /// @file
18a515125bSLeila Ghaffari /// Euler traveling vortex initial condition and operator for Navier-Stokes
19a515125bSLeila Ghaffari /// example using PETSc
20a515125bSLeila Ghaffari 
21a515125bSLeila Ghaffari // Model from:
22a515125bSLeila Ghaffari //   On the Order of Accuracy and Numerical Performance of Two Classes of
23a515125bSLeila Ghaffari //   Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011).
24a515125bSLeila Ghaffari 
25a515125bSLeila Ghaffari #ifndef eulervortex_h
26a515125bSLeila Ghaffari #define eulervortex_h
27a515125bSLeila Ghaffari 
28a515125bSLeila Ghaffari #include <math.h>
29a515125bSLeila Ghaffari 
30a515125bSLeila Ghaffari #ifndef M_PI
31a515125bSLeila Ghaffari #define M_PI    3.14159265358979323846
32a515125bSLeila Ghaffari #endif
33a515125bSLeila Ghaffari 
34a515125bSLeila Ghaffari #ifndef euler_context_struct
35a515125bSLeila Ghaffari #define euler_context_struct
36a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext;
37a515125bSLeila Ghaffari struct EulerContext_ {
38a515125bSLeila Ghaffari   CeedScalar center[3];
39a515125bSLeila Ghaffari   CeedScalar curr_time;
40a515125bSLeila Ghaffari   CeedScalar vortex_strength;
41*d8a22b9eSJed Brown   CeedScalar c_tau;
42a515125bSLeila Ghaffari   CeedScalar mean_velocity[3];
43a515125bSLeila Ghaffari   bool implicit;
44139613f2SLeila Ghaffari   int euler_test;
45139613f2SLeila Ghaffari   int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG
46a515125bSLeila Ghaffari };
47a515125bSLeila Ghaffari #endif
48a515125bSLeila Ghaffari 
49a515125bSLeila Ghaffari // *****************************************************************************
50a515125bSLeila Ghaffari // This function sets the initial conditions
51a515125bSLeila Ghaffari //
52a515125bSLeila Ghaffari //   Temperature:
53a515125bSLeila Ghaffari //     T   = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
54a515125bSLeila Ghaffari //   Density:
55a515125bSLeila Ghaffari //     rho = (T/S_vortex)^(1 / (gamma - 1))
56a515125bSLeila Ghaffari //   Pressure:
57a515125bSLeila Ghaffari //     P   = rho * T
58a515125bSLeila Ghaffari //   Velocity:
59a515125bSLeila Ghaffari //     ui  = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
60a515125bSLeila Ghaffari //     r   = sqrt( (x - xc)**2 + (y - yc)**2 )
61a515125bSLeila Ghaffari //   Velocity/Momentum Density:
62a515125bSLeila Ghaffari //     Ui  = rho ui
63a515125bSLeila Ghaffari //   Total Energy:
64a515125bSLeila Ghaffari //     E   = P / (gamma - 1) + rho (u u)/2
65a515125bSLeila Ghaffari //
66a515125bSLeila Ghaffari // Constants:
67a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
68a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
69a515125bSLeila Ghaffari //   vortex_strength ,  Strength of vortex
70a515125bSLeila Ghaffari //   center          ,  Location of bubble center
71a515125bSLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
72a515125bSLeila Ghaffari //
73a515125bSLeila Ghaffari // *****************************************************************************
74a515125bSLeila Ghaffari 
75a515125bSLeila Ghaffari // *****************************************************************************
76a515125bSLeila Ghaffari // This helper function provides support for the exact, time-dependent solution
77a515125bSLeila Ghaffari //   (currently not implemented) and IC formulation for Euler traveling vortex
78a515125bSLeila Ghaffari // *****************************************************************************
79a515125bSLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time,
80a515125bSLeila Ghaffari                                       const CeedScalar X[], CeedInt Nf, CeedScalar q[],
81a515125bSLeila Ghaffari                                       void *ctx) {
82a515125bSLeila Ghaffari   // Context
83a515125bSLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
84a515125bSLeila Ghaffari   const CeedScalar vortex_strength    = context->vortex_strength;
85a515125bSLeila Ghaffari   const CeedScalar *center            = context->center; // Center of the domain
86a515125bSLeila Ghaffari   const CeedScalar *mean_velocity = context->mean_velocity;
87a515125bSLeila Ghaffari 
88a515125bSLeila Ghaffari   // Setup
89a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
90a515125bSLeila Ghaffari   const CeedScalar cv    = 2.5;
91a515125bSLeila Ghaffari   const CeedScalar R     = 1.;
92a515125bSLeila Ghaffari   const CeedScalar x     = X[0], y = X[1]; // Coordinates
93a515125bSLeila Ghaffari   // Vortex center
94a515125bSLeila Ghaffari   const CeedScalar xc = center[0] + mean_velocity[0] * time;
95a515125bSLeila Ghaffari   const CeedScalar yc = center[1] + mean_velocity[1] * time;
96a515125bSLeila Ghaffari 
97a515125bSLeila Ghaffari   const CeedScalar x0       = x - xc;
98a515125bSLeila Ghaffari   const CeedScalar y0       = y - yc;
99a515125bSLeila Ghaffari   const CeedScalar r        = sqrt( x0*x0 + y0*y0 );
100a515125bSLeila Ghaffari   const CeedScalar C        = vortex_strength * exp((1. - r*r)/2.) / (2. * M_PI);
101139613f2SLeila Ghaffari   const CeedScalar delta_T  = - (gamma - 1.) * vortex_strength * vortex_strength *
102139613f2SLeila Ghaffari                               exp(1 - r*r) / (8. * gamma * M_PI * M_PI);
103a515125bSLeila Ghaffari   const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma
104a515125bSLeila Ghaffari   const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength /
105a515125bSLeila Ghaffari                               (8.*gamma*M_PI*M_PI);
106a515125bSLeila Ghaffari   CeedScalar rho, P, T, E, u[3] = {0.};
107a515125bSLeila Ghaffari 
108a515125bSLeila Ghaffari   // Initial Conditions
109a515125bSLeila Ghaffari   switch (context->euler_test) {
110a515125bSLeila Ghaffari   case 0: // Traveling vortex
111a515125bSLeila Ghaffari     T = 1 + delta_T;
112a515125bSLeila Ghaffari     // P = rho * T
113a515125bSLeila Ghaffari     // P = S * rho^gamma
114a515125bSLeila Ghaffari     // Solve for rho, then substitute for P
115139613f2SLeila Ghaffari     rho  = pow(T/S_vortex, 1 / (gamma - 1.));
116a515125bSLeila Ghaffari     P    = rho * T;
117a515125bSLeila Ghaffari     u[0] = mean_velocity[0] - C*y0;
118a515125bSLeila Ghaffari     u[1] = mean_velocity[1] + C*x0;
119a515125bSLeila Ghaffari 
120a515125bSLeila Ghaffari     // Assign exact solution
121a515125bSLeila Ghaffari     q[0] = rho;
122a515125bSLeila Ghaffari     q[1] = rho * u[0];
123a515125bSLeila Ghaffari     q[2] = rho * u[1];
124a515125bSLeila Ghaffari     q[3] = rho * u[2];
125a515125bSLeila Ghaffari     q[4] = P / (gamma - 1.) + rho * (u[0]*u[0] + u[1]*u[1]) / 2.;
126a515125bSLeila Ghaffari     break;
127a515125bSLeila Ghaffari   case 1: // Constant zero velocity, density constant, total energy constant
128a515125bSLeila Ghaffari     rho  = 1.;
129a515125bSLeila Ghaffari     E    = 2.;
130a515125bSLeila Ghaffari 
131a515125bSLeila Ghaffari     // Assign exact solution
132a515125bSLeila Ghaffari     q[0] = rho;
133a515125bSLeila Ghaffari     q[1] = rho * u[0];
134a515125bSLeila Ghaffari     q[2] = rho * u[1];
135a515125bSLeila Ghaffari     q[3] = rho * u[2];
136a515125bSLeila Ghaffari     q[4] = E;
137a515125bSLeila Ghaffari     break;
138a515125bSLeila Ghaffari   case 2: // Constant nonzero velocity, density constant, total energy constant
139a515125bSLeila Ghaffari     rho  = 1.;
140a515125bSLeila Ghaffari     E    = 2.;
141a515125bSLeila Ghaffari     u[0] = mean_velocity[0];
142a515125bSLeila Ghaffari     u[1] = mean_velocity[1];
143a515125bSLeila Ghaffari 
144a515125bSLeila Ghaffari     // Assign exact solution
145a515125bSLeila Ghaffari     q[0] = rho;
146a515125bSLeila Ghaffari     q[1] = rho * u[0];
147a515125bSLeila Ghaffari     q[2] = rho * u[1];
148a515125bSLeila Ghaffari     q[3] = rho * u[2];
149a515125bSLeila Ghaffari     q[4] = E;
150a515125bSLeila Ghaffari     break;
151a515125bSLeila Ghaffari   case 3: // Velocity zero, pressure constant
152a515125bSLeila Ghaffari     // (so density and internal energy will be non-constant),
153a515125bSLeila Ghaffari     // but the velocity should stay zero and the bubble won't diffuse
154a515125bSLeila Ghaffari     // (for Euler, where there is no thermal conductivity)
155a515125bSLeila Ghaffari     P    = 1.;
156a515125bSLeila Ghaffari     T    = 1. - S_bubble * exp(1. - r*r);
157a515125bSLeila Ghaffari     rho  = P / (R*T);
158a515125bSLeila Ghaffari 
159a515125bSLeila Ghaffari     // Assign exact solution
160a515125bSLeila Ghaffari     q[0] = rho;
161a515125bSLeila Ghaffari     q[1] = rho * u[0];
162a515125bSLeila Ghaffari     q[2] = rho * u[1];
163a515125bSLeila Ghaffari     q[3] = rho * u[2];
164a515125bSLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
165a515125bSLeila Ghaffari     break;
166a515125bSLeila Ghaffari   case 4: // Constant nonzero velocity, pressure constant
167a515125bSLeila Ghaffari     // (so density and internal energy will be non-constant),
168a515125bSLeila Ghaffari     // it should be transported across the domain, but velocity stays constant
169a515125bSLeila Ghaffari     P    = 1.;
170a515125bSLeila Ghaffari     T    = 1. - S_bubble * exp(1. - r*r);
171a515125bSLeila Ghaffari     rho  = P / (R*T);
172a515125bSLeila Ghaffari     u[0] = mean_velocity[0];
173a515125bSLeila Ghaffari     u[1] = mean_velocity[1];
174a515125bSLeila Ghaffari 
175a515125bSLeila Ghaffari     // Assign exact solution
176a515125bSLeila Ghaffari     q[0] = rho;
177a515125bSLeila Ghaffari     q[1] = rho * u[0];
178a515125bSLeila Ghaffari     q[2] = rho * u[1];
179a515125bSLeila Ghaffari     q[3] = rho * u[2];
180a515125bSLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
181a515125bSLeila Ghaffari     break;
182a515125bSLeila Ghaffari   }
183a515125bSLeila Ghaffari   // Return
184a515125bSLeila Ghaffari   return 0;
185a515125bSLeila Ghaffari }
186a515125bSLeila Ghaffari 
187a515125bSLeila Ghaffari // *****************************************************************************
188139613f2SLeila Ghaffari // Helper function for computing flux Jacobian
189139613f2SLeila Ghaffari // *****************************************************************************
190*d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5],
191139613f2SLeila Ghaffari     const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
192139613f2SLeila Ghaffari     const CeedScalar gamma) {
193139613f2SLeila Ghaffari   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square
194139613f2SLeila Ghaffari   for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions
195139613f2SLeila Ghaffari     for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix
196139613f2SLeila Ghaffari       dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j];
197139613f2SLeila Ghaffari       for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix
198139613f2SLeila Ghaffari         dF[i][0][k+1]   = ((i==k) ? 1. : 0.);
199139613f2SLeila Ghaffari         dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) +
200139613f2SLeila Ghaffari                           ((i==k) ? u[j] : 0.) -
201139613f2SLeila Ghaffari                           ((i==j) ? u[k] : 0.) * (gamma-1.);
202139613f2SLeila Ghaffari         dF[i][4][k+1]   = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) -
203139613f2SLeila Ghaffari                           (gamma-1.)*u[i]*u[k];
204139613f2SLeila Ghaffari       }
205139613f2SLeila Ghaffari       dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.);
206139613f2SLeila Ghaffari     }
207139613f2SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho);
208139613f2SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
209139613f2SLeila Ghaffari   }
210139613f2SLeila Ghaffari }
211139613f2SLeila Ghaffari 
212139613f2SLeila Ghaffari // *****************************************************************************
213*d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant)
214*d8a22b9eSJed Brown //   Model from:
215*d8a22b9eSJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
216*d8a22b9eSJed Brown //
217*d8a22b9eSJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
218*d8a22b9eSJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
219*d8a22b9eSJed Brown //
220*d8a22b9eSJed Brown // Where
221*d8a22b9eSJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
222*d8a22b9eSJed Brown //   h[i]      = 2 length(dxdX[i])
223*d8a22b9eSJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
224*d8a22b9eSJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
225*d8a22b9eSJed Brown //   rho(A[i]) = spectral radius of the convective flux Jacobian i,
226*d8a22b9eSJed Brown //               wave speed in direction i
227*d8a22b9eSJed Brown // *****************************************************************************
228*d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3],
229*d8a22b9eSJed Brown                                        const CeedScalar dXdx[3][3], const CeedScalar u[3],
230*d8a22b9eSJed Brown                                        const CeedScalar sound_speed, const CeedScalar c_tau) {
231*d8a22b9eSJed Brown   for (int i=0; i<3; i++) {
232*d8a22b9eSJed Brown     // length of element in direction i
233*d8a22b9eSJed Brown     CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] +
234*d8a22b9eSJed Brown                             dXdx[2][i]*dXdx[2][i]);
235*d8a22b9eSJed Brown     // fastest wave in direction i
236*d8a22b9eSJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
237*d8a22b9eSJed Brown     Tau_x[i] = c_tau * h / fastest_wave;
238*d8a22b9eSJed Brown   }
239*d8a22b9eSJed Brown }
240*d8a22b9eSJed Brown 
241*d8a22b9eSJed Brown // *****************************************************************************
242a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
243a515125bSLeila Ghaffari // *****************************************************************************
244a515125bSLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q,
245a515125bSLeila Ghaffari                          const CeedScalar *const *in, CeedScalar *const *out) {
246a515125bSLeila Ghaffari   // Inputs
247a515125bSLeila Ghaffari   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
248a515125bSLeila Ghaffari 
249a515125bSLeila Ghaffari   // Outputs
250a515125bSLeila Ghaffari   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
251a515125bSLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
252a515125bSLeila Ghaffari 
253a515125bSLeila Ghaffari   CeedPragmaSIMD
254a515125bSLeila Ghaffari   // Quadrature Point Loop
255a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
256a515125bSLeila Ghaffari     const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
257139613f2SLeila Ghaffari     CeedScalar q[5] = {0.};
258a515125bSLeila Ghaffari 
259a515125bSLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
260a515125bSLeila Ghaffari 
261a515125bSLeila Ghaffari     for (CeedInt j=0; j<5; j++)
262a515125bSLeila Ghaffari       q0[j][i] = q[j];
263a515125bSLeila Ghaffari   } // End of Quadrature Point Loop
264a515125bSLeila Ghaffari 
265a515125bSLeila Ghaffari   // Return
266a515125bSLeila Ghaffari   return 0;
267a515125bSLeila Ghaffari }
268a515125bSLeila Ghaffari 
269a515125bSLeila Ghaffari // *****************************************************************************
270a515125bSLeila Ghaffari // This QFunction implements the following formulation of Euler equations
271a515125bSLeila Ghaffari //   with explicit time stepping method
272a515125bSLeila Ghaffari //
273a515125bSLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation
274a515125bSLeila Ghaffari //   form with state variables of density, momentum density, and total
275a515125bSLeila Ghaffari //   energy density.
276a515125bSLeila Ghaffari //
277a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
278a515125bSLeila Ghaffari //   rho - Mass Density
279a515125bSLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
280a515125bSLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
281a515125bSLeila Ghaffari //
282a515125bSLeila Ghaffari // Euler Equations:
283a515125bSLeila Ghaffari //   drho/dt + div( U )                   = 0
284a515125bSLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
285a515125bSLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
286a515125bSLeila Ghaffari //
287a515125bSLeila Ghaffari // Equation of State:
288a515125bSLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
289a515125bSLeila Ghaffari //
290a515125bSLeila Ghaffari // Constants:
291a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
292a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
293a515125bSLeila Ghaffari //   g               ,  Gravity
294a515125bSLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
295a515125bSLeila Ghaffari // *****************************************************************************
296a515125bSLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q,
297a515125bSLeila Ghaffari                       const CeedScalar *const *in, CeedScalar *const *out) {
298a515125bSLeila Ghaffari   // *INDENT-OFF*
299a515125bSLeila Ghaffari   // Inputs
300a515125bSLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
301139613f2SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
302a515125bSLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
303a515125bSLeila Ghaffari   // Outputs
304a515125bSLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
305a515125bSLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
306a515125bSLeila Ghaffari 
307139613f2SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
308*d8a22b9eSJed Brown   const CeedScalar c_tau = context->c_tau;
309a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
310a515125bSLeila Ghaffari 
311a515125bSLeila Ghaffari   CeedPragmaSIMD
312a515125bSLeila Ghaffari   // Quadrature Point Loop
313a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
314a515125bSLeila Ghaffari     // *INDENT-OFF*
315a515125bSLeila Ghaffari     // Setup
316a515125bSLeila Ghaffari     // -- Interp in
317a515125bSLeila Ghaffari     const CeedScalar rho        =   q[0][i];
318a515125bSLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
319a515125bSLeila Ghaffari                                     q[2][i] / rho,
320a515125bSLeila Ghaffari                                     q[3][i] / rho
321a515125bSLeila Ghaffari                                    };
322a515125bSLeila Ghaffari     const CeedScalar E          =   q[4][i];
323139613f2SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
324139613f2SLeila Ghaffari                                     dq[1][0][i],
325139613f2SLeila Ghaffari                                     dq[2][0][i]
326139613f2SLeila Ghaffari                                    };
327139613f2SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
328139613f2SLeila Ghaffari                                     dq[1][1][i],
329139613f2SLeila Ghaffari                                     dq[2][1][i]},
330139613f2SLeila Ghaffari                                    {dq[0][2][i],
331139613f2SLeila Ghaffari                                     dq[1][2][i],
332139613f2SLeila Ghaffari                                     dq[2][2][i]},
333139613f2SLeila Ghaffari                                    {dq[0][3][i],
334139613f2SLeila Ghaffari                                     dq[1][3][i],
335139613f2SLeila Ghaffari                                     dq[2][3][i]}
336139613f2SLeila Ghaffari                                   };
337139613f2SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
338139613f2SLeila Ghaffari                                     dq[1][4][i],
339139613f2SLeila Ghaffari                                     dq[2][4][i]
340139613f2SLeila Ghaffari                                    };
341a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
342a515125bSLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
343a515125bSLeila Ghaffari     // -- Interp-to-Grad q_data
344a515125bSLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
345a515125bSLeila Ghaffari     // *INDENT-OFF*
346a515125bSLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
347a515125bSLeila Ghaffari                                     q_data[2][i],
348a515125bSLeila Ghaffari                                     q_data[3][i]},
349a515125bSLeila Ghaffari                                    {q_data[4][i],
350a515125bSLeila Ghaffari                                     q_data[5][i],
351a515125bSLeila Ghaffari                                     q_data[6][i]},
352a515125bSLeila Ghaffari                                    {q_data[7][i],
353a515125bSLeila Ghaffari                                     q_data[8][i],
354a515125bSLeila Ghaffari                                     q_data[9][i]}
355a515125bSLeila Ghaffari                                   };
356a515125bSLeila Ghaffari     // *INDENT-ON*
357139613f2SLeila Ghaffari     // dU/dx
358139613f2SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
359139613f2SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
360139613f2SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
361139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
362139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
363139613f2SLeila Ghaffari       for (int k=0; k<3; k++) {
364139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
365139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
366139613f2SLeila Ghaffari         for (int l=0; l<3; l++) {
367139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
368139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
369139613f2SLeila Ghaffari         }
370139613f2SLeila Ghaffari       }
371139613f2SLeila Ghaffari     }
372139613f2SLeila Ghaffari     // Pressure
373a515125bSLeila Ghaffari     const CeedScalar
374a515125bSLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
375a515125bSLeila Ghaffari     E_internal = E - E_kinetic,
376139613f2SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
377a515125bSLeila Ghaffari 
378a515125bSLeila Ghaffari     // The Physics
379a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
380a515125bSLeila Ghaffari     for (int j=0; j<5; j++) {
381139613f2SLeila Ghaffari       v[j][i] = 0.;
382a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
383139613f2SLeila Ghaffari         dv[k][j][i] = 0.;
384a515125bSLeila Ghaffari     }
385a515125bSLeila Ghaffari 
386a515125bSLeila Ghaffari     // -- Density
387a515125bSLeila Ghaffari     // ---- u rho
388a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
389a515125bSLeila Ghaffari       dv[j][0][i]  += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
390a515125bSLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
391a515125bSLeila Ghaffari     // -- Momentum
392a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
393a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
394a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
395139613f2SLeila Ghaffari         dv[k][j+1][i]  += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
396139613f2SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
397139613f2SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
398a515125bSLeila Ghaffari     // -- Total Energy Density
399a515125bSLeila Ghaffari     // ---- (E + P) u
400a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
401a515125bSLeila Ghaffari       dv[j][4][i]  += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
402a515125bSLeila Ghaffari                                          u[2]*dXdx[j][2]);
403139613f2SLeila Ghaffari 
404139613f2SLeila Ghaffari     // --Stabilization terms
405139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
406139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
407*d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
408139613f2SLeila Ghaffari 
409139613f2SLeila Ghaffari     // ---- Transpose of the Jacobian
410139613f2SLeila Ghaffari     CeedScalar jacob_F_conv_T[3][5][5];
411139613f2SLeila Ghaffari     for (int j=0; j<3; j++)
412139613f2SLeila Ghaffari       for (int k=0; k<5; k++)
413139613f2SLeila Ghaffari         for (int l=0; l<5; l++)
414139613f2SLeila Ghaffari           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
415139613f2SLeila Ghaffari 
416139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
417139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
418139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
419139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
420139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
421139613f2SLeila Ghaffari       for (int k=0; k<3; k++)
422139613f2SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
423139613f2SLeila Ghaffari     }
424139613f2SLeila Ghaffari 
425139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
426139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
427139613f2SLeila Ghaffari     for (int j=0; j<3; j++)
428139613f2SLeila Ghaffari       for (int k=0; k<5; k++)
429139613f2SLeila Ghaffari         for (int l=0; l<5; l++)
430139613f2SLeila Ghaffari           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
431139613f2SLeila Ghaffari 
432*d8a22b9eSJed Brown     // Stabilization
433*d8a22b9eSJed Brown     // -- Tau elements
434*d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
435*d8a22b9eSJed Brown     CeedScalar Tau_x[3] = {0.};
436*d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
437139613f2SLeila Ghaffari 
438*d8a22b9eSJed Brown     // -- Stabilization method: none or SU
439139613f2SLeila Ghaffari     CeedScalar stab[5][3];
440139613f2SLeila Ghaffari     switch (context->stabilization) {
441139613f2SLeila Ghaffari     case 0:        // Galerkin
442139613f2SLeila Ghaffari       break;
443139613f2SLeila Ghaffari     case 1:        // SU
444139613f2SLeila Ghaffari       for (int j=0; j<3; j++)
445139613f2SLeila Ghaffari         for (int k=0; k<5; k++)
446139613f2SLeila Ghaffari           for (int l=0; l<5; l++)
447*d8a22b9eSJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
448139613f2SLeila Ghaffari 
449139613f2SLeila Ghaffari       for (int j=0; j<5; j++)
450139613f2SLeila Ghaffari         for (int k=0; k<3; k++)
451139613f2SLeila Ghaffari           dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] +
452139613f2SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
453139613f2SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
454139613f2SLeila Ghaffari       break;
455139613f2SLeila Ghaffari     case 2:        // SUPG is not implemented for explicit scheme
456139613f2SLeila Ghaffari       break;
457139613f2SLeila Ghaffari     }
458139613f2SLeila Ghaffari 
459a515125bSLeila Ghaffari   } // End Quadrature Point Loop
460a515125bSLeila Ghaffari 
461a515125bSLeila Ghaffari   // Return
462a515125bSLeila Ghaffari   return 0;
463a515125bSLeila Ghaffari }
464a515125bSLeila Ghaffari 
465a515125bSLeila Ghaffari // *****************************************************************************
466a515125bSLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above)
467a515125bSLeila Ghaffari //   with implicit time stepping method
468a515125bSLeila Ghaffari //
469a515125bSLeila Ghaffari // *****************************************************************************
470a515125bSLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q,
471a515125bSLeila Ghaffari                                 const CeedScalar *const *in, CeedScalar *const *out) {
472a515125bSLeila Ghaffari   // *INDENT-OFF*
473a515125bSLeila Ghaffari   // Inputs
474a515125bSLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
475139613f2SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
476a515125bSLeila Ghaffari                    (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
477a515125bSLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
478a515125bSLeila Ghaffari   // Outputs
479a515125bSLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
480a515125bSLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
481a515125bSLeila Ghaffari 
482139613f2SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
483*d8a22b9eSJed Brown   const CeedScalar c_tau = context->c_tau;
484a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
485a515125bSLeila Ghaffari 
486a515125bSLeila Ghaffari   CeedPragmaSIMD
487a515125bSLeila Ghaffari   // Quadrature Point Loop
488a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
489a515125bSLeila Ghaffari     // *INDENT-OFF*
490a515125bSLeila Ghaffari     // Setup
491a515125bSLeila Ghaffari     // -- Interp in
492a515125bSLeila Ghaffari     const CeedScalar rho        =   q[0][i];
493a515125bSLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
494a515125bSLeila Ghaffari                                     q[2][i] / rho,
495a515125bSLeila Ghaffari                                     q[3][i] / rho
496a515125bSLeila Ghaffari                                    };
497a515125bSLeila Ghaffari     const CeedScalar E          =   q[4][i];
498139613f2SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
499139613f2SLeila Ghaffari                                     dq[1][0][i],
500139613f2SLeila Ghaffari                                     dq[2][0][i]
501139613f2SLeila Ghaffari                                    };
502139613f2SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
503139613f2SLeila Ghaffari                                     dq[1][1][i],
504139613f2SLeila Ghaffari                                     dq[2][1][i]},
505139613f2SLeila Ghaffari                                    {dq[0][2][i],
506139613f2SLeila Ghaffari                                     dq[1][2][i],
507139613f2SLeila Ghaffari                                     dq[2][2][i]},
508139613f2SLeila Ghaffari                                    {dq[0][3][i],
509139613f2SLeila Ghaffari                                     dq[1][3][i],
510139613f2SLeila Ghaffari                                     dq[2][3][i]}
511139613f2SLeila Ghaffari                                   };
512139613f2SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
513139613f2SLeila Ghaffari                                     dq[1][4][i],
514139613f2SLeila Ghaffari                                     dq[2][4][i]
515139613f2SLeila Ghaffari                                    };
516a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
517a515125bSLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
518a515125bSLeila Ghaffari     // -- Interp-to-Grad q_data
519a515125bSLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
520a515125bSLeila Ghaffari     // *INDENT-OFF*
521a515125bSLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
522a515125bSLeila Ghaffari                                     q_data[2][i],
523a515125bSLeila Ghaffari                                     q_data[3][i]},
524a515125bSLeila Ghaffari                                    {q_data[4][i],
525a515125bSLeila Ghaffari                                     q_data[5][i],
526a515125bSLeila Ghaffari                                     q_data[6][i]},
527a515125bSLeila Ghaffari                                    {q_data[7][i],
528a515125bSLeila Ghaffari                                     q_data[8][i],
529a515125bSLeila Ghaffari                                     q_data[9][i]}
530a515125bSLeila Ghaffari                                   };
531a515125bSLeila Ghaffari     // *INDENT-ON*
532139613f2SLeila Ghaffari     // dU/dx
533139613f2SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
534139613f2SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
535139613f2SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
536139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
537139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
538139613f2SLeila Ghaffari       for (int k=0; k<3; k++) {
539139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
540139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
541139613f2SLeila Ghaffari         for (int l=0; l<3; l++) {
542139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
543139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
544139613f2SLeila Ghaffari         }
545139613f2SLeila Ghaffari       }
546139613f2SLeila Ghaffari     }
547a515125bSLeila Ghaffari     const CeedScalar
548a515125bSLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
549a515125bSLeila Ghaffari     E_internal = E - E_kinetic,
550139613f2SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
551a515125bSLeila Ghaffari 
552a515125bSLeila Ghaffari     // The Physics
553a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
554a515125bSLeila Ghaffari     for (int j=0; j<5; j++) {
555139613f2SLeila Ghaffari       v[j][i] = 0.;
556a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
557139613f2SLeila Ghaffari         dv[k][j][i] = 0.;
558a515125bSLeila Ghaffari     }
559a515125bSLeila Ghaffari     //-----mass matrix
560a515125bSLeila Ghaffari     for (int j=0; j<5; j++)
561a515125bSLeila Ghaffari       v[j][i] += wdetJ*q_dot[j][i];
562a515125bSLeila Ghaffari 
563a515125bSLeila Ghaffari     // -- Density
564a515125bSLeila Ghaffari     // ---- u rho
565a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
566a515125bSLeila Ghaffari       dv[j][0][i]  -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
567a515125bSLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
568a515125bSLeila Ghaffari     // -- Momentum
569a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
570a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
571a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
572139613f2SLeila Ghaffari         dv[k][j+1][i]  -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
573139613f2SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
574139613f2SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
575a515125bSLeila Ghaffari     // -- Total Energy Density
576a515125bSLeila Ghaffari     // ---- (E + P) u
577a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
578a515125bSLeila Ghaffari       dv[j][4][i]  -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
579a515125bSLeila Ghaffari                                          u[2]*dXdx[j][2]);
580139613f2SLeila Ghaffari 
581139613f2SLeila Ghaffari     // -- Stabilization terms
582139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
583139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
584*d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
585139613f2SLeila Ghaffari 
586139613f2SLeila Ghaffari     // ---- Transpose of the Jacobian
587139613f2SLeila Ghaffari     CeedScalar jacob_F_conv_T[3][5][5];
588139613f2SLeila Ghaffari     for (int j=0; j<3; j++)
589139613f2SLeila Ghaffari       for (int k=0; k<5; k++)
590139613f2SLeila Ghaffari         for (int l=0; l<5; l++)
591139613f2SLeila Ghaffari           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
592139613f2SLeila Ghaffari 
593139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
594139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
595139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
596139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
597139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
598139613f2SLeila Ghaffari       for (int k=0; k<3; k++)
599139613f2SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
600139613f2SLeila Ghaffari     }
601139613f2SLeila Ghaffari 
602139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
603139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
604139613f2SLeila Ghaffari     for (int j=0; j<3; j++)
605139613f2SLeila Ghaffari       for (int k=0; k<5; k++)
606139613f2SLeila Ghaffari         for (int l=0; l<5; l++)
607139613f2SLeila Ghaffari           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
608139613f2SLeila Ghaffari 
609139613f2SLeila Ghaffari     // ---- Strong residual
610139613f2SLeila Ghaffari     CeedScalar strong_res[5];
611139613f2SLeila Ghaffari     for (int j=0; j<5; j++)
612139613f2SLeila Ghaffari       strong_res[j] = q_dot[j][i] + strong_conv[j];
613139613f2SLeila Ghaffari 
614*d8a22b9eSJed Brown     // Stabilization
615*d8a22b9eSJed Brown     // -- Tau elements
616*d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
617*d8a22b9eSJed Brown     CeedScalar Tau_x[3] = {0.};
618*d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
619139613f2SLeila Ghaffari 
620*d8a22b9eSJed Brown     // -- Stabilization method: none, SU, or SUPG
621139613f2SLeila Ghaffari     CeedScalar stab[5][3];
622139613f2SLeila Ghaffari     switch (context->stabilization) {
623139613f2SLeila Ghaffari     case 0:        // Galerkin
624139613f2SLeila Ghaffari       break;
625139613f2SLeila Ghaffari     case 1:        // SU
626139613f2SLeila Ghaffari       for (int j=0; j<3; j++)
627139613f2SLeila Ghaffari         for (int k=0; k<5; k++)
628139613f2SLeila Ghaffari           for (int l=0; l<5; l++)
629*d8a22b9eSJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
630139613f2SLeila Ghaffari 
631139613f2SLeila Ghaffari       for (int j=0; j<5; j++)
632139613f2SLeila Ghaffari         for (int k=0; k<3; k++)
633139613f2SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
634139613f2SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
635139613f2SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
636139613f2SLeila Ghaffari       break;
637139613f2SLeila Ghaffari     case 2:        // SUPG
638139613f2SLeila Ghaffari       for (int j=0; j<3; j++)
639139613f2SLeila Ghaffari         for (int k=0; k<5; k++)
640139613f2SLeila Ghaffari           for (int l=0; l<5; l++)
641*d8a22b9eSJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_res[l];
642139613f2SLeila Ghaffari 
643139613f2SLeila Ghaffari       for (int j=0; j<5; j++)
644139613f2SLeila Ghaffari         for (int k=0; k<3; k++)
645139613f2SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
646139613f2SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
647139613f2SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
648139613f2SLeila Ghaffari       break;
649139613f2SLeila Ghaffari     }
650a515125bSLeila Ghaffari   } // End Quadrature Point Loop
651a515125bSLeila Ghaffari 
652a515125bSLeila Ghaffari   // Return
653a515125bSLeila Ghaffari   return 0;
654a515125bSLeila Ghaffari }
655a515125bSLeila Ghaffari // *****************************************************************************
656a515125bSLeila Ghaffari // This QFunction sets the boundary conditions
657a515125bSLeila Ghaffari //   In this problem, only in/outflow BCs are implemented
658a515125bSLeila Ghaffari //
659a515125bSLeila Ghaffari //  Inflow and outflow faces are determined based on
660a515125bSLeila Ghaffari //    sign(dot(mean_velocity, normal)):
661a515125bSLeila Ghaffari //      sign(dot(mean_velocity, normal)) > 0 : outflow BCs
662a515125bSLeila Ghaffari //      sign(dot(mean_velocity, normal)) < 0 : inflow BCs
663a515125bSLeila Ghaffari //
664a515125bSLeila Ghaffari //  Outflow BCs:
665a515125bSLeila Ghaffari //    The validity of the weak form of the governing equations is
666a515125bSLeila Ghaffari //      extended to the outflow.
667a515125bSLeila Ghaffari //
668a515125bSLeila Ghaffari //  Inflow BCs:
669a515125bSLeila Ghaffari //    Prescribed T_inlet and P_inlet are converted to conservative variables
670a515125bSLeila Ghaffari //      and applied weakly.
671a515125bSLeila Ghaffari //
672a515125bSLeila Ghaffari // *****************************************************************************
673a515125bSLeila Ghaffari CEED_QFUNCTION(Euler_Sur)(void *ctx, CeedInt Q,
674a515125bSLeila Ghaffari                           const CeedScalar *const *in,
675a515125bSLeila Ghaffari                           CeedScalar *const *out) {
676a515125bSLeila Ghaffari   // *INDENT-OFF*
677a515125bSLeila Ghaffari   // Inputs
678a515125bSLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
679a515125bSLeila Ghaffari                    (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
680a515125bSLeila Ghaffari   // Outputs
681a515125bSLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
682a515125bSLeila Ghaffari   // *INDENT-ON*
683a515125bSLeila Ghaffari   EulerContext context = (EulerContext)ctx;
684a515125bSLeila Ghaffari   const int euler_test      = context->euler_test;
685a515125bSLeila Ghaffari   const bool implicit       = context->implicit;
686a515125bSLeila Ghaffari   CeedScalar *mean_velocity = context->mean_velocity;
687a515125bSLeila Ghaffari 
688a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
689a515125bSLeila Ghaffari   const CeedScalar cv    = 2.5;
690a515125bSLeila Ghaffari   const CeedScalar R     = 1.;
691a515125bSLeila Ghaffari   CeedScalar T_inlet;
692a515125bSLeila Ghaffari   CeedScalar P_inlet;
693a515125bSLeila Ghaffari 
694a515125bSLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
695a515125bSLeila Ghaffari   if (euler_test == 1 || euler_test == 3)
696a515125bSLeila Ghaffari     for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.;
697a515125bSLeila Ghaffari 
698a515125bSLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
699a515125bSLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
700a515125bSLeila Ghaffari   else T_inlet = P_inlet = 1.;
701a515125bSLeila Ghaffari 
702a515125bSLeila Ghaffari   CeedPragmaSIMD
703a515125bSLeila Ghaffari   // Quadrature Point Loop
704a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
705a515125bSLeila Ghaffari     // Setup
706a515125bSLeila Ghaffari     // -- Interp in
707a515125bSLeila Ghaffari     const CeedScalar rho      =  q[0][i];
708a515125bSLeila Ghaffari     const CeedScalar u[3]     = {q[1][i] / rho,
709a515125bSLeila Ghaffari                                  q[2][i] / rho,
710a515125bSLeila Ghaffari                                  q[3][i] / rho
711a515125bSLeila Ghaffari                                 };
712a515125bSLeila Ghaffari     const CeedScalar E        =  q[4][i];
713a515125bSLeila Ghaffari 
714a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
715a515125bSLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
716a515125bSLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
717a515125bSLeila Ghaffari     // We can effect this by swapping the sign on this weight
718a515125bSLeila Ghaffari     const CeedScalar wdetJb     =   (implicit ? -1. : 1.) * q_data_sur[0][i];
719a515125bSLeila Ghaffari     // ---- Normal vectors
720a515125bSLeila Ghaffari     const CeedScalar norm[3]    =   {q_data_sur[1][i],
721a515125bSLeila Ghaffari                                      q_data_sur[2][i],
722a515125bSLeila Ghaffari                                      q_data_sur[3][i]
723a515125bSLeila Ghaffari                                     };
724a515125bSLeila Ghaffari 
725a515125bSLeila Ghaffari     // face_normal = Normal vector of the face
726a515125bSLeila Ghaffari     const CeedScalar face_normal = norm[0]*mean_velocity[0] +
727a515125bSLeila Ghaffari                                    norm[1]*mean_velocity[1] +
728a515125bSLeila Ghaffari                                    norm[2]*mean_velocity[2];
729a515125bSLeila Ghaffari     // The Physics
730a515125bSLeila Ghaffari     // Zero v so all future terms can safely sum into it
731139613f2SLeila Ghaffari     for (int j=0; j<5; j++) v[j][i] = 0.;
732a515125bSLeila Ghaffari 
733a515125bSLeila Ghaffari     // Implementing in/outflow BCs
734a515125bSLeila Ghaffari     if (face_normal > 0) { // outflow
735a515125bSLeila Ghaffari       const CeedScalar E_kinetic = (u[0]*u[0] + u[1]*u[1]) / 2.;
736a515125bSLeila Ghaffari       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.); // pressure
737a515125bSLeila Ghaffari       const CeedScalar u_normal  = norm[0]*u[0] + norm[1]*u[1] +
738a515125bSLeila Ghaffari                                    norm[2]*u[2]; // Normal velocity
739a515125bSLeila Ghaffari       // The Physics
740a515125bSLeila Ghaffari       // -- Density
741a515125bSLeila Ghaffari       v[0][i] -= wdetJb * rho * u_normal;
742a515125bSLeila Ghaffari 
743a515125bSLeila Ghaffari       // -- Momentum
744a515125bSLeila Ghaffari       for (int j=0; j<3; j++)
745a515125bSLeila Ghaffari         v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + norm[j] * P);
746a515125bSLeila Ghaffari 
747a515125bSLeila Ghaffari       // -- Total Energy Density
748a515125bSLeila Ghaffari       v[4][i] -= wdetJb * u_normal * (E + P);
749a515125bSLeila Ghaffari 
750a515125bSLeila Ghaffari     } else { // inflow
751a515125bSLeila Ghaffari       const CeedScalar rho_inlet = P_inlet/(R*T_inlet);
752a515125bSLeila Ghaffari       const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] +
753a515125bSLeila Ghaffari                                           mean_velocity[1]*mean_velocity[1]) / 2.;
754a515125bSLeila Ghaffari       // incoming total energy
755a515125bSLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
756a515125bSLeila Ghaffari 
757a515125bSLeila Ghaffari       // The Physics
758a515125bSLeila Ghaffari       // -- Density
759a515125bSLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
760a515125bSLeila Ghaffari 
761a515125bSLeila Ghaffari       // -- Momentum
762a515125bSLeila Ghaffari       for (int j=0; j<3; j++)
763a515125bSLeila Ghaffari         v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] +
764a515125bSLeila Ghaffari                               norm[j] * P_inlet);
765a515125bSLeila Ghaffari 
766a515125bSLeila Ghaffari       // -- Total Energy Density
767a515125bSLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
768a515125bSLeila Ghaffari     }
769a515125bSLeila Ghaffari 
770a515125bSLeila Ghaffari   } // End Quadrature Point Loop
771a515125bSLeila Ghaffari   return 0;
772a515125bSLeila Ghaffari }
773a515125bSLeila Ghaffari 
774a515125bSLeila Ghaffari // *****************************************************************************
775a515125bSLeila Ghaffari 
776a515125bSLeila Ghaffari #endif // eulervortex_h
777