xref: /honee/qfunctions/eulervortex.h (revision 3d65b166467e11bf054ebd404213121dbd593270)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Euler traveling vortex initial condition and operator for Navier-Stokes
10a515125bSLeila Ghaffari /// example using PETSc
11a515125bSLeila Ghaffari 
12a515125bSLeila Ghaffari // Model from:
13a515125bSLeila Ghaffari //   On the Order of Accuracy and Numerical Performance of Two Classes of
14a515125bSLeila Ghaffari //   Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011).
15a515125bSLeila Ghaffari 
16a515125bSLeila Ghaffari #ifndef eulervortex_h
17a515125bSLeila Ghaffari #define eulervortex_h
18a515125bSLeila Ghaffari 
193a8779fbSJames Wright #include <ceed.h>
20d0cce58aSJeremy L Thompson #include <math.h>
212b916ea7SJeremy L Thompson 
22704b8bbeSJames Wright #include "utils.h"
23a515125bSLeila Ghaffari 
24a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext;
25a515125bSLeila Ghaffari struct EulerContext_ {
26a515125bSLeila Ghaffari   CeedScalar center[3];
27a515125bSLeila Ghaffari   CeedScalar curr_time;
28a515125bSLeila Ghaffari   CeedScalar vortex_strength;
29d8a22b9eSJed Brown   CeedScalar c_tau;
30a515125bSLeila Ghaffari   CeedScalar mean_velocity[3];
31a515125bSLeila Ghaffari   bool       implicit;
32139613f2SLeila Ghaffari   int        euler_test;
33139613f2SLeila Ghaffari   int        stabilization;  // See StabilizationType: 0=none, 1=SU, 2=SUPG
34a515125bSLeila Ghaffari };
35a515125bSLeila Ghaffari 
36a515125bSLeila Ghaffari // *****************************************************************************
37a515125bSLeila Ghaffari // This function sets the initial conditions
38a515125bSLeila Ghaffari //
39a515125bSLeila Ghaffari //   Temperature:
40a515125bSLeila Ghaffari //     T   = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
41a515125bSLeila Ghaffari //   Density:
42a515125bSLeila Ghaffari //     rho = (T/S_vortex)^(1 / (gamma - 1))
43a515125bSLeila Ghaffari //   Pressure:
44a515125bSLeila Ghaffari //     P   = rho * T
45a515125bSLeila Ghaffari //   Velocity:
46a515125bSLeila Ghaffari //     ui  = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
47a515125bSLeila Ghaffari //     r   = sqrt( (x - xc)**2 + (y - yc)**2 )
48a515125bSLeila Ghaffari //   Velocity/Momentum Density:
49a515125bSLeila Ghaffari //     Ui  = rho ui
50a515125bSLeila Ghaffari //   Total Energy:
51a515125bSLeila Ghaffari //     E   = P / (gamma - 1) + rho (u u)/2
52a515125bSLeila Ghaffari //
53a515125bSLeila Ghaffari // Constants:
54a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
55a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
56a515125bSLeila Ghaffari //   vortex_strength ,  Strength of vortex
57a515125bSLeila Ghaffari //   center          ,  Location of bubble center
58a515125bSLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
59a515125bSLeila Ghaffari //
60a515125bSLeila Ghaffari // *****************************************************************************
61a515125bSLeila Ghaffari 
62a515125bSLeila Ghaffari // *****************************************************************************
63a515125bSLeila Ghaffari // This helper function provides support for the exact, time-dependent solution
64a515125bSLeila Ghaffari //   (currently not implemented) and IC formulation for Euler traveling vortex
65a515125bSLeila Ghaffari // *****************************************************************************
662b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
67a515125bSLeila Ghaffari   // Context
68a515125bSLeila Ghaffari   const EulerContext context         = (EulerContext)ctx;
69a515125bSLeila Ghaffari   const CeedScalar   vortex_strength = context->vortex_strength;
70a515125bSLeila Ghaffari   const CeedScalar  *center          = context->center;  // Center of the domain
71a515125bSLeila Ghaffari   const CeedScalar  *mean_velocity   = context->mean_velocity;
72a515125bSLeila Ghaffari 
73a515125bSLeila Ghaffari   // Setup
74a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
75a515125bSLeila Ghaffari   const CeedScalar cv    = 2.5;
76a515125bSLeila Ghaffari   const CeedScalar R     = 1.;
77a515125bSLeila Ghaffari   const CeedScalar x = X[0], y = X[1];  // Coordinates
78a515125bSLeila Ghaffari   // Vortex center
79a515125bSLeila Ghaffari   const CeedScalar xc = center[0] + mean_velocity[0] * time;
80a515125bSLeila Ghaffari   const CeedScalar yc = center[1] + mean_velocity[1] * time;
81a515125bSLeila Ghaffari 
82a515125bSLeila Ghaffari   const CeedScalar x0       = x - xc;
83a515125bSLeila Ghaffari   const CeedScalar y0       = y - yc;
84a515125bSLeila Ghaffari   const CeedScalar r        = sqrt(x0 * x0 + y0 * y0);
85a515125bSLeila Ghaffari   const CeedScalar C        = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI);
862b916ea7SJeremy L Thompson   const CeedScalar delta_T  = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI);
87a515125bSLeila Ghaffari   const CeedScalar S_vortex = 1;  // no perturbation in the entropy P / rho^gamma
882b916ea7SJeremy L Thompson   const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI);
89a515125bSLeila Ghaffari   CeedScalar       rho, P, T, E, u[3] = {0.};
90a515125bSLeila Ghaffari 
91a515125bSLeila Ghaffari   // Initial Conditions
92a515125bSLeila Ghaffari   switch (context->euler_test) {
93a515125bSLeila Ghaffari     case 0:  // Traveling vortex
94a515125bSLeila Ghaffari       T = 1 + delta_T;
95a515125bSLeila Ghaffari       // P = rho * T
96a515125bSLeila Ghaffari       // P = S * rho^gamma
97a515125bSLeila Ghaffari       // Solve for rho, then substitute for P
98139613f2SLeila Ghaffari       rho  = pow(T / S_vortex, 1 / (gamma - 1.));
99a515125bSLeila Ghaffari       P    = rho * T;
100a515125bSLeila Ghaffari       u[0] = mean_velocity[0] - C * y0;
101a515125bSLeila Ghaffari       u[1] = mean_velocity[1] + C * x0;
102a515125bSLeila Ghaffari 
103a515125bSLeila Ghaffari       // Assign exact solution
104a515125bSLeila Ghaffari       q[0] = rho;
105a515125bSLeila Ghaffari       q[1] = rho * u[0];
106a515125bSLeila Ghaffari       q[2] = rho * u[1];
107a515125bSLeila Ghaffari       q[3] = rho * u[2];
108a515125bSLeila Ghaffari       q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.;
109a515125bSLeila Ghaffari       break;
110a515125bSLeila Ghaffari     case 1:  // Constant zero velocity, density constant, total energy constant
111a515125bSLeila Ghaffari       rho = 1.;
112a515125bSLeila Ghaffari       E   = 2.;
113a515125bSLeila Ghaffari 
114a515125bSLeila Ghaffari       // Assign exact solution
115a515125bSLeila Ghaffari       q[0] = rho;
116a515125bSLeila Ghaffari       q[1] = rho * u[0];
117a515125bSLeila Ghaffari       q[2] = rho * u[1];
118a515125bSLeila Ghaffari       q[3] = rho * u[2];
119a515125bSLeila Ghaffari       q[4] = E;
120a515125bSLeila Ghaffari       break;
121a515125bSLeila Ghaffari     case 2:  // Constant nonzero velocity, density constant, total energy constant
122a515125bSLeila Ghaffari       rho  = 1.;
123a515125bSLeila Ghaffari       E    = 2.;
124a515125bSLeila Ghaffari       u[0] = mean_velocity[0];
125a515125bSLeila Ghaffari       u[1] = mean_velocity[1];
126a515125bSLeila Ghaffari 
127a515125bSLeila Ghaffari       // Assign exact solution
128a515125bSLeila Ghaffari       q[0] = rho;
129a515125bSLeila Ghaffari       q[1] = rho * u[0];
130a515125bSLeila Ghaffari       q[2] = rho * u[1];
131a515125bSLeila Ghaffari       q[3] = rho * u[2];
132a515125bSLeila Ghaffari       q[4] = E;
133a515125bSLeila Ghaffari       break;
134a515125bSLeila Ghaffari     case 3:  // Velocity zero, pressure constant
135a515125bSLeila Ghaffari       // (so density and internal energy will be non-constant),
136a515125bSLeila Ghaffari       // but the velocity should stay zero and the bubble won't diffuse
137a515125bSLeila Ghaffari       // (for Euler, where there is no thermal conductivity)
138a515125bSLeila Ghaffari       P   = 1.;
139a515125bSLeila Ghaffari       T   = 1. - S_bubble * exp(1. - r * r);
140a515125bSLeila Ghaffari       rho = P / (R * T);
141a515125bSLeila Ghaffari 
142a515125bSLeila Ghaffari       // Assign exact solution
143a515125bSLeila Ghaffari       q[0] = rho;
144a515125bSLeila Ghaffari       q[1] = rho * u[0];
145a515125bSLeila Ghaffari       q[2] = rho * u[1];
146a515125bSLeila Ghaffari       q[3] = rho * u[2];
147a515125bSLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
148a515125bSLeila Ghaffari       break;
149a515125bSLeila Ghaffari     case 4:  // Constant nonzero velocity, pressure constant
150a515125bSLeila Ghaffari       // (so density and internal energy will be non-constant),
151a515125bSLeila Ghaffari       // it should be transported across the domain, but velocity stays constant
152a515125bSLeila Ghaffari       P    = 1.;
153a515125bSLeila Ghaffari       T    = 1. - S_bubble * exp(1. - r * r);
154a515125bSLeila Ghaffari       rho  = P / (R * T);
155a515125bSLeila Ghaffari       u[0] = mean_velocity[0];
156a515125bSLeila Ghaffari       u[1] = mean_velocity[1];
157a515125bSLeila Ghaffari 
158a515125bSLeila Ghaffari       // Assign exact solution
159a515125bSLeila Ghaffari       q[0] = rho;
160a515125bSLeila Ghaffari       q[1] = rho * u[0];
161a515125bSLeila Ghaffari       q[2] = rho * u[1];
162a515125bSLeila Ghaffari       q[3] = rho * u[2];
163a515125bSLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
164a515125bSLeila Ghaffari       break;
1650df2634dSLeila Ghaffari     case 5:  // non-smooth thermal bubble - cylinder
1660df2634dSLeila Ghaffari       P    = 1.;
1670df2634dSLeila Ghaffari       T    = 1. - (r < 1. ? S_bubble : 0.);
1680df2634dSLeila Ghaffari       rho  = P / (R * T);
1690df2634dSLeila Ghaffari       u[0] = mean_velocity[0];
1700df2634dSLeila Ghaffari       u[1] = mean_velocity[1];
1710df2634dSLeila Ghaffari 
1720df2634dSLeila Ghaffari       // Assign exact solution
1730df2634dSLeila Ghaffari       q[0] = rho;
1740df2634dSLeila Ghaffari       q[1] = rho * u[0];
1750df2634dSLeila Ghaffari       q[2] = rho * u[1];
1760df2634dSLeila Ghaffari       q[3] = rho * u[2];
1770df2634dSLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
1780df2634dSLeila Ghaffari       break;
179a515125bSLeila Ghaffari   }
180a515125bSLeila Ghaffari   // Return
181a515125bSLeila Ghaffari   return 0;
182a515125bSLeila Ghaffari }
183a515125bSLeila Ghaffari 
184a515125bSLeila Ghaffari // *****************************************************************************
185139613f2SLeila Ghaffari // Helper function for computing flux Jacobian
186139613f2SLeila Ghaffari // *****************************************************************************
1872b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
188139613f2SLeila Ghaffari                                                         const CeedScalar gamma) {
189139613f2SLeila Ghaffari   CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2];  // Velocity square
190139613f2SLeila Ghaffari   for (CeedInt i = 0; i < 3; i++) {                           // Jacobian matrices for 3 directions
191139613f2SLeila Ghaffari     for (CeedInt j = 0; j < 3; j++) {                         // Rows of each Jacobian matrix
192139613f2SLeila Ghaffari       dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j];
193139613f2SLeila Ghaffari       for (CeedInt k = 0; k < 3; k++) {  // Columns of each Jacobian matrix
194139613f2SLeila Ghaffari         dF[i][0][k + 1]     = ((i == k) ? 1. : 0.);
1952b916ea7SJeremy L Thompson         dF[i][j + 1][k + 1] = ((j == k) ? u[i] : 0.) + ((i == k) ? u[j] : 0.) - ((i == j) ? u[k] : 0.) * (gamma - 1.);
1962b916ea7SJeremy L Thompson         dF[i][4][k + 1]     = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k];
197139613f2SLeila Ghaffari       }
198139613f2SLeila Ghaffari       dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.);
199139613f2SLeila Ghaffari     }
200139613f2SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho);
201139613f2SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
202139613f2SLeila Ghaffari   }
203139613f2SLeila Ghaffari }
204139613f2SLeila Ghaffari 
205139613f2SLeila Ghaffari // *****************************************************************************
206d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant)
207d8a22b9eSJed Brown //   Model from:
208d8a22b9eSJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
209d8a22b9eSJed Brown //
210d8a22b9eSJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
211d8a22b9eSJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
212d8a22b9eSJed Brown //
213d8a22b9eSJed Brown // Where
214d8a22b9eSJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
215d8a22b9eSJed Brown //   h[i]      = 2 length(dxdX[i])
216d8a22b9eSJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
217d8a22b9eSJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
218d8a22b9eSJed Brown //   rho(A[i]) = spectral radius of the convective flux Jacobian i,
219d8a22b9eSJed Brown //               wave speed in direction i
220d8a22b9eSJed Brown // *****************************************************************************
2212b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], const CeedScalar dXdx[3][3], const CeedScalar u[3], const CeedScalar sound_speed,
2222b916ea7SJeremy L Thompson                                        const CeedScalar c_tau) {
223493642f1SJames Wright   for (CeedInt i = 0; i < 3; i++) {
224d8a22b9eSJed Brown     // length of element in direction i
2252b916ea7SJeremy L Thompson     CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]);
226d8a22b9eSJed Brown     // fastest wave in direction i
227d8a22b9eSJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
228d8a22b9eSJed Brown     Tau_x[i]                = c_tau * h / fastest_wave;
229d8a22b9eSJed Brown   }
230d8a22b9eSJed Brown }
231d8a22b9eSJed Brown 
232d8a22b9eSJed Brown // *****************************************************************************
233a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
234a515125bSLeila Ghaffari // *****************************************************************************
2352b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
236a515125bSLeila Ghaffari   // Inputs
237a515125bSLeila Ghaffari   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
238a515125bSLeila Ghaffari 
239a515125bSLeila Ghaffari   // Outputs
240a515125bSLeila Ghaffari   CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
241a515125bSLeila Ghaffari   const EulerContext context  = (EulerContext)ctx;
242a515125bSLeila Ghaffari 
243a515125bSLeila Ghaffari   // Quadrature Point Loop
244*3d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
245a515125bSLeila Ghaffari     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
246139613f2SLeila Ghaffari     CeedScalar       q[5] = {0.};
247a515125bSLeila Ghaffari 
248a515125bSLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
249a515125bSLeila Ghaffari 
2502b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
251a515125bSLeila Ghaffari   }  // End of Quadrature Point Loop
252a515125bSLeila Ghaffari 
253a515125bSLeila Ghaffari   // Return
254a515125bSLeila Ghaffari   return 0;
255a515125bSLeila Ghaffari }
256a515125bSLeila Ghaffari 
257a515125bSLeila Ghaffari // *****************************************************************************
258a515125bSLeila Ghaffari // This QFunction implements the following formulation of Euler equations
259a515125bSLeila Ghaffari //   with explicit time stepping method
260a515125bSLeila Ghaffari //
261a515125bSLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation
262a515125bSLeila Ghaffari //   form with state variables of density, momentum density, and total
263a515125bSLeila Ghaffari //   energy density.
264a515125bSLeila Ghaffari //
265a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
266a515125bSLeila Ghaffari //   rho - Mass Density
267a515125bSLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
268a515125bSLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
269a515125bSLeila Ghaffari //
270a515125bSLeila Ghaffari // Euler Equations:
271a515125bSLeila Ghaffari //   drho/dt + div( U )                   = 0
272a515125bSLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
273a515125bSLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
274a515125bSLeila Ghaffari //
275a515125bSLeila Ghaffari // Equation of State:
276a515125bSLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
277a515125bSLeila Ghaffari //
278a515125bSLeila Ghaffari // Constants:
279a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
280a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
281a515125bSLeila Ghaffari //   g               ,  Gravity
282a515125bSLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
283a515125bSLeila Ghaffari // *****************************************************************************
2842b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
285a515125bSLeila Ghaffari   // Inputs
286*3d65b166SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]      = (const CeedScalar(*)[CEED_Q_VLA])in[0];
287*3d65b166SJames Wright   const CeedScalar(*dq)[5][CEED_Q_VLA]  = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
288*3d65b166SJames Wright   const CeedScalar(*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
289*3d65b166SJames Wright 
290a515125bSLeila Ghaffari   // Outputs
291*3d65b166SJames Wright   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
292*3d65b166SJames Wright   CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
293a515125bSLeila Ghaffari 
294139613f2SLeila Ghaffari   EulerContext     context = (EulerContext)ctx;
295d8a22b9eSJed Brown   const CeedScalar c_tau   = context->c_tau;
296a515125bSLeila Ghaffari   const CeedScalar gamma   = 1.4;
297a515125bSLeila Ghaffari 
298a515125bSLeila Ghaffari   // Quadrature Point Loop
299*3d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
300a515125bSLeila Ghaffari     // Setup
301a515125bSLeila Ghaffari     // -- Interp in
302a515125bSLeila Ghaffari     const CeedScalar rho      = q[0][i];
3032b916ea7SJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
304a515125bSLeila Ghaffari     const CeedScalar E        = q[4][i];
3052b916ea7SJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
3062b916ea7SJeremy L Thompson     const CeedScalar dU[3][3] = {
3072b916ea7SJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
3082b916ea7SJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
3092b916ea7SJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
310139613f2SLeila Ghaffari     };
3112b916ea7SJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
312a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
313a515125bSLeila Ghaffari     const CeedScalar wdetJ = q_data[0][i];
314a515125bSLeila Ghaffari     // -- Interp-to-Grad q_data
315a515125bSLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
3162b916ea7SJeremy L Thompson     const CeedScalar dXdx[3][3] = {
3172b916ea7SJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
3182b916ea7SJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
3192b916ea7SJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
320a515125bSLeila Ghaffari     };
321139613f2SLeila Ghaffari     // dU/dx
322139613f2SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
323139613f2SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
324139613f2SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
325139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
326493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
327493642f1SJames Wright       for (CeedInt k = 0; k < 3; k++) {
328139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
329139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
330493642f1SJames Wright         for (CeedInt l = 0; l < 3; l++) {
331139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
332139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
333139613f2SLeila Ghaffari         }
334139613f2SLeila Ghaffari       }
335139613f2SLeila Ghaffari     }
336139613f2SLeila Ghaffari     // Pressure
3372b916ea7SJeremy L Thompson     const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic,
338139613f2SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
339a515125bSLeila Ghaffari 
340a515125bSLeila Ghaffari     // The Physics
341a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
342493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) {
343139613f2SLeila Ghaffari       v[j][i] = 0.;
3442b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
345a515125bSLeila Ghaffari     }
346a515125bSLeila Ghaffari 
347a515125bSLeila Ghaffari     // -- Density
348a515125bSLeila Ghaffari     // ---- u rho
3492b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) dv[j][0][i] += wdetJ * (rho * u[0] * dXdx[j][0] + rho * u[1] * dXdx[j][1] + rho * u[2] * dXdx[j][2]);
350a515125bSLeila Ghaffari     // -- Momentum
351a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
3522b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3532b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
3542b916ea7SJeremy L Thompson         dv[k][j + 1][i] += wdetJ * ((rho * u[j] * u[0] + (j == 0 ? P : 0.)) * dXdx[k][0] + (rho * u[j] * u[1] + (j == 1 ? P : 0.)) * dXdx[k][1] +
355139613f2SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
3562b916ea7SJeremy L Thompson       }
3572b916ea7SJeremy L Thompson     }
358a515125bSLeila Ghaffari     // -- Total Energy Density
359a515125bSLeila Ghaffari     // ---- (E + P) u
3602b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) dv[j][4][i] += wdetJ * (E + P) * (u[0] * dXdx[j][0] + u[1] * dXdx[j][1] + u[2] * dXdx[j][2]);
361139613f2SLeila Ghaffari 
362139613f2SLeila Ghaffari     // --Stabilization terms
363139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
364139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
365d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
366139613f2SLeila Ghaffari 
367139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
368139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
369493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
370139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
371139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
3722b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
373139613f2SLeila Ghaffari     }
374139613f2SLeila Ghaffari 
375139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
376139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
3772b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3782b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
3792b916ea7SJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
3802b916ea7SJeremy L Thompson       }
3812b916ea7SJeremy L Thompson     }
382139613f2SLeila Ghaffari 
383d8a22b9eSJed Brown     // Stabilization
384d8a22b9eSJed Brown     // -- Tau elements
385d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
386d8a22b9eSJed Brown     CeedScalar       Tau_x[3]    = {0.};
387d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
388139613f2SLeila Ghaffari 
389d8a22b9eSJed Brown     // -- Stabilization method: none or SU
390bb8a0c61SJames Wright     CeedScalar stab[5][3] = {{0.}};
391139613f2SLeila Ghaffari     switch (context->stabilization) {
392139613f2SLeila Ghaffari       case 0:  // Galerkin
393139613f2SLeila Ghaffari         break;
394139613f2SLeila Ghaffari       case 1:  // SU
3952b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
3962b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
3972b916ea7SJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
3982b916ea7SJeremy L Thompson           }
3992b916ea7SJeremy L Thompson         }
400139613f2SLeila Ghaffari 
4012b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
4022b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 3; k++) dv[k][j][i] -= wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]);
4032b916ea7SJeremy L Thompson         }
404139613f2SLeila Ghaffari         break;
405139613f2SLeila Ghaffari       case 2:  // SUPG is not implemented for explicit scheme
406139613f2SLeila Ghaffari         break;
407139613f2SLeila Ghaffari     }
408139613f2SLeila Ghaffari 
409a515125bSLeila Ghaffari   }  // End Quadrature Point Loop
410a515125bSLeila Ghaffari 
411a515125bSLeila Ghaffari   // Return
412a515125bSLeila Ghaffari   return 0;
413a515125bSLeila Ghaffari }
414a515125bSLeila Ghaffari 
415a515125bSLeila Ghaffari // *****************************************************************************
416a515125bSLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above)
417a515125bSLeila Ghaffari //   with implicit time stepping method
418a515125bSLeila Ghaffari //
419a515125bSLeila Ghaffari // *****************************************************************************
4202b916ea7SJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
421a515125bSLeila Ghaffari   // Inputs
422*3d65b166SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]      = (const CeedScalar(*)[CEED_Q_VLA])in[0];
423*3d65b166SJames Wright   const CeedScalar(*dq)[5][CEED_Q_VLA]  = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
424*3d65b166SJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA]  = (const CeedScalar(*)[CEED_Q_VLA])in[2];
425*3d65b166SJames Wright   const CeedScalar(*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
426*3d65b166SJames Wright 
427a515125bSLeila Ghaffari   // Outputs
428*3d65b166SJames Wright   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
429*3d65b166SJames Wright   CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
430a515125bSLeila Ghaffari 
431139613f2SLeila Ghaffari   EulerContext     context = (EulerContext)ctx;
432d8a22b9eSJed Brown   const CeedScalar c_tau   = context->c_tau;
433a515125bSLeila Ghaffari   const CeedScalar gamma   = 1.4;
434a515125bSLeila Ghaffari 
435a515125bSLeila Ghaffari   // Quadrature Point Loop
436*3d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
437a515125bSLeila Ghaffari     // Setup
438a515125bSLeila Ghaffari     // -- Interp in
439a515125bSLeila Ghaffari     const CeedScalar rho      = q[0][i];
4402b916ea7SJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
441a515125bSLeila Ghaffari     const CeedScalar E        = q[4][i];
4422b916ea7SJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
4432b916ea7SJeremy L Thompson     const CeedScalar dU[3][3] = {
4442b916ea7SJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
4452b916ea7SJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
4462b916ea7SJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
447139613f2SLeila Ghaffari     };
4482b916ea7SJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
449a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
450a515125bSLeila Ghaffari     const CeedScalar wdetJ = q_data[0][i];
451a515125bSLeila Ghaffari     // -- Interp-to-Grad q_data
452a515125bSLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
4532b916ea7SJeremy L Thompson     const CeedScalar dXdx[3][3] = {
4542b916ea7SJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
4552b916ea7SJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
4562b916ea7SJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
457a515125bSLeila Ghaffari     };
458139613f2SLeila Ghaffari     // dU/dx
459139613f2SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
460139613f2SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
461139613f2SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
462139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
463493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
464493642f1SJames Wright       for (CeedInt k = 0; k < 3; k++) {
465139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
466139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
467493642f1SJames Wright         for (CeedInt l = 0; l < 3; l++) {
468139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
469139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
470139613f2SLeila Ghaffari         }
471139613f2SLeila Ghaffari       }
472139613f2SLeila Ghaffari     }
4732b916ea7SJeremy L Thompson     const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic,
474139613f2SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
475a515125bSLeila Ghaffari 
476a515125bSLeila Ghaffari     // The Physics
477a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
478493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) {
479139613f2SLeila Ghaffari       v[j][i] = 0.;
4802b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
481a515125bSLeila Ghaffari     }
482a515125bSLeila Ghaffari     //-----mass matrix
4832b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i];
484a515125bSLeila Ghaffari 
485a515125bSLeila Ghaffari     // -- Density
486a515125bSLeila Ghaffari     // ---- u rho
4872b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) dv[j][0][i] -= wdetJ * (rho * u[0] * dXdx[j][0] + rho * u[1] * dXdx[j][1] + rho * u[2] * dXdx[j][2]);
488a515125bSLeila Ghaffari     // -- Momentum
489a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
4902b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
4912b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
4922b916ea7SJeremy L Thompson         dv[k][j + 1][i] -= wdetJ * ((rho * u[j] * u[0] + (j == 0 ? P : 0.)) * dXdx[k][0] + (rho * u[j] * u[1] + (j == 1 ? P : 0.)) * dXdx[k][1] +
493139613f2SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
4942b916ea7SJeremy L Thompson       }
4952b916ea7SJeremy L Thompson     }
496a515125bSLeila Ghaffari     // -- Total Energy Density
497a515125bSLeila Ghaffari     // ---- (E + P) u
4982b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) dv[j][4][i] -= wdetJ * (E + P) * (u[0] * dXdx[j][0] + u[1] * dXdx[j][1] + u[2] * dXdx[j][2]);
499139613f2SLeila Ghaffari 
500139613f2SLeila Ghaffari     // -- Stabilization terms
501139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
502139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
503d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
504139613f2SLeila Ghaffari 
505139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
506139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
507493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
508139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
509139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
5102b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
511139613f2SLeila Ghaffari     }
512139613f2SLeila Ghaffari 
513139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
514139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
5152b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
5162b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
5172b916ea7SJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
5182b916ea7SJeremy L Thompson       }
5192b916ea7SJeremy L Thompson     }
520139613f2SLeila Ghaffari 
521139613f2SLeila Ghaffari     // ---- Strong residual
522139613f2SLeila Ghaffari     CeedScalar strong_res[5];
5232b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j];
524139613f2SLeila Ghaffari 
525d8a22b9eSJed Brown     // Stabilization
526d8a22b9eSJed Brown     // -- Tau elements
527d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
528d8a22b9eSJed Brown     CeedScalar       Tau_x[3]    = {0.};
529d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
530139613f2SLeila Ghaffari 
531d8a22b9eSJed Brown     // -- Stabilization method: none, SU, or SUPG
532bb8a0c61SJames Wright     CeedScalar stab[5][3] = {{0.}};
533139613f2SLeila Ghaffari     switch (context->stabilization) {
534139613f2SLeila Ghaffari       case 0:  // Galerkin
535139613f2SLeila Ghaffari         break;
536139613f2SLeila Ghaffari       case 1:  // SU
5372b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
5382b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
5392b916ea7SJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
5402b916ea7SJeremy L Thompson           }
5412b916ea7SJeremy L Thompson         }
542139613f2SLeila Ghaffari 
5432b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
5442b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 3; k++) dv[k][j][i] += wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]);
5452b916ea7SJeremy L Thompson         }
546139613f2SLeila Ghaffari         break;
547139613f2SLeila Ghaffari       case 2:  // SUPG
5482b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
5492b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
5502b916ea7SJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
5512b916ea7SJeremy L Thompson           }
5522b916ea7SJeremy L Thompson         }
553139613f2SLeila Ghaffari 
5542b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
5552b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 3; k++) dv[k][j][i] += wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]);
5562b916ea7SJeremy L Thompson         }
557139613f2SLeila Ghaffari         break;
558139613f2SLeila Ghaffari     }
559a515125bSLeila Ghaffari   }  // End Quadrature Point Loop
560a515125bSLeila Ghaffari 
561a515125bSLeila Ghaffari   // Return
562a515125bSLeila Ghaffari   return 0;
563a515125bSLeila Ghaffari }
564a515125bSLeila Ghaffari // *****************************************************************************
565002797a3SLeila Ghaffari // This QFunction sets the inflow boundary conditions for
566002797a3SLeila Ghaffari //   the traveling vortex problem.
567a515125bSLeila Ghaffari //
568a515125bSLeila Ghaffari //  Prescribed T_inlet and P_inlet are converted to conservative variables
569a515125bSLeila Ghaffari //      and applied weakly.
570a515125bSLeila Ghaffari //
571a515125bSLeila Ghaffari // *****************************************************************************
5722b916ea7SJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
573a515125bSLeila Ghaffari   // Inputs
574dd64951cSJames Wright   const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
575a515125bSLeila Ghaffari   // Outputs
576a515125bSLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
577a515125bSLeila Ghaffari   EulerContext     context       = (EulerContext)ctx;
578a515125bSLeila Ghaffari   const int        euler_test    = context->euler_test;
579a515125bSLeila Ghaffari   const bool       implicit      = context->implicit;
580a515125bSLeila Ghaffari   CeedScalar      *mean_velocity = context->mean_velocity;
581a515125bSLeila Ghaffari   const CeedScalar cv            = 2.5;
582a515125bSLeila Ghaffari   const CeedScalar R             = 1.;
583a515125bSLeila Ghaffari   CeedScalar       T_inlet;
584a515125bSLeila Ghaffari   CeedScalar       P_inlet;
585a515125bSLeila Ghaffari 
586a515125bSLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
5872b916ea7SJeremy L Thompson   if (euler_test == 1 || euler_test == 3) {
588a515125bSLeila Ghaffari     for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.;
5892b916ea7SJeremy L Thompson   }
590a515125bSLeila Ghaffari 
591a515125bSLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
592a515125bSLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
593a515125bSLeila Ghaffari   else T_inlet = P_inlet = 1.;
594a515125bSLeila Ghaffari 
595a515125bSLeila Ghaffari   // Quadrature Point Loop
596*3d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
597a515125bSLeila Ghaffari     // Setup
598a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
599a515125bSLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
600a515125bSLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
601a515125bSLeila Ghaffari     // We can effect this by swapping the sign on this weight
602a515125bSLeila Ghaffari     const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i];
603002797a3SLeila Ghaffari     // ---- Normal vect
6042b916ea7SJeremy L Thompson     const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
605a515125bSLeila Ghaffari 
606a515125bSLeila Ghaffari     // face_normal = Normal vector of the face
6072b916ea7SJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
608a515125bSLeila Ghaffari     // The Physics
609a515125bSLeila Ghaffari     // Zero v so all future terms can safely sum into it
610493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
611a515125bSLeila Ghaffari 
612a515125bSLeila Ghaffari     // Implementing in/outflow BCs
613002797a3SLeila Ghaffari     if (face_normal > 0) {
614a515125bSLeila Ghaffari     } else {  // inflow
615a515125bSLeila Ghaffari       const CeedScalar rho_inlet       = P_inlet / (R * T_inlet);
6162b916ea7SJeremy L Thompson       const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.;
617a515125bSLeila Ghaffari       // incoming total energy
618a515125bSLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
619a515125bSLeila Ghaffari 
620a515125bSLeila Ghaffari       // The Physics
621a515125bSLeila Ghaffari       // -- Density
622a515125bSLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
623a515125bSLeila Ghaffari 
624a515125bSLeila Ghaffari       // -- Momentum
6252b916ea7SJeremy L Thompson       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho_inlet * face_normal * mean_velocity[j] + norm[j] * P_inlet);
626a515125bSLeila Ghaffari 
627a515125bSLeila Ghaffari       // -- Total Energy Density
628a515125bSLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
629a515125bSLeila Ghaffari     }
630a515125bSLeila Ghaffari 
631a515125bSLeila Ghaffari   }  // End Quadrature Point Loop
632a515125bSLeila Ghaffari   return 0;
633a515125bSLeila Ghaffari }
634a515125bSLeila Ghaffari 
635a515125bSLeila Ghaffari // *****************************************************************************
63668ef3d20SLeila Ghaffari // This QFunction sets the outflow boundary conditions for
63768ef3d20SLeila Ghaffari //   the Euler solver.
63868ef3d20SLeila Ghaffari //
63968ef3d20SLeila Ghaffari //  Outflow BCs:
64068ef3d20SLeila Ghaffari //    The validity of the weak form of the governing equations is
64168ef3d20SLeila Ghaffari //      extended to the outflow.
64268ef3d20SLeila Ghaffari //
64368ef3d20SLeila Ghaffari // *****************************************************************************
6442b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
64568ef3d20SLeila Ghaffari   // Inputs
646*3d65b166SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA])in[0];
647*3d65b166SJames Wright   const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
648*3d65b166SJames Wright 
64968ef3d20SLeila Ghaffari   // Outputs
65068ef3d20SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
65168ef3d20SLeila Ghaffari   EulerContext context       = (EulerContext)ctx;
65268ef3d20SLeila Ghaffari   const bool   implicit      = context->implicit;
65368ef3d20SLeila Ghaffari   CeedScalar  *mean_velocity = context->mean_velocity;
65468ef3d20SLeila Ghaffari 
65568ef3d20SLeila Ghaffari   const CeedScalar gamma = 1.4;
65668ef3d20SLeila Ghaffari 
65768ef3d20SLeila Ghaffari   // Quadrature Point Loop
658*3d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
65968ef3d20SLeila Ghaffari     // Setup
66068ef3d20SLeila Ghaffari     // -- Interp in
66168ef3d20SLeila Ghaffari     const CeedScalar rho  = q[0][i];
6622b916ea7SJeremy L Thompson     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
66368ef3d20SLeila Ghaffari     const CeedScalar E    = q[4][i];
66468ef3d20SLeila Ghaffari 
66568ef3d20SLeila Ghaffari     // -- Interp-to-Interp q_data
66668ef3d20SLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
66768ef3d20SLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
66868ef3d20SLeila Ghaffari     // We can effect this by swapping the sign on this weight
66968ef3d20SLeila Ghaffari     const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i];
67068ef3d20SLeila Ghaffari     // ---- Normal vectors
6712b916ea7SJeremy L Thompson     const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
67268ef3d20SLeila Ghaffari 
67368ef3d20SLeila Ghaffari     // face_normal = Normal vector of the face
6742b916ea7SJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
67568ef3d20SLeila Ghaffari     // The Physics
67668ef3d20SLeila Ghaffari     // Zero v so all future terms can safely sum into it
677493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0;
67868ef3d20SLeila Ghaffari 
67968ef3d20SLeila Ghaffari     // Implementing in/outflow BCs
68068ef3d20SLeila Ghaffari     if (face_normal > 0) {  // outflow
68168ef3d20SLeila Ghaffari       const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.;
68268ef3d20SLeila Ghaffari       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.);              // pressure
6832b916ea7SJeremy L Thompson       const CeedScalar u_normal  = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2];  // Normal velocity
68468ef3d20SLeila Ghaffari       // The Physics
68568ef3d20SLeila Ghaffari       // -- Density
68668ef3d20SLeila Ghaffari       v[0][i] -= wdetJb * rho * u_normal;
68768ef3d20SLeila Ghaffari 
68868ef3d20SLeila Ghaffari       // -- Momentum
6892b916ea7SJeremy L Thompson       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P);
69068ef3d20SLeila Ghaffari 
69168ef3d20SLeila Ghaffari       // -- Total Energy Density
69268ef3d20SLeila Ghaffari       v[4][i] -= wdetJb * u_normal * (E + P);
69368ef3d20SLeila Ghaffari     }
69468ef3d20SLeila Ghaffari   }  // End Quadrature Point Loop
69568ef3d20SLeila Ghaffari   return 0;
69668ef3d20SLeila Ghaffari }
69768ef3d20SLeila Ghaffari 
69868ef3d20SLeila Ghaffari // *****************************************************************************
699a515125bSLeila Ghaffari 
700a515125bSLeila Ghaffari #endif  // eulervortex_h
701