xref: /honee/qfunctions/eulervortex.h (revision ade4951188ee005b66ffd0a9de0bd2ad8f48d7f6)
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:
1304e40bb6SJeremy L Thompson //   On the Order of Accuracy and Numerical Performance of Two Classes of Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011).
14a515125bSLeila Ghaffari 
15a515125bSLeila Ghaffari #ifndef eulervortex_h
16a515125bSLeila Ghaffari #define eulervortex_h
17a515125bSLeila Ghaffari 
183a8779fbSJames Wright #include <ceed.h>
19d0cce58aSJeremy L Thompson #include <math.h>
202b916ea7SJeremy L Thompson 
21704b8bbeSJames Wright #include "utils.h"
22a515125bSLeila Ghaffari 
23a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext;
24a515125bSLeila Ghaffari struct EulerContext_ {
25a515125bSLeila Ghaffari   CeedScalar center[3];
26a515125bSLeila Ghaffari   CeedScalar curr_time;
27a515125bSLeila Ghaffari   CeedScalar vortex_strength;
28d8a22b9eSJed Brown   CeedScalar c_tau;
29a515125bSLeila Ghaffari   CeedScalar mean_velocity[3];
30a515125bSLeila Ghaffari   bool       implicit;
31139613f2SLeila Ghaffari   int        euler_test;
32139613f2SLeila Ghaffari   int        stabilization;  // See StabilizationType: 0=none, 1=SU, 2=SUPG
33a515125bSLeila Ghaffari };
34a515125bSLeila Ghaffari 
35a515125bSLeila Ghaffari // *****************************************************************************
36a515125bSLeila Ghaffari // This function sets the initial conditions
37a515125bSLeila Ghaffari //
38a515125bSLeila Ghaffari //   Temperature:
39a515125bSLeila Ghaffari //     T   = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
40a515125bSLeila Ghaffari //   Density:
41a515125bSLeila Ghaffari //     rho = (T/S_vortex)^(1 / (gamma - 1))
42a515125bSLeila Ghaffari //   Pressure:
43a515125bSLeila Ghaffari //     P   = rho * T
44a515125bSLeila Ghaffari //   Velocity:
45a515125bSLeila Ghaffari //     ui  = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
46a515125bSLeila Ghaffari //     r   = sqrt( (x - xc)**2 + (y - yc)**2 )
47a515125bSLeila Ghaffari //   Velocity/Momentum Density:
48a515125bSLeila Ghaffari //     Ui  = rho ui
49a515125bSLeila Ghaffari //   Total Energy:
50a515125bSLeila Ghaffari //     E   = P / (gamma - 1) + rho (u u)/2
51a515125bSLeila Ghaffari //
52a515125bSLeila Ghaffari // Constants:
53a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
54a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
55a515125bSLeila Ghaffari //   vortex_strength ,  Strength of vortex
56a515125bSLeila Ghaffari //   center          ,  Location of bubble center
57a515125bSLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
58a515125bSLeila Ghaffari //
59a515125bSLeila Ghaffari // *****************************************************************************
60a515125bSLeila Ghaffari 
61a515125bSLeila Ghaffari // *****************************************************************************
6204e40bb6SJeremy L Thompson // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling
6304e40bb6SJeremy L Thompson // vortex
64a515125bSLeila Ghaffari // *****************************************************************************
652b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
66a515125bSLeila Ghaffari   // Context
67a515125bSLeila Ghaffari   const EulerContext context         = (EulerContext)ctx;
68a515125bSLeila Ghaffari   const CeedScalar   vortex_strength = context->vortex_strength;
69a515125bSLeila Ghaffari   const CeedScalar  *center          = context->center;  // Center of the domain
70a515125bSLeila Ghaffari   const CeedScalar  *mean_velocity   = context->mean_velocity;
71a515125bSLeila Ghaffari 
72a515125bSLeila Ghaffari   // Setup
73a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
74a515125bSLeila Ghaffari   const CeedScalar cv    = 2.5;
75a515125bSLeila Ghaffari   const CeedScalar R     = 1.;
76a515125bSLeila Ghaffari   const CeedScalar x = X[0], y = X[1];  // Coordinates
77a515125bSLeila Ghaffari   // Vortex center
78a515125bSLeila Ghaffari   const CeedScalar xc = center[0] + mean_velocity[0] * time;
79a515125bSLeila Ghaffari   const CeedScalar yc = center[1] + mean_velocity[1] * time;
80a515125bSLeila Ghaffari 
81a515125bSLeila Ghaffari   const CeedScalar x0       = x - xc;
82a515125bSLeila Ghaffari   const CeedScalar y0       = y - yc;
83a515125bSLeila Ghaffari   const CeedScalar r        = sqrt(x0 * x0 + y0 * y0);
84a515125bSLeila Ghaffari   const CeedScalar C        = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI);
852b916ea7SJeremy L Thompson   const CeedScalar delta_T  = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI);
86a515125bSLeila Ghaffari   const CeedScalar S_vortex = 1;  // no perturbation in the entropy P / rho^gamma
872b916ea7SJeremy L Thompson   const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI);
88a515125bSLeila Ghaffari   CeedScalar       rho, P, T, E, u[3] = {0.};
89a515125bSLeila Ghaffari 
90a515125bSLeila Ghaffari   // Initial Conditions
91a515125bSLeila Ghaffari   switch (context->euler_test) {
92a515125bSLeila Ghaffari     case 0:  // Traveling vortex
93a515125bSLeila Ghaffari       T = 1 + delta_T;
94a515125bSLeila Ghaffari       // P = rho * T
95a515125bSLeila Ghaffari       // P = S * rho^gamma
96a515125bSLeila Ghaffari       // Solve for rho, then substitute for P
97139613f2SLeila Ghaffari       rho  = pow(T / S_vortex, 1 / (gamma - 1.));
98a515125bSLeila Ghaffari       P    = rho * T;
99a515125bSLeila Ghaffari       u[0] = mean_velocity[0] - C * y0;
100a515125bSLeila Ghaffari       u[1] = mean_velocity[1] + C * x0;
101a515125bSLeila Ghaffari 
102a515125bSLeila Ghaffari       // Assign exact solution
103a515125bSLeila Ghaffari       q[0] = rho;
104a515125bSLeila Ghaffari       q[1] = rho * u[0];
105a515125bSLeila Ghaffari       q[2] = rho * u[1];
106a515125bSLeila Ghaffari       q[3] = rho * u[2];
107a515125bSLeila Ghaffari       q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.;
108a515125bSLeila Ghaffari       break;
109a515125bSLeila Ghaffari     case 1:  // Constant zero velocity, density constant, total energy constant
110a515125bSLeila Ghaffari       rho = 1.;
111a515125bSLeila Ghaffari       E   = 2.;
112a515125bSLeila Ghaffari 
113a515125bSLeila Ghaffari       // Assign exact solution
114a515125bSLeila Ghaffari       q[0] = rho;
115a515125bSLeila Ghaffari       q[1] = rho * u[0];
116a515125bSLeila Ghaffari       q[2] = rho * u[1];
117a515125bSLeila Ghaffari       q[3] = rho * u[2];
118a515125bSLeila Ghaffari       q[4] = E;
119a515125bSLeila Ghaffari       break;
120a515125bSLeila Ghaffari     case 2:  // Constant nonzero velocity, density constant, total energy constant
121a515125bSLeila Ghaffari       rho  = 1.;
122a515125bSLeila Ghaffari       E    = 2.;
123a515125bSLeila Ghaffari       u[0] = mean_velocity[0];
124a515125bSLeila Ghaffari       u[1] = mean_velocity[1];
125a515125bSLeila Ghaffari 
126a515125bSLeila Ghaffari       // Assign exact solution
127a515125bSLeila Ghaffari       q[0] = rho;
128a515125bSLeila Ghaffari       q[1] = rho * u[0];
129a515125bSLeila Ghaffari       q[2] = rho * u[1];
130a515125bSLeila Ghaffari       q[3] = rho * u[2];
131a515125bSLeila Ghaffari       q[4] = E;
132a515125bSLeila Ghaffari       break;
13304e40bb6SJeremy L Thompson     case 3:  // Velocity zero, pressure constant (so density and internal energy will be non-constant), but the velocity should stay zero and the
13404e40bb6SJeremy L Thompson              // bubble won't diffuse
135a515125bSLeila Ghaffari       // (for Euler, where there is no thermal conductivity)
136a515125bSLeila Ghaffari       P   = 1.;
137a515125bSLeila Ghaffari       T   = 1. - S_bubble * exp(1. - r * r);
138a515125bSLeila Ghaffari       rho = P / (R * T);
139a515125bSLeila Ghaffari 
140a515125bSLeila Ghaffari       // Assign exact solution
141a515125bSLeila Ghaffari       q[0] = rho;
142a515125bSLeila Ghaffari       q[1] = rho * u[0];
143a515125bSLeila Ghaffari       q[2] = rho * u[1];
144a515125bSLeila Ghaffari       q[3] = rho * u[2];
145a515125bSLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
146a515125bSLeila Ghaffari       break;
14704e40bb6SJeremy L Thompson     case 4:  // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant),
14804e40bb6SJeremy L Thompson       // It should be transported across the domain, but velocity stays constant
149a515125bSLeila Ghaffari       P    = 1.;
150a515125bSLeila Ghaffari       T    = 1. - S_bubble * exp(1. - r * r);
151a515125bSLeila Ghaffari       rho  = P / (R * T);
152a515125bSLeila Ghaffari       u[0] = mean_velocity[0];
153a515125bSLeila Ghaffari       u[1] = mean_velocity[1];
154a515125bSLeila Ghaffari 
155a515125bSLeila Ghaffari       // Assign exact solution
156a515125bSLeila Ghaffari       q[0] = rho;
157a515125bSLeila Ghaffari       q[1] = rho * u[0];
158a515125bSLeila Ghaffari       q[2] = rho * u[1];
159a515125bSLeila Ghaffari       q[3] = rho * u[2];
160a515125bSLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
161a515125bSLeila Ghaffari       break;
1620df2634dSLeila Ghaffari     case 5:  // non-smooth thermal bubble - cylinder
1630df2634dSLeila Ghaffari       P    = 1.;
1640df2634dSLeila Ghaffari       T    = 1. - (r < 1. ? S_bubble : 0.);
1650df2634dSLeila Ghaffari       rho  = P / (R * T);
1660df2634dSLeila Ghaffari       u[0] = mean_velocity[0];
1670df2634dSLeila Ghaffari       u[1] = mean_velocity[1];
1680df2634dSLeila Ghaffari 
1690df2634dSLeila Ghaffari       // Assign exact solution
1700df2634dSLeila Ghaffari       q[0] = rho;
1710df2634dSLeila Ghaffari       q[1] = rho * u[0];
1720df2634dSLeila Ghaffari       q[2] = rho * u[1];
1730df2634dSLeila Ghaffari       q[3] = rho * u[2];
1740df2634dSLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
1750df2634dSLeila Ghaffari       break;
176a515125bSLeila Ghaffari   }
177a515125bSLeila Ghaffari   // Return
178a515125bSLeila Ghaffari   return 0;
179a515125bSLeila Ghaffari }
180a515125bSLeila Ghaffari 
181a515125bSLeila Ghaffari // *****************************************************************************
182139613f2SLeila Ghaffari // Helper function for computing flux Jacobian
183139613f2SLeila Ghaffari // *****************************************************************************
1842b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
185139613f2SLeila Ghaffari                                                         const CeedScalar gamma) {
186139613f2SLeila Ghaffari   CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2];  // Velocity square
187139613f2SLeila Ghaffari   for (CeedInt i = 0; i < 3; i++) {                           // Jacobian matrices for 3 directions
188139613f2SLeila Ghaffari     for (CeedInt j = 0; j < 3; j++) {                         // Rows of each Jacobian matrix
189139613f2SLeila Ghaffari       dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j];
190139613f2SLeila Ghaffari       for (CeedInt k = 0; k < 3; k++) {  // Columns of each Jacobian matrix
191139613f2SLeila Ghaffari         dF[i][0][k + 1]     = ((i == k) ? 1. : 0.);
1922b916ea7SJeremy 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.);
1932b916ea7SJeremy L Thompson         dF[i][4][k + 1]     = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k];
194139613f2SLeila Ghaffari       }
195139613f2SLeila Ghaffari       dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.);
196139613f2SLeila Ghaffari     }
197139613f2SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho);
198139613f2SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
199139613f2SLeila Ghaffari   }
200139613f2SLeila Ghaffari }
201139613f2SLeila Ghaffari 
202139613f2SLeila Ghaffari // *****************************************************************************
203d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant)
204d8a22b9eSJed Brown //   Model from:
205d8a22b9eSJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
206d8a22b9eSJed Brown //
207d8a22b9eSJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
208d8a22b9eSJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
209d8a22b9eSJed Brown //
210d8a22b9eSJed Brown // Where
211d8a22b9eSJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
212d8a22b9eSJed Brown //   h[i]      = 2 length(dxdX[i])
213d8a22b9eSJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
214d8a22b9eSJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
21504e40bb6SJeremy L Thompson //   rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i
216d8a22b9eSJed Brown // *****************************************************************************
2172b916ea7SJeremy 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,
2182b916ea7SJeremy L Thompson                                        const CeedScalar c_tau) {
219493642f1SJames Wright   for (CeedInt i = 0; i < 3; i++) {
220d8a22b9eSJed Brown     // length of element in direction i
2212b916ea7SJeremy L Thompson     CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]);
222d8a22b9eSJed Brown     // fastest wave in direction i
223d8a22b9eSJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
224d8a22b9eSJed Brown     Tau_x[i]                = c_tau * h / fastest_wave;
225d8a22b9eSJed Brown   }
226d8a22b9eSJed Brown }
227d8a22b9eSJed Brown 
228d8a22b9eSJed Brown // *****************************************************************************
229a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
230a515125bSLeila Ghaffari // *****************************************************************************
2312b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
232a515125bSLeila Ghaffari   // Inputs
233a515125bSLeila Ghaffari   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
234a515125bSLeila Ghaffari 
235a515125bSLeila Ghaffari   // Outputs
236a515125bSLeila Ghaffari   CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
237a515125bSLeila Ghaffari   const EulerContext context  = (EulerContext)ctx;
238a515125bSLeila Ghaffari 
239a515125bSLeila Ghaffari   // Quadrature Point Loop
2403d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
241a515125bSLeila Ghaffari     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
242139613f2SLeila Ghaffari     CeedScalar       q[5] = {0.};
243a515125bSLeila Ghaffari 
244a515125bSLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
245a515125bSLeila Ghaffari 
2462b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
247a515125bSLeila Ghaffari   }  // End of Quadrature Point Loop
248a515125bSLeila Ghaffari 
249a515125bSLeila Ghaffari   // Return
250a515125bSLeila Ghaffari   return 0;
251a515125bSLeila Ghaffari }
252a515125bSLeila Ghaffari 
253a515125bSLeila Ghaffari // *****************************************************************************
25404e40bb6SJeremy L Thompson // This QFunction implements the following formulation of Euler equations with explicit time stepping method
255a515125bSLeila Ghaffari //
25604e40bb6SJeremy L Thompson // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density.
257a515125bSLeila Ghaffari //
258a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
259a515125bSLeila Ghaffari //   rho - Mass Density
260a515125bSLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
261a515125bSLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
262a515125bSLeila Ghaffari //
263a515125bSLeila Ghaffari // Euler Equations:
264a515125bSLeila Ghaffari //   drho/dt + div( U )                   = 0
265a515125bSLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
266a515125bSLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
267a515125bSLeila Ghaffari //
268a515125bSLeila Ghaffari // Equation of State:
269a515125bSLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
270a515125bSLeila Ghaffari //
271a515125bSLeila Ghaffari // Constants:
272a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
273a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
274a515125bSLeila Ghaffari //   g               ,  Gravity
275a515125bSLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
276a515125bSLeila Ghaffari // *****************************************************************************
2772b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
278a515125bSLeila Ghaffari   // Inputs
2793d65b166SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
2803d65b166SJames Wright   const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
281*ade49511SJames Wright   const CeedScalar(*q_data)            = in[2];
2823d65b166SJames Wright 
283a515125bSLeila Ghaffari   // Outputs
2843d65b166SJames Wright   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
2853d65b166SJames Wright   CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
286a515125bSLeila Ghaffari 
287139613f2SLeila Ghaffari   EulerContext     context = (EulerContext)ctx;
288d8a22b9eSJed Brown   const CeedScalar c_tau   = context->c_tau;
289a515125bSLeila Ghaffari   const CeedScalar gamma   = 1.4;
290a515125bSLeila Ghaffari 
291a515125bSLeila Ghaffari   // Quadrature Point Loop
2923d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
293a515125bSLeila Ghaffari     // Setup
294a515125bSLeila Ghaffari     // -- Interp in
295a515125bSLeila Ghaffari     const CeedScalar rho      = q[0][i];
2962b916ea7SJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
297a515125bSLeila Ghaffari     const CeedScalar E        = q[4][i];
2982b916ea7SJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
2992b916ea7SJeremy L Thompson     const CeedScalar dU[3][3] = {
3002b916ea7SJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
3012b916ea7SJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
3022b916ea7SJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
303139613f2SLeila Ghaffari     };
3042b916ea7SJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
305*ade49511SJames Wright     CeedScalar       wdetJ, dXdx[3][3];
306*ade49511SJames Wright     QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
307139613f2SLeila Ghaffari     // dU/dx
308139613f2SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
309139613f2SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
310139613f2SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
311139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
312493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
313493642f1SJames Wright       for (CeedInt k = 0; k < 3; k++) {
314139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
315139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
316493642f1SJames Wright         for (CeedInt l = 0; l < 3; l++) {
317139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
318139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
319139613f2SLeila Ghaffari         }
320139613f2SLeila Ghaffari       }
321139613f2SLeila Ghaffari     }
322139613f2SLeila Ghaffari     // Pressure
3232b916ea7SJeremy 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,
324139613f2SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
325a515125bSLeila Ghaffari 
326a515125bSLeila Ghaffari     // The Physics
327a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
328493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) {
329139613f2SLeila Ghaffari       v[j][i] = 0.;
3302b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
331a515125bSLeila Ghaffari     }
332a515125bSLeila Ghaffari 
333a515125bSLeila Ghaffari     // -- Density
334a515125bSLeila Ghaffari     // ---- u rho
3352b916ea7SJeremy 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]);
336a515125bSLeila Ghaffari     // -- Momentum
337a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
3382b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3392b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
3402b916ea7SJeremy 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] +
341139613f2SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
3422b916ea7SJeremy L Thompson       }
3432b916ea7SJeremy L Thompson     }
344a515125bSLeila Ghaffari     // -- Total Energy Density
345a515125bSLeila Ghaffari     // ---- (E + P) u
3462b916ea7SJeremy 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]);
347139613f2SLeila Ghaffari 
348139613f2SLeila Ghaffari     // --Stabilization terms
349139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
350139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
351d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
352139613f2SLeila Ghaffari 
353139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
354139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
355493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
356139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
357139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
3582b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
359139613f2SLeila Ghaffari     }
360139613f2SLeila Ghaffari 
361139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
362139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
3632b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3642b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
3652b916ea7SJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
3662b916ea7SJeremy L Thompson       }
3672b916ea7SJeremy L Thompson     }
368139613f2SLeila Ghaffari 
369d8a22b9eSJed Brown     // Stabilization
370d8a22b9eSJed Brown     // -- Tau elements
371d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
372d8a22b9eSJed Brown     CeedScalar       Tau_x[3]    = {0.};
373d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
374139613f2SLeila Ghaffari 
375d8a22b9eSJed Brown     // -- Stabilization method: none or SU
376bb8a0c61SJames Wright     CeedScalar stab[5][3] = {{0.}};
377139613f2SLeila Ghaffari     switch (context->stabilization) {
378139613f2SLeila Ghaffari       case 0:  // Galerkin
379139613f2SLeila Ghaffari         break;
380139613f2SLeila Ghaffari       case 1:  // SU
3812b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
3822b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
3832b916ea7SJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
3842b916ea7SJeremy L Thompson           }
3852b916ea7SJeremy L Thompson         }
386139613f2SLeila Ghaffari 
3872b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
3882b916ea7SJeremy 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]);
3892b916ea7SJeremy L Thompson         }
390139613f2SLeila Ghaffari         break;
391139613f2SLeila Ghaffari       case 2:  // SUPG is not implemented for explicit scheme
392139613f2SLeila Ghaffari         break;
393139613f2SLeila Ghaffari     }
394139613f2SLeila Ghaffari 
395a515125bSLeila Ghaffari   }  // End Quadrature Point Loop
396a515125bSLeila Ghaffari 
397a515125bSLeila Ghaffari   // Return
398a515125bSLeila Ghaffari   return 0;
399a515125bSLeila Ghaffari }
400a515125bSLeila Ghaffari 
401a515125bSLeila Ghaffari // *****************************************************************************
40204e40bb6SJeremy L Thompson // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method
403a515125bSLeila Ghaffari // *****************************************************************************
4042b916ea7SJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
405a515125bSLeila Ghaffari   // Inputs
4063d65b166SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
4073d65b166SJames Wright   const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
4083d65b166SJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
409*ade49511SJames Wright   const CeedScalar(*q_data)            = in[3];
4103d65b166SJames Wright 
411a515125bSLeila Ghaffari   // Outputs
4123d65b166SJames Wright   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
4133d65b166SJames Wright   CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
414a515125bSLeila Ghaffari 
415139613f2SLeila Ghaffari   EulerContext     context = (EulerContext)ctx;
416d8a22b9eSJed Brown   const CeedScalar c_tau   = context->c_tau;
417a515125bSLeila Ghaffari   const CeedScalar gamma   = 1.4;
418a515125bSLeila Ghaffari 
419a515125bSLeila Ghaffari   // Quadrature Point Loop
4203d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
421a515125bSLeila Ghaffari     // Setup
422a515125bSLeila Ghaffari     // -- Interp in
423a515125bSLeila Ghaffari     const CeedScalar rho      = q[0][i];
4242b916ea7SJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
425a515125bSLeila Ghaffari     const CeedScalar E        = q[4][i];
4262b916ea7SJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
4272b916ea7SJeremy L Thompson     const CeedScalar dU[3][3] = {
4282b916ea7SJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
4292b916ea7SJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
4302b916ea7SJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
431139613f2SLeila Ghaffari     };
4322b916ea7SJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
433*ade49511SJames Wright     CeedScalar       wdetJ, dXdx[3][3];
434*ade49511SJames Wright     QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
435139613f2SLeila Ghaffari     // dU/dx
436139613f2SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
437139613f2SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
438139613f2SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
439139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
440493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
441493642f1SJames Wright       for (CeedInt k = 0; k < 3; k++) {
442139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
443139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
444493642f1SJames Wright         for (CeedInt l = 0; l < 3; l++) {
445139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
446139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
447139613f2SLeila Ghaffari         }
448139613f2SLeila Ghaffari       }
449139613f2SLeila Ghaffari     }
4502b916ea7SJeremy 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,
451139613f2SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
452a515125bSLeila Ghaffari 
453a515125bSLeila Ghaffari     // The Physics
454a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
455493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) {
456139613f2SLeila Ghaffari       v[j][i] = 0.;
4572b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
458a515125bSLeila Ghaffari     }
459a515125bSLeila Ghaffari     //-----mass matrix
4602b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i];
461a515125bSLeila Ghaffari 
462a515125bSLeila Ghaffari     // -- Density
463a515125bSLeila Ghaffari     // ---- u rho
4642b916ea7SJeremy 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]);
465a515125bSLeila Ghaffari     // -- Momentum
466a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
4672b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
4682b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
4692b916ea7SJeremy 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] +
470139613f2SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
4712b916ea7SJeremy L Thompson       }
4722b916ea7SJeremy L Thompson     }
473a515125bSLeila Ghaffari     // -- Total Energy Density
474a515125bSLeila Ghaffari     // ---- (E + P) u
4752b916ea7SJeremy 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]);
476139613f2SLeila Ghaffari 
477139613f2SLeila Ghaffari     // -- Stabilization terms
478139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
479139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
480d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
481139613f2SLeila Ghaffari 
482139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
483139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
484493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
485139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
486139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
4872b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
488139613f2SLeila Ghaffari     }
489139613f2SLeila Ghaffari 
490139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
491139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
4922b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
4932b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
4942b916ea7SJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
4952b916ea7SJeremy L Thompson       }
4962b916ea7SJeremy L Thompson     }
497139613f2SLeila Ghaffari 
498139613f2SLeila Ghaffari     // ---- Strong residual
499139613f2SLeila Ghaffari     CeedScalar strong_res[5];
5002b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j];
501139613f2SLeila Ghaffari 
502d8a22b9eSJed Brown     // Stabilization
503d8a22b9eSJed Brown     // -- Tau elements
504d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
505d8a22b9eSJed Brown     CeedScalar       Tau_x[3]    = {0.};
506d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
507139613f2SLeila Ghaffari 
508d8a22b9eSJed Brown     // -- Stabilization method: none, SU, or SUPG
509bb8a0c61SJames Wright     CeedScalar stab[5][3] = {{0.}};
510139613f2SLeila Ghaffari     switch (context->stabilization) {
511139613f2SLeila Ghaffari       case 0:  // Galerkin
512139613f2SLeila Ghaffari         break;
513139613f2SLeila Ghaffari       case 1:  // SU
5142b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
5152b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
5162b916ea7SJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
5172b916ea7SJeremy L Thompson           }
5182b916ea7SJeremy L Thompson         }
519139613f2SLeila Ghaffari 
5202b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
5212b916ea7SJeremy 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]);
5222b916ea7SJeremy L Thompson         }
523139613f2SLeila Ghaffari         break;
524139613f2SLeila Ghaffari       case 2:  // SUPG
5252b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
5262b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
5272b916ea7SJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
5282b916ea7SJeremy L Thompson           }
5292b916ea7SJeremy L Thompson         }
530139613f2SLeila Ghaffari 
5312b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
5322b916ea7SJeremy 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]);
5332b916ea7SJeremy L Thompson         }
534139613f2SLeila Ghaffari         break;
535139613f2SLeila Ghaffari     }
536a515125bSLeila Ghaffari   }  // End Quadrature Point Loop
537a515125bSLeila Ghaffari 
538a515125bSLeila Ghaffari   // Return
539a515125bSLeila Ghaffari   return 0;
540a515125bSLeila Ghaffari }
541a515125bSLeila Ghaffari // *****************************************************************************
54204e40bb6SJeremy L Thompson // This QFunction sets the inflow boundary conditions for the traveling vortex problem.
543a515125bSLeila Ghaffari //
54404e40bb6SJeremy L Thompson //  Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly.
545a515125bSLeila Ghaffari // *****************************************************************************
5462b916ea7SJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
547a515125bSLeila Ghaffari   // Inputs
548*ade49511SJames Wright   const CeedScalar(*q_data_sur) = in[2];
549a515125bSLeila Ghaffari   // Outputs
550a515125bSLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
551a515125bSLeila Ghaffari   EulerContext     context       = (EulerContext)ctx;
552a515125bSLeila Ghaffari   const int        euler_test    = context->euler_test;
553*ade49511SJames Wright   const bool       is_implicit   = context->implicit;
554a515125bSLeila Ghaffari   CeedScalar      *mean_velocity = context->mean_velocity;
555a515125bSLeila Ghaffari   const CeedScalar cv            = 2.5;
556a515125bSLeila Ghaffari   const CeedScalar R             = 1.;
557a515125bSLeila Ghaffari   CeedScalar       T_inlet;
558a515125bSLeila Ghaffari   CeedScalar       P_inlet;
559a515125bSLeila Ghaffari 
560a515125bSLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
5612b916ea7SJeremy L Thompson   if (euler_test == 1 || euler_test == 3) {
562a515125bSLeila Ghaffari     for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.;
5632b916ea7SJeremy L Thompson   }
564a515125bSLeila Ghaffari 
565a515125bSLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
566a515125bSLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
567a515125bSLeila Ghaffari   else T_inlet = P_inlet = 1.;
568a515125bSLeila Ghaffari 
569a515125bSLeila Ghaffari   // Quadrature Point Loop
5703d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
571*ade49511SJames Wright     CeedScalar wdetJb, norm[3];
572*ade49511SJames Wright     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
573*ade49511SJames Wright     wdetJb *= is_implicit ? -1. : 1.;
574a515125bSLeila Ghaffari 
575a515125bSLeila Ghaffari     // face_normal = Normal vector of the face
5762b916ea7SJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
577a515125bSLeila Ghaffari     // The Physics
578a515125bSLeila Ghaffari     // Zero v so all future terms can safely sum into it
579493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
580a515125bSLeila Ghaffari 
581a515125bSLeila Ghaffari     // Implementing in/outflow BCs
582002797a3SLeila Ghaffari     if (face_normal > 0) {
583a515125bSLeila Ghaffari     } else {  // inflow
584a515125bSLeila Ghaffari       const CeedScalar rho_inlet       = P_inlet / (R * T_inlet);
5852b916ea7SJeremy L Thompson       const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.;
586a515125bSLeila Ghaffari       // incoming total energy
587a515125bSLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
588a515125bSLeila Ghaffari 
589a515125bSLeila Ghaffari       // The Physics
590a515125bSLeila Ghaffari       // -- Density
591a515125bSLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
592a515125bSLeila Ghaffari 
593a515125bSLeila Ghaffari       // -- Momentum
5942b916ea7SJeremy 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);
595a515125bSLeila Ghaffari 
596a515125bSLeila Ghaffari       // -- Total Energy Density
597a515125bSLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
598a515125bSLeila Ghaffari     }
599a515125bSLeila Ghaffari 
600a515125bSLeila Ghaffari   }  // End Quadrature Point Loop
601a515125bSLeila Ghaffari   return 0;
602a515125bSLeila Ghaffari }
603a515125bSLeila Ghaffari 
604a515125bSLeila Ghaffari // *****************************************************************************
60504e40bb6SJeremy L Thompson // This QFunction sets the outflow boundary conditions for the Euler solver.
60668ef3d20SLeila Ghaffari //
60768ef3d20SLeila Ghaffari //  Outflow BCs:
60804e40bb6SJeremy L Thompson //    The validity of the weak form of the governing equations is extended to the outflow.
60968ef3d20SLeila Ghaffari // *****************************************************************************
6102b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
61168ef3d20SLeila Ghaffari   // Inputs
6123d65b166SJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
613*ade49511SJames Wright   const CeedScalar(*q_data_sur)    = in[2];
6143d65b166SJames Wright 
61568ef3d20SLeila Ghaffari   // Outputs
61668ef3d20SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
61768ef3d20SLeila Ghaffari   EulerContext context       = (EulerContext)ctx;
618*ade49511SJames Wright   const bool   is_implicit   = context->implicit;
61968ef3d20SLeila Ghaffari   CeedScalar  *mean_velocity = context->mean_velocity;
62068ef3d20SLeila Ghaffari 
62168ef3d20SLeila Ghaffari   const CeedScalar gamma = 1.4;
62268ef3d20SLeila Ghaffari 
62368ef3d20SLeila Ghaffari   // Quadrature Point Loop
6243d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
62568ef3d20SLeila Ghaffari     // Setup
62668ef3d20SLeila Ghaffari     // -- Interp in
62768ef3d20SLeila Ghaffari     const CeedScalar rho  = q[0][i];
6282b916ea7SJeremy L Thompson     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
62968ef3d20SLeila Ghaffari     const CeedScalar E    = q[4][i];
63068ef3d20SLeila Ghaffari 
631*ade49511SJames Wright     CeedScalar wdetJb, norm[3];
632*ade49511SJames Wright     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
633*ade49511SJames Wright     wdetJb *= is_implicit ? -1. : 1.;
63468ef3d20SLeila Ghaffari 
63568ef3d20SLeila Ghaffari     // face_normal = Normal vector of the face
6362b916ea7SJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
63768ef3d20SLeila Ghaffari     // The Physics
63868ef3d20SLeila Ghaffari     // Zero v so all future terms can safely sum into it
639493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0;
64068ef3d20SLeila Ghaffari 
64168ef3d20SLeila Ghaffari     // Implementing in/outflow BCs
64268ef3d20SLeila Ghaffari     if (face_normal > 0) {  // outflow
64368ef3d20SLeila Ghaffari       const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.;
64468ef3d20SLeila Ghaffari       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.);              // pressure
6452b916ea7SJeremy L Thompson       const CeedScalar u_normal  = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2];  // Normal velocity
64668ef3d20SLeila Ghaffari       // The Physics
64768ef3d20SLeila Ghaffari       // -- Density
64868ef3d20SLeila Ghaffari       v[0][i] -= wdetJb * rho * u_normal;
64968ef3d20SLeila Ghaffari 
65068ef3d20SLeila Ghaffari       // -- Momentum
6512b916ea7SJeremy L Thompson       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P);
65268ef3d20SLeila Ghaffari 
65368ef3d20SLeila Ghaffari       // -- Total Energy Density
65468ef3d20SLeila Ghaffari       v[4][i] -= wdetJb * u_normal * (E + P);
65568ef3d20SLeila Ghaffari     }
65668ef3d20SLeila Ghaffari   }  // End Quadrature Point Loop
65768ef3d20SLeila Ghaffari   return 0;
65868ef3d20SLeila Ghaffari }
65968ef3d20SLeila Ghaffari 
66068ef3d20SLeila Ghaffari // *****************************************************************************
661a515125bSLeila Ghaffari 
662a515125bSLeila Ghaffari #endif  // eulervortex_h
663