xref: /honee/qfunctions/eulervortex.h (revision 002797a32e7e6bab43fe78a1f2bfae476ddb23d4)
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;
41d8a22b9eSJed 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;
1820df2634dSLeila Ghaffari   case 5: // non-smooth thermal bubble - cylinder
1830df2634dSLeila Ghaffari     P    = 1.;
1840df2634dSLeila Ghaffari     T = 1. - (r < 1. ? S_bubble : 0.);
1850df2634dSLeila Ghaffari     rho  = P / (R*T);
1860df2634dSLeila Ghaffari     u[0] = mean_velocity[0];
1870df2634dSLeila Ghaffari     u[1] = mean_velocity[1];
1880df2634dSLeila Ghaffari 
1890df2634dSLeila Ghaffari     // Assign exact solution
1900df2634dSLeila Ghaffari     q[0] = rho;
1910df2634dSLeila Ghaffari     q[1] = rho * u[0];
1920df2634dSLeila Ghaffari     q[2] = rho * u[1];
1930df2634dSLeila Ghaffari     q[3] = rho * u[2];
1940df2634dSLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
1950df2634dSLeila Ghaffari     break;
196a515125bSLeila Ghaffari   }
197a515125bSLeila Ghaffari   // Return
198a515125bSLeila Ghaffari   return 0;
199a515125bSLeila Ghaffari }
200a515125bSLeila Ghaffari 
201a515125bSLeila Ghaffari // *****************************************************************************
202139613f2SLeila Ghaffari // Helper function for computing flux Jacobian
203139613f2SLeila Ghaffari // *****************************************************************************
204d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5],
205139613f2SLeila Ghaffari     const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
206139613f2SLeila Ghaffari     const CeedScalar gamma) {
207139613f2SLeila Ghaffari   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square
208139613f2SLeila Ghaffari   for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions
209139613f2SLeila Ghaffari     for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix
210139613f2SLeila Ghaffari       dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j];
211139613f2SLeila Ghaffari       for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix
212139613f2SLeila Ghaffari         dF[i][0][k+1]   = ((i==k) ? 1. : 0.);
213139613f2SLeila Ghaffari         dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) +
214139613f2SLeila Ghaffari                           ((i==k) ? u[j] : 0.) -
215139613f2SLeila Ghaffari                           ((i==j) ? u[k] : 0.) * (gamma-1.);
216139613f2SLeila Ghaffari         dF[i][4][k+1]   = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) -
217139613f2SLeila Ghaffari                           (gamma-1.)*u[i]*u[k];
218139613f2SLeila Ghaffari       }
219139613f2SLeila Ghaffari       dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.);
220139613f2SLeila Ghaffari     }
221139613f2SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho);
222139613f2SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
223139613f2SLeila Ghaffari   }
224139613f2SLeila Ghaffari }
225139613f2SLeila Ghaffari 
226139613f2SLeila Ghaffari // *****************************************************************************
227d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant)
228d8a22b9eSJed Brown //   Model from:
229d8a22b9eSJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
230d8a22b9eSJed Brown //
231d8a22b9eSJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
232d8a22b9eSJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
233d8a22b9eSJed Brown //
234d8a22b9eSJed Brown // Where
235d8a22b9eSJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
236d8a22b9eSJed Brown //   h[i]      = 2 length(dxdX[i])
237d8a22b9eSJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
238d8a22b9eSJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
239d8a22b9eSJed Brown //   rho(A[i]) = spectral radius of the convective flux Jacobian i,
240d8a22b9eSJed Brown //               wave speed in direction i
241d8a22b9eSJed Brown // *****************************************************************************
242d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3],
243d8a22b9eSJed Brown                                        const CeedScalar dXdx[3][3], const CeedScalar u[3],
244d8a22b9eSJed Brown                                        const CeedScalar sound_speed, const CeedScalar c_tau) {
245d8a22b9eSJed Brown   for (int i=0; i<3; i++) {
246d8a22b9eSJed Brown     // length of element in direction i
247d8a22b9eSJed Brown     CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] +
248d8a22b9eSJed Brown                             dXdx[2][i]*dXdx[2][i]);
249d8a22b9eSJed Brown     // fastest wave in direction i
250d8a22b9eSJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
251d8a22b9eSJed Brown     Tau_x[i] = c_tau * h / fastest_wave;
252d8a22b9eSJed Brown   }
253d8a22b9eSJed Brown }
254d8a22b9eSJed Brown 
255d8a22b9eSJed Brown // *****************************************************************************
256a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
257a515125bSLeila Ghaffari // *****************************************************************************
258a515125bSLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q,
259a515125bSLeila Ghaffari                          const CeedScalar *const *in, CeedScalar *const *out) {
260a515125bSLeila Ghaffari   // Inputs
261a515125bSLeila Ghaffari   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
262a515125bSLeila Ghaffari 
263a515125bSLeila Ghaffari   // Outputs
264a515125bSLeila Ghaffari   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
265a515125bSLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
266a515125bSLeila Ghaffari 
267a515125bSLeila Ghaffari   CeedPragmaSIMD
268a515125bSLeila Ghaffari   // Quadrature Point Loop
269a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
270a515125bSLeila Ghaffari     const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
271139613f2SLeila Ghaffari     CeedScalar q[5] = {0.};
272a515125bSLeila Ghaffari 
273a515125bSLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
274a515125bSLeila Ghaffari 
275a515125bSLeila Ghaffari     for (CeedInt j=0; j<5; j++)
276a515125bSLeila Ghaffari       q0[j][i] = q[j];
277a515125bSLeila Ghaffari   } // End of Quadrature Point Loop
278a515125bSLeila Ghaffari 
279a515125bSLeila Ghaffari   // Return
280a515125bSLeila Ghaffari   return 0;
281a515125bSLeila Ghaffari }
282a515125bSLeila Ghaffari 
283a515125bSLeila Ghaffari // *****************************************************************************
284a515125bSLeila Ghaffari // This QFunction implements the following formulation of Euler equations
285a515125bSLeila Ghaffari //   with explicit time stepping method
286a515125bSLeila Ghaffari //
287a515125bSLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation
288a515125bSLeila Ghaffari //   form with state variables of density, momentum density, and total
289a515125bSLeila Ghaffari //   energy density.
290a515125bSLeila Ghaffari //
291a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
292a515125bSLeila Ghaffari //   rho - Mass Density
293a515125bSLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
294a515125bSLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
295a515125bSLeila Ghaffari //
296a515125bSLeila Ghaffari // Euler Equations:
297a515125bSLeila Ghaffari //   drho/dt + div( U )                   = 0
298a515125bSLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
299a515125bSLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
300a515125bSLeila Ghaffari //
301a515125bSLeila Ghaffari // Equation of State:
302a515125bSLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
303a515125bSLeila Ghaffari //
304a515125bSLeila Ghaffari // Constants:
305a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
306a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
307a515125bSLeila Ghaffari //   g               ,  Gravity
308a515125bSLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
309a515125bSLeila Ghaffari // *****************************************************************************
310a515125bSLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q,
311a515125bSLeila Ghaffari                       const CeedScalar *const *in, CeedScalar *const *out) {
312a515125bSLeila Ghaffari   // *INDENT-OFF*
313a515125bSLeila Ghaffari   // Inputs
314a515125bSLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
315139613f2SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
316a515125bSLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
317a515125bSLeila Ghaffari   // Outputs
318a515125bSLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
319a515125bSLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
320a515125bSLeila Ghaffari 
321139613f2SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
322d8a22b9eSJed Brown   const CeedScalar c_tau = context->c_tau;
323a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
324a515125bSLeila Ghaffari 
325a515125bSLeila Ghaffari   CeedPragmaSIMD
326a515125bSLeila Ghaffari   // Quadrature Point Loop
327a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
328a515125bSLeila Ghaffari     // *INDENT-OFF*
329a515125bSLeila Ghaffari     // Setup
330a515125bSLeila Ghaffari     // -- Interp in
331a515125bSLeila Ghaffari     const CeedScalar rho        =   q[0][i];
332a515125bSLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
333a515125bSLeila Ghaffari                                     q[2][i] / rho,
334a515125bSLeila Ghaffari                                     q[3][i] / rho
335a515125bSLeila Ghaffari                                    };
336a515125bSLeila Ghaffari     const CeedScalar E          =   q[4][i];
337139613f2SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
338139613f2SLeila Ghaffari                                     dq[1][0][i],
339139613f2SLeila Ghaffari                                     dq[2][0][i]
340139613f2SLeila Ghaffari                                    };
341139613f2SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
342139613f2SLeila Ghaffari                                     dq[1][1][i],
343139613f2SLeila Ghaffari                                     dq[2][1][i]},
344139613f2SLeila Ghaffari                                    {dq[0][2][i],
345139613f2SLeila Ghaffari                                     dq[1][2][i],
346139613f2SLeila Ghaffari                                     dq[2][2][i]},
347139613f2SLeila Ghaffari                                    {dq[0][3][i],
348139613f2SLeila Ghaffari                                     dq[1][3][i],
349139613f2SLeila Ghaffari                                     dq[2][3][i]}
350139613f2SLeila Ghaffari                                   };
351139613f2SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
352139613f2SLeila Ghaffari                                     dq[1][4][i],
353139613f2SLeila Ghaffari                                     dq[2][4][i]
354139613f2SLeila Ghaffari                                    };
355a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
356a515125bSLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
357a515125bSLeila Ghaffari     // -- Interp-to-Grad q_data
358a515125bSLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
359a515125bSLeila Ghaffari     // *INDENT-OFF*
360a515125bSLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
361a515125bSLeila Ghaffari                                     q_data[2][i],
362a515125bSLeila Ghaffari                                     q_data[3][i]},
363a515125bSLeila Ghaffari                                    {q_data[4][i],
364a515125bSLeila Ghaffari                                     q_data[5][i],
365a515125bSLeila Ghaffari                                     q_data[6][i]},
366a515125bSLeila Ghaffari                                    {q_data[7][i],
367a515125bSLeila Ghaffari                                     q_data[8][i],
368a515125bSLeila Ghaffari                                     q_data[9][i]}
369a515125bSLeila Ghaffari                                   };
370a515125bSLeila Ghaffari     // *INDENT-ON*
371139613f2SLeila Ghaffari     // dU/dx
372139613f2SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
373139613f2SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
374139613f2SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
375139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
376139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
377139613f2SLeila Ghaffari       for (int k=0; k<3; k++) {
378139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
379139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
380139613f2SLeila Ghaffari         for (int l=0; l<3; l++) {
381139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
382139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
383139613f2SLeila Ghaffari         }
384139613f2SLeila Ghaffari       }
385139613f2SLeila Ghaffari     }
386139613f2SLeila Ghaffari     // Pressure
387a515125bSLeila Ghaffari     const CeedScalar
388a515125bSLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
389a515125bSLeila Ghaffari     E_internal = E - E_kinetic,
390139613f2SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
391a515125bSLeila Ghaffari 
392a515125bSLeila Ghaffari     // The Physics
393a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
394a515125bSLeila Ghaffari     for (int j=0; j<5; j++) {
395139613f2SLeila Ghaffari       v[j][i] = 0.;
396a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
397139613f2SLeila Ghaffari         dv[k][j][i] = 0.;
398a515125bSLeila Ghaffari     }
399a515125bSLeila Ghaffari 
400a515125bSLeila Ghaffari     // -- Density
401a515125bSLeila Ghaffari     // ---- u rho
402a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
403a515125bSLeila Ghaffari       dv[j][0][i]  += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
404a515125bSLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
405a515125bSLeila Ghaffari     // -- Momentum
406a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
407a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
408a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
409139613f2SLeila Ghaffari         dv[k][j+1][i]  += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
410139613f2SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
411139613f2SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
412a515125bSLeila Ghaffari     // -- Total Energy Density
413a515125bSLeila Ghaffari     // ---- (E + P) u
414a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
415a515125bSLeila Ghaffari       dv[j][4][i]  += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
416a515125bSLeila Ghaffari                                          u[2]*dXdx[j][2]);
417139613f2SLeila Ghaffari 
418139613f2SLeila Ghaffari     // --Stabilization terms
419139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
420139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
421d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
422139613f2SLeila Ghaffari 
423139613f2SLeila Ghaffari     // ---- Transpose of the Jacobian
424139613f2SLeila Ghaffari     CeedScalar jacob_F_conv_T[3][5][5];
425139613f2SLeila Ghaffari     for (int j=0; j<3; j++)
426139613f2SLeila Ghaffari       for (int k=0; k<5; k++)
427139613f2SLeila Ghaffari         for (int l=0; l<5; l++)
428139613f2SLeila Ghaffari           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
429139613f2SLeila Ghaffari 
430139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
431139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
432139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
433139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
434139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
435139613f2SLeila Ghaffari       for (int k=0; k<3; k++)
436139613f2SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
437139613f2SLeila Ghaffari     }
438139613f2SLeila Ghaffari 
439139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
440139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
441139613f2SLeila Ghaffari     for (int j=0; j<3; j++)
442139613f2SLeila Ghaffari       for (int k=0; k<5; k++)
443139613f2SLeila Ghaffari         for (int l=0; l<5; l++)
444139613f2SLeila Ghaffari           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
445139613f2SLeila Ghaffari 
446d8a22b9eSJed Brown     // Stabilization
447d8a22b9eSJed Brown     // -- Tau elements
448d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
449d8a22b9eSJed Brown     CeedScalar Tau_x[3] = {0.};
450d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
451139613f2SLeila Ghaffari 
452d8a22b9eSJed Brown     // -- Stabilization method: none or SU
453139613f2SLeila Ghaffari     CeedScalar stab[5][3];
454139613f2SLeila Ghaffari     switch (context->stabilization) {
455139613f2SLeila Ghaffari     case 0:        // Galerkin
456139613f2SLeila Ghaffari       break;
457139613f2SLeila Ghaffari     case 1:        // SU
458139613f2SLeila Ghaffari       for (int j=0; j<3; j++)
459139613f2SLeila Ghaffari         for (int k=0; k<5; k++)
460139613f2SLeila Ghaffari           for (int l=0; l<5; l++)
461d8a22b9eSJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
462139613f2SLeila Ghaffari 
463139613f2SLeila Ghaffari       for (int j=0; j<5; j++)
464139613f2SLeila Ghaffari         for (int k=0; k<3; k++)
465139613f2SLeila Ghaffari           dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] +
466139613f2SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
467139613f2SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
468139613f2SLeila Ghaffari       break;
469139613f2SLeila Ghaffari     case 2:        // SUPG is not implemented for explicit scheme
470139613f2SLeila Ghaffari       break;
471139613f2SLeila Ghaffari     }
472139613f2SLeila Ghaffari 
473a515125bSLeila Ghaffari   } // End Quadrature Point Loop
474a515125bSLeila Ghaffari 
475a515125bSLeila Ghaffari   // Return
476a515125bSLeila Ghaffari   return 0;
477a515125bSLeila Ghaffari }
478a515125bSLeila Ghaffari 
479a515125bSLeila Ghaffari // *****************************************************************************
480a515125bSLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above)
481a515125bSLeila Ghaffari //   with implicit time stepping method
482a515125bSLeila Ghaffari //
483a515125bSLeila Ghaffari // *****************************************************************************
484a515125bSLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q,
485a515125bSLeila Ghaffari                                 const CeedScalar *const *in, CeedScalar *const *out) {
486a515125bSLeila Ghaffari   // *INDENT-OFF*
487a515125bSLeila Ghaffari   // Inputs
488a515125bSLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
489139613f2SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
490a515125bSLeila Ghaffari                    (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
491a515125bSLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
492a515125bSLeila Ghaffari   // Outputs
493a515125bSLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
494a515125bSLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
495a515125bSLeila Ghaffari 
496139613f2SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
497d8a22b9eSJed Brown   const CeedScalar c_tau = context->c_tau;
498a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
499a515125bSLeila Ghaffari 
500a515125bSLeila Ghaffari   CeedPragmaSIMD
501a515125bSLeila Ghaffari   // Quadrature Point Loop
502a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
503a515125bSLeila Ghaffari     // *INDENT-OFF*
504a515125bSLeila Ghaffari     // Setup
505a515125bSLeila Ghaffari     // -- Interp in
506a515125bSLeila Ghaffari     const CeedScalar rho        =   q[0][i];
507a515125bSLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
508a515125bSLeila Ghaffari                                     q[2][i] / rho,
509a515125bSLeila Ghaffari                                     q[3][i] / rho
510a515125bSLeila Ghaffari                                    };
511a515125bSLeila Ghaffari     const CeedScalar E          =   q[4][i];
512139613f2SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
513139613f2SLeila Ghaffari                                     dq[1][0][i],
514139613f2SLeila Ghaffari                                     dq[2][0][i]
515139613f2SLeila Ghaffari                                    };
516139613f2SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
517139613f2SLeila Ghaffari                                     dq[1][1][i],
518139613f2SLeila Ghaffari                                     dq[2][1][i]},
519139613f2SLeila Ghaffari                                    {dq[0][2][i],
520139613f2SLeila Ghaffari                                     dq[1][2][i],
521139613f2SLeila Ghaffari                                     dq[2][2][i]},
522139613f2SLeila Ghaffari                                    {dq[0][3][i],
523139613f2SLeila Ghaffari                                     dq[1][3][i],
524139613f2SLeila Ghaffari                                     dq[2][3][i]}
525139613f2SLeila Ghaffari                                   };
526139613f2SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
527139613f2SLeila Ghaffari                                     dq[1][4][i],
528139613f2SLeila Ghaffari                                     dq[2][4][i]
529139613f2SLeila Ghaffari                                    };
530a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
531a515125bSLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
532a515125bSLeila Ghaffari     // -- Interp-to-Grad q_data
533a515125bSLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
534a515125bSLeila Ghaffari     // *INDENT-OFF*
535a515125bSLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
536a515125bSLeila Ghaffari                                     q_data[2][i],
537a515125bSLeila Ghaffari                                     q_data[3][i]},
538a515125bSLeila Ghaffari                                    {q_data[4][i],
539a515125bSLeila Ghaffari                                     q_data[5][i],
540a515125bSLeila Ghaffari                                     q_data[6][i]},
541a515125bSLeila Ghaffari                                    {q_data[7][i],
542a515125bSLeila Ghaffari                                     q_data[8][i],
543a515125bSLeila Ghaffari                                     q_data[9][i]}
544a515125bSLeila Ghaffari                                   };
545a515125bSLeila Ghaffari     // *INDENT-ON*
546139613f2SLeila Ghaffari     // dU/dx
547139613f2SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
548139613f2SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
549139613f2SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
550139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
551139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
552139613f2SLeila Ghaffari       for (int k=0; k<3; k++) {
553139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
554139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
555139613f2SLeila Ghaffari         for (int l=0; l<3; l++) {
556139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
557139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
558139613f2SLeila Ghaffari         }
559139613f2SLeila Ghaffari       }
560139613f2SLeila Ghaffari     }
561a515125bSLeila Ghaffari     const CeedScalar
562a515125bSLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
563a515125bSLeila Ghaffari     E_internal = E - E_kinetic,
564139613f2SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
565a515125bSLeila Ghaffari 
566a515125bSLeila Ghaffari     // The Physics
567a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
568a515125bSLeila Ghaffari     for (int j=0; j<5; j++) {
569139613f2SLeila Ghaffari       v[j][i] = 0.;
570a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
571139613f2SLeila Ghaffari         dv[k][j][i] = 0.;
572a515125bSLeila Ghaffari     }
573a515125bSLeila Ghaffari     //-----mass matrix
574a515125bSLeila Ghaffari     for (int j=0; j<5; j++)
575a515125bSLeila Ghaffari       v[j][i] += wdetJ*q_dot[j][i];
576a515125bSLeila Ghaffari 
577a515125bSLeila Ghaffari     // -- Density
578a515125bSLeila Ghaffari     // ---- u rho
579a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
580a515125bSLeila Ghaffari       dv[j][0][i]  -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
581a515125bSLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
582a515125bSLeila Ghaffari     // -- Momentum
583a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
584a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
585a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
586139613f2SLeila Ghaffari         dv[k][j+1][i]  -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
587139613f2SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
588139613f2SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
589a515125bSLeila Ghaffari     // -- Total Energy Density
590a515125bSLeila Ghaffari     // ---- (E + P) u
591a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
592a515125bSLeila Ghaffari       dv[j][4][i]  -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
593a515125bSLeila Ghaffari                                          u[2]*dXdx[j][2]);
594139613f2SLeila Ghaffari 
595139613f2SLeila Ghaffari     // -- Stabilization terms
596139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
597139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
598d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
599139613f2SLeila Ghaffari 
600139613f2SLeila Ghaffari     // ---- Transpose of the Jacobian
601139613f2SLeila Ghaffari     CeedScalar jacob_F_conv_T[3][5][5];
602139613f2SLeila Ghaffari     for (int j=0; j<3; j++)
603139613f2SLeila Ghaffari       for (int k=0; k<5; k++)
604139613f2SLeila Ghaffari         for (int l=0; l<5; l++)
605139613f2SLeila Ghaffari           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
606139613f2SLeila Ghaffari 
607139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
608139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
609139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
610139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
611139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
612139613f2SLeila Ghaffari       for (int k=0; k<3; k++)
613139613f2SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
614139613f2SLeila Ghaffari     }
615139613f2SLeila Ghaffari 
616139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
617139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
618139613f2SLeila Ghaffari     for (int j=0; j<3; j++)
619139613f2SLeila Ghaffari       for (int k=0; k<5; k++)
620139613f2SLeila Ghaffari         for (int l=0; l<5; l++)
621139613f2SLeila Ghaffari           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
622139613f2SLeila Ghaffari 
623139613f2SLeila Ghaffari     // ---- Strong residual
624139613f2SLeila Ghaffari     CeedScalar strong_res[5];
625139613f2SLeila Ghaffari     for (int j=0; j<5; j++)
626139613f2SLeila Ghaffari       strong_res[j] = q_dot[j][i] + strong_conv[j];
627139613f2SLeila Ghaffari 
628d8a22b9eSJed Brown     // Stabilization
629d8a22b9eSJed Brown     // -- Tau elements
630d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
631d8a22b9eSJed Brown     CeedScalar Tau_x[3] = {0.};
632d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
633139613f2SLeila Ghaffari 
634d8a22b9eSJed Brown     // -- Stabilization method: none, SU, or SUPG
635139613f2SLeila Ghaffari     CeedScalar stab[5][3];
636139613f2SLeila Ghaffari     switch (context->stabilization) {
637139613f2SLeila Ghaffari     case 0:        // Galerkin
638139613f2SLeila Ghaffari       break;
639139613f2SLeila Ghaffari     case 1:        // SU
640139613f2SLeila Ghaffari       for (int j=0; j<3; j++)
641139613f2SLeila Ghaffari         for (int k=0; k<5; k++)
642139613f2SLeila Ghaffari           for (int l=0; l<5; l++)
643d8a22b9eSJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
644139613f2SLeila Ghaffari 
645139613f2SLeila Ghaffari       for (int j=0; j<5; j++)
646139613f2SLeila Ghaffari         for (int k=0; k<3; k++)
647139613f2SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
648139613f2SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
649139613f2SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
650139613f2SLeila Ghaffari       break;
651139613f2SLeila Ghaffari     case 2:        // SUPG
652139613f2SLeila Ghaffari       for (int j=0; j<3; j++)
653139613f2SLeila Ghaffari         for (int k=0; k<5; k++)
654139613f2SLeila Ghaffari           for (int l=0; l<5; l++)
655d8a22b9eSJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_res[l];
656139613f2SLeila Ghaffari 
657139613f2SLeila Ghaffari       for (int j=0; j<5; j++)
658139613f2SLeila Ghaffari         for (int k=0; k<3; k++)
659139613f2SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
660139613f2SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
661139613f2SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
662139613f2SLeila Ghaffari       break;
663139613f2SLeila Ghaffari     }
664a515125bSLeila Ghaffari   } // End Quadrature Point Loop
665a515125bSLeila Ghaffari 
666a515125bSLeila Ghaffari   // Return
667a515125bSLeila Ghaffari   return 0;
668a515125bSLeila Ghaffari }
669a515125bSLeila Ghaffari // *****************************************************************************
670*002797a3SLeila Ghaffari // This QFunction sets the inflow boundary conditions for
671*002797a3SLeila Ghaffari //   the traveling vortex problem.
672a515125bSLeila Ghaffari //
673a515125bSLeila Ghaffari //  Prescribed T_inlet and P_inlet are converted to conservative variables
674a515125bSLeila Ghaffari //      and applied weakly.
675a515125bSLeila Ghaffari //
676a515125bSLeila Ghaffari // *****************************************************************************
677*002797a3SLeila Ghaffari CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q,
678a515125bSLeila Ghaffari                                        const CeedScalar *const *in,
679a515125bSLeila Ghaffari                                        CeedScalar *const *out) {
680a515125bSLeila Ghaffari   // *INDENT-OFF*
681a515125bSLeila Ghaffari   // Inputs
682*002797a3SLeila Ghaffari   const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
683a515125bSLeila Ghaffari   // Outputs
684a515125bSLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
685a515125bSLeila Ghaffari   // *INDENT-ON*
686a515125bSLeila Ghaffari   EulerContext context = (EulerContext)ctx;
687a515125bSLeila Ghaffari   const int euler_test      = context->euler_test;
688a515125bSLeila Ghaffari   const bool implicit       = context->implicit;
689a515125bSLeila Ghaffari   CeedScalar *mean_velocity = context->mean_velocity;
690a515125bSLeila Ghaffari   const CeedScalar cv    = 2.5;
691a515125bSLeila Ghaffari   const CeedScalar R     = 1.;
692a515125bSLeila Ghaffari   CeedScalar T_inlet;
693a515125bSLeila Ghaffari   CeedScalar P_inlet;
694a515125bSLeila Ghaffari 
695a515125bSLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
696a515125bSLeila Ghaffari   if (euler_test == 1 || euler_test == 3)
697a515125bSLeila Ghaffari     for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.;
698a515125bSLeila Ghaffari 
699a515125bSLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
700a515125bSLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
701a515125bSLeila Ghaffari   else T_inlet = P_inlet = 1.;
702a515125bSLeila Ghaffari 
703a515125bSLeila Ghaffari   CeedPragmaSIMD
704a515125bSLeila Ghaffari   // Quadrature Point Loop
705a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
706a515125bSLeila Ghaffari     // Setup
707a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
708a515125bSLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
709a515125bSLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
710a515125bSLeila Ghaffari     // We can effect this by swapping the sign on this weight
711a515125bSLeila Ghaffari     const CeedScalar wdetJb  = (implicit ? -1. : 1.) * q_data_sur[0][i];
712*002797a3SLeila Ghaffari     // ---- Normal vect
713a515125bSLeila Ghaffari     const CeedScalar norm[3] = {q_data_sur[1][i],
714a515125bSLeila Ghaffari                                 q_data_sur[2][i],
715a515125bSLeila Ghaffari                                 q_data_sur[3][i]
716a515125bSLeila Ghaffari                                };
717a515125bSLeila Ghaffari 
718a515125bSLeila Ghaffari     // face_normal = Normal vector of the face
719a515125bSLeila Ghaffari     const CeedScalar face_normal = norm[0]*mean_velocity[0] +
720a515125bSLeila Ghaffari                                    norm[1]*mean_velocity[1] +
721a515125bSLeila Ghaffari                                    norm[2]*mean_velocity[2];
722a515125bSLeila Ghaffari     // The Physics
723a515125bSLeila Ghaffari     // Zero v so all future terms can safely sum into it
724139613f2SLeila Ghaffari     for (int j=0; j<5; j++) v[j][i] = 0.;
725a515125bSLeila Ghaffari 
726a515125bSLeila Ghaffari     // Implementing in/outflow BCs
727*002797a3SLeila Ghaffari     if (face_normal > 0) {
728a515125bSLeila Ghaffari     } else { // inflow
729a515125bSLeila Ghaffari       const CeedScalar rho_inlet = P_inlet/(R*T_inlet);
730a515125bSLeila Ghaffari       const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] +
731a515125bSLeila Ghaffari                                           mean_velocity[1]*mean_velocity[1]) / 2.;
732a515125bSLeila Ghaffari       // incoming total energy
733a515125bSLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
734a515125bSLeila Ghaffari 
735a515125bSLeila Ghaffari       // The Physics
736a515125bSLeila Ghaffari       // -- Density
737a515125bSLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
738a515125bSLeila Ghaffari 
739a515125bSLeila Ghaffari       // -- Momentum
740a515125bSLeila Ghaffari       for (int j=0; j<3; j++)
741a515125bSLeila Ghaffari         v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] +
742a515125bSLeila Ghaffari                               norm[j] * P_inlet);
743a515125bSLeila Ghaffari 
744a515125bSLeila Ghaffari       // -- Total Energy Density
745a515125bSLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
746a515125bSLeila Ghaffari     }
747a515125bSLeila Ghaffari 
748a515125bSLeila Ghaffari   } // End Quadrature Point Loop
749a515125bSLeila Ghaffari   return 0;
750a515125bSLeila Ghaffari }
751a515125bSLeila Ghaffari 
752a515125bSLeila Ghaffari // *****************************************************************************
753a515125bSLeila Ghaffari 
754a515125bSLeila Ghaffari #endif // eulervortex_h
755