xref: /honee/qfunctions/eulervortex.h (revision 3e17a7a1849b799c24669fa5a24aea7b05055bba)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5a515125bSLeila Ghaffari /// Euler traveling vortex initial condition and operator for Navier-Stokes
6a515125bSLeila Ghaffari /// example using PETSc
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari // Model from:
904e40bb6SJeremy L Thompson //   On the Order of Accuracy and Numerical Performance of Two Classes of Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011).
10*3e17a7a1SJames Wright #include <ceed/types.h>
11*3e17a7a1SJames Wright #ifndef CEED_RUNNING_JIT_PASS
12*3e17a7a1SJames Wright #include <stdbool.h>
13*3e17a7a1SJames Wright #endif
142b916ea7SJeremy L Thompson 
15704b8bbeSJames Wright #include "utils.h"
16a515125bSLeila Ghaffari 
17a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext;
18a515125bSLeila Ghaffari struct EulerContext_ {
19a515125bSLeila Ghaffari   CeedScalar center[3];
20a515125bSLeila Ghaffari   CeedScalar curr_time;
21a515125bSLeila Ghaffari   CeedScalar vortex_strength;
22d8a22b9eSJed Brown   CeedScalar c_tau;
23a515125bSLeila Ghaffari   CeedScalar mean_velocity[3];
24a515125bSLeila Ghaffari   bool       implicit;
25139613f2SLeila Ghaffari   int        euler_test;
26139613f2SLeila Ghaffari   int        stabilization;  // See StabilizationType: 0=none, 1=SU, 2=SUPG
27a515125bSLeila Ghaffari };
28a515125bSLeila Ghaffari 
29a515125bSLeila Ghaffari // *****************************************************************************
30a515125bSLeila Ghaffari // This function sets the initial conditions
31a515125bSLeila Ghaffari //
32a515125bSLeila Ghaffari //   Temperature:
33a515125bSLeila Ghaffari //     T   = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
34a515125bSLeila Ghaffari //   Density:
35a515125bSLeila Ghaffari //     rho = (T/S_vortex)^(1 / (gamma - 1))
36a515125bSLeila Ghaffari //   Pressure:
37a515125bSLeila Ghaffari //     P   = rho * T
38a515125bSLeila Ghaffari //   Velocity:
39a515125bSLeila Ghaffari //     ui  = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
40a515125bSLeila Ghaffari //     r   = sqrt( (x - xc)**2 + (y - yc)**2 )
41a515125bSLeila Ghaffari //   Velocity/Momentum Density:
42a515125bSLeila Ghaffari //     Ui  = rho ui
43a515125bSLeila Ghaffari //   Total Energy:
44a515125bSLeila Ghaffari //     E   = P / (gamma - 1) + rho (u u)/2
45a515125bSLeila Ghaffari //
46a515125bSLeila Ghaffari // Constants:
47a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
48a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
49a515125bSLeila Ghaffari //   vortex_strength ,  Strength of vortex
50a515125bSLeila Ghaffari //   center          ,  Location of bubble center
51a515125bSLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
52a515125bSLeila Ghaffari //
53a515125bSLeila Ghaffari // *****************************************************************************
54a515125bSLeila Ghaffari 
55a515125bSLeila Ghaffari // *****************************************************************************
5604e40bb6SJeremy L Thompson // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling
5704e40bb6SJeremy L Thompson // vortex
58a515125bSLeila Ghaffari // *****************************************************************************
592b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
60a515125bSLeila Ghaffari   // Context
61a515125bSLeila Ghaffari   const EulerContext context         = (EulerContext)ctx;
62a515125bSLeila Ghaffari   const CeedScalar   vortex_strength = context->vortex_strength;
63a515125bSLeila Ghaffari   const CeedScalar  *center          = context->center;  // Center of the domain
64a515125bSLeila Ghaffari   const CeedScalar  *mean_velocity   = context->mean_velocity;
65a515125bSLeila Ghaffari 
66a515125bSLeila Ghaffari   // Setup
67a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
68a515125bSLeila Ghaffari   const CeedScalar cv    = 2.5;
69a515125bSLeila Ghaffari   const CeedScalar R     = 1.;
70a515125bSLeila Ghaffari   const CeedScalar x = X[0], y = X[1];  // Coordinates
71a515125bSLeila Ghaffari   // Vortex center
72a515125bSLeila Ghaffari   const CeedScalar xc = center[0] + mean_velocity[0] * time;
73a515125bSLeila Ghaffari   const CeedScalar yc = center[1] + mean_velocity[1] * time;
74a515125bSLeila Ghaffari 
75a515125bSLeila Ghaffari   const CeedScalar x0       = x - xc;
76a515125bSLeila Ghaffari   const CeedScalar y0       = y - yc;
7774960ff5SJames Wright   const CeedScalar r        = sqrt(Square(x0) + Square(y0));
78a515125bSLeila Ghaffari   const CeedScalar C        = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI);
792b916ea7SJeremy L Thompson   const CeedScalar delta_T  = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI);
80a515125bSLeila Ghaffari   const CeedScalar S_vortex = 1;  // no perturbation in the entropy P / rho^gamma
812b916ea7SJeremy L Thompson   const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI);
82a515125bSLeila Ghaffari   CeedScalar       rho, P, T, E, u[3] = {0.};
83a515125bSLeila Ghaffari 
84a515125bSLeila Ghaffari   // Initial Conditions
85a515125bSLeila Ghaffari   switch (context->euler_test) {
86a515125bSLeila Ghaffari     case 0:  // Traveling vortex
87a515125bSLeila Ghaffari       T = 1 + delta_T;
88a515125bSLeila Ghaffari       // P = rho * T
89a515125bSLeila Ghaffari       // P = S * rho^gamma
90a515125bSLeila Ghaffari       // Solve for rho, then substitute for P
91139613f2SLeila Ghaffari       rho  = pow(T / S_vortex, 1 / (gamma - 1.));
92a515125bSLeila Ghaffari       P    = rho * T;
93a515125bSLeila Ghaffari       u[0] = mean_velocity[0] - C * y0;
94a515125bSLeila Ghaffari       u[1] = mean_velocity[1] + C * x0;
95a515125bSLeila Ghaffari 
96a515125bSLeila Ghaffari       // Assign exact solution
97a515125bSLeila Ghaffari       q[0] = rho;
98a515125bSLeila Ghaffari       q[1] = rho * u[0];
99a515125bSLeila Ghaffari       q[2] = rho * u[1];
100a515125bSLeila Ghaffari       q[3] = rho * u[2];
101a515125bSLeila Ghaffari       q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.;
102a515125bSLeila Ghaffari       break;
103a515125bSLeila Ghaffari     case 1:  // Constant zero velocity, density constant, total energy constant
104a515125bSLeila Ghaffari       rho = 1.;
105a515125bSLeila Ghaffari       E   = 2.;
106a515125bSLeila Ghaffari 
107a515125bSLeila Ghaffari       // Assign exact solution
108a515125bSLeila Ghaffari       q[0] = rho;
109a515125bSLeila Ghaffari       q[1] = rho * u[0];
110a515125bSLeila Ghaffari       q[2] = rho * u[1];
111a515125bSLeila Ghaffari       q[3] = rho * u[2];
112a515125bSLeila Ghaffari       q[4] = E;
113a515125bSLeila Ghaffari       break;
114a515125bSLeila Ghaffari     case 2:  // Constant nonzero velocity, density constant, total energy constant
115a515125bSLeila Ghaffari       rho  = 1.;
116a515125bSLeila Ghaffari       E    = 2.;
117a515125bSLeila Ghaffari       u[0] = mean_velocity[0];
118a515125bSLeila Ghaffari       u[1] = mean_velocity[1];
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] = E;
126a515125bSLeila Ghaffari       break;
12704e40bb6SJeremy 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
12804e40bb6SJeremy L Thompson              // bubble won't diffuse
129a515125bSLeila Ghaffari       // (for Euler, where there is no thermal conductivity)
130a515125bSLeila Ghaffari       P   = 1.;
131a515125bSLeila Ghaffari       T   = 1. - S_bubble * exp(1. - r * r);
132a515125bSLeila Ghaffari       rho = P / (R * T);
133a515125bSLeila Ghaffari 
134a515125bSLeila Ghaffari       // Assign exact solution
135a515125bSLeila Ghaffari       q[0] = rho;
136a515125bSLeila Ghaffari       q[1] = rho * u[0];
137a515125bSLeila Ghaffari       q[2] = rho * u[1];
138a515125bSLeila Ghaffari       q[3] = rho * u[2];
139a515125bSLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
140a515125bSLeila Ghaffari       break;
14104e40bb6SJeremy L Thompson     case 4:  // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant),
14204e40bb6SJeremy L Thompson       // It should be transported across the domain, but velocity stays constant
143a515125bSLeila Ghaffari       P    = 1.;
144a515125bSLeila Ghaffari       T    = 1. - S_bubble * exp(1. - r * r);
145a515125bSLeila Ghaffari       rho  = P / (R * T);
146a515125bSLeila Ghaffari       u[0] = mean_velocity[0];
147a515125bSLeila Ghaffari       u[1] = mean_velocity[1];
148a515125bSLeila Ghaffari 
149a515125bSLeila Ghaffari       // Assign exact solution
150a515125bSLeila Ghaffari       q[0] = rho;
151a515125bSLeila Ghaffari       q[1] = rho * u[0];
152a515125bSLeila Ghaffari       q[2] = rho * u[1];
153a515125bSLeila Ghaffari       q[3] = rho * u[2];
154a515125bSLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
155a515125bSLeila Ghaffari       break;
1560df2634dSLeila Ghaffari     case 5:  // non-smooth thermal bubble - cylinder
1570df2634dSLeila Ghaffari       P    = 1.;
1580df2634dSLeila Ghaffari       T    = 1. - (r < 1. ? S_bubble : 0.);
1590df2634dSLeila Ghaffari       rho  = P / (R * T);
1600df2634dSLeila Ghaffari       u[0] = mean_velocity[0];
1610df2634dSLeila Ghaffari       u[1] = mean_velocity[1];
1620df2634dSLeila Ghaffari 
1630df2634dSLeila Ghaffari       // Assign exact solution
1640df2634dSLeila Ghaffari       q[0] = rho;
1650df2634dSLeila Ghaffari       q[1] = rho * u[0];
1660df2634dSLeila Ghaffari       q[2] = rho * u[1];
1670df2634dSLeila Ghaffari       q[3] = rho * u[2];
1680df2634dSLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
1690df2634dSLeila Ghaffari       break;
170a515125bSLeila Ghaffari   }
171a515125bSLeila Ghaffari   return 0;
172a515125bSLeila Ghaffari }
173a515125bSLeila Ghaffari 
174a515125bSLeila Ghaffari // *****************************************************************************
175139613f2SLeila Ghaffari // Helper function for computing flux Jacobian
176139613f2SLeila Ghaffari // *****************************************************************************
1772b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
178139613f2SLeila Ghaffari                                                         const CeedScalar gamma) {
179139613f2SLeila Ghaffari   CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2];  // Velocity square
180139613f2SLeila Ghaffari   for (CeedInt i = 0; i < 3; i++) {                           // Jacobian matrices for 3 directions
181139613f2SLeila Ghaffari     for (CeedInt j = 0; j < 3; j++) {                         // Rows of each Jacobian matrix
182139613f2SLeila Ghaffari       dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j];
183139613f2SLeila Ghaffari       for (CeedInt k = 0; k < 3; k++) {  // Columns of each Jacobian matrix
184139613f2SLeila Ghaffari         dF[i][0][k + 1]     = ((i == k) ? 1. : 0.);
1852b916ea7SJeremy 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.);
1862b916ea7SJeremy L Thompson         dF[i][4][k + 1]     = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k];
187139613f2SLeila Ghaffari       }
188139613f2SLeila Ghaffari       dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.);
189139613f2SLeila Ghaffari     }
190139613f2SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho);
191139613f2SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
192139613f2SLeila Ghaffari   }
193139613f2SLeila Ghaffari }
194139613f2SLeila Ghaffari 
195139613f2SLeila Ghaffari // *****************************************************************************
196d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant)
197d8a22b9eSJed Brown //   Model from:
198d8a22b9eSJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
199d8a22b9eSJed Brown //
200d8a22b9eSJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
201d8a22b9eSJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
202d8a22b9eSJed Brown //
203d8a22b9eSJed Brown // Where
204d8a22b9eSJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
205d8a22b9eSJed Brown //   h[i]      = 2 length(dxdX[i])
206d8a22b9eSJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
207d8a22b9eSJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
20804e40bb6SJeremy L Thompson //   rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i
209d8a22b9eSJed Brown // *****************************************************************************
2102b916ea7SJeremy 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,
2112b916ea7SJeremy L Thompson                                        const CeedScalar c_tau) {
212493642f1SJames Wright   for (CeedInt i = 0; i < 3; i++) {
213d8a22b9eSJed Brown     // length of element in direction i
21474960ff5SJames Wright     CeedScalar h = 2 / sqrt(Square(dXdx[0][i]) + Square(dXdx[1][i]) + Square(dXdx[2][i]));
215d8a22b9eSJed Brown     // fastest wave in direction i
216d8a22b9eSJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
217d8a22b9eSJed Brown     Tau_x[i]                = c_tau * h / fastest_wave;
218d8a22b9eSJed Brown   }
219d8a22b9eSJed Brown }
220d8a22b9eSJed Brown 
221d8a22b9eSJed Brown // *****************************************************************************
222a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
223a515125bSLeila Ghaffari // *****************************************************************************
2242b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
225a515125bSLeila Ghaffari   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
226a515125bSLeila Ghaffari   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
227b193fadcSJames Wright 
228a515125bSLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
229a515125bSLeila Ghaffari 
2303d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
231a515125bSLeila Ghaffari     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
232139613f2SLeila Ghaffari     CeedScalar       q[5] = {0.};
233a515125bSLeila Ghaffari 
234a515125bSLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
2352b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
236b193fadcSJames Wright   }
237a515125bSLeila Ghaffari   return 0;
238a515125bSLeila Ghaffari }
239a515125bSLeila Ghaffari 
240a515125bSLeila Ghaffari // *****************************************************************************
24104e40bb6SJeremy L Thompson // This QFunction implements the following formulation of Euler equations with explicit time stepping method
242a515125bSLeila Ghaffari //
24304e40bb6SJeremy L Thompson // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density.
244a515125bSLeila Ghaffari //
245a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
246a515125bSLeila Ghaffari //   rho - Mass Density
247a515125bSLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
248a515125bSLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
249a515125bSLeila Ghaffari //
250a515125bSLeila Ghaffari // Euler Equations:
251a515125bSLeila Ghaffari //   drho/dt + div( U )                   = 0
252a515125bSLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
253a515125bSLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
254a515125bSLeila Ghaffari //
255a515125bSLeila Ghaffari // Equation of State:
256a515125bSLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
257a515125bSLeila Ghaffari //
258a515125bSLeila Ghaffari // Constants:
259a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
260a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
261a515125bSLeila Ghaffari //   g               ,  Gravity
262a515125bSLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
263a515125bSLeila Ghaffari // *****************************************************************************
2642b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
2653d65b166SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
2663d65b166SJames Wright   const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
267ade49511SJames Wright   const CeedScalar(*q_data)            = in[2];
2683d65b166SJames Wright   CeedScalar(*v)[CEED_Q_VLA]           = (CeedScalar(*)[CEED_Q_VLA])out[0];
2693d65b166SJames Wright   CeedScalar(*dv)[5][CEED_Q_VLA]       = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
270a515125bSLeila Ghaffari 
271139613f2SLeila Ghaffari   EulerContext     context = (EulerContext)ctx;
272d8a22b9eSJed Brown   const CeedScalar c_tau   = context->c_tau;
273a515125bSLeila Ghaffari   const CeedScalar gamma   = 1.4;
274a515125bSLeila Ghaffari 
2753d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
276a515125bSLeila Ghaffari     // Setup
277a515125bSLeila Ghaffari     // -- Interp in
278a515125bSLeila Ghaffari     const CeedScalar rho      = q[0][i];
2792b916ea7SJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
280a515125bSLeila Ghaffari     const CeedScalar E        = q[4][i];
2812b916ea7SJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
2822b916ea7SJeremy L Thompson     const CeedScalar dU[3][3] = {
2832b916ea7SJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
2842b916ea7SJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
2852b916ea7SJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
286139613f2SLeila Ghaffari     };
2872b916ea7SJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
288ade49511SJames Wright     CeedScalar       wdetJ, dXdx[3][3];
289ade49511SJames Wright     QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
290139613f2SLeila Ghaffari     // dU/dx
291139613f2SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
292139613f2SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
293139613f2SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
294139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
295493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
296493642f1SJames Wright       for (CeedInt k = 0; k < 3; k++) {
297139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
298139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
299493642f1SJames Wright         for (CeedInt l = 0; l < 3; l++) {
300139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
301139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
302139613f2SLeila Ghaffari         }
303139613f2SLeila Ghaffari       }
304139613f2SLeila Ghaffari     }
305139613f2SLeila Ghaffari     // Pressure
3062b916ea7SJeremy 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,
307139613f2SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
308a515125bSLeila Ghaffari 
309a515125bSLeila Ghaffari     // The Physics
310a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
311493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) {
312139613f2SLeila Ghaffari       v[j][i] = 0.;
3132b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
314a515125bSLeila Ghaffari     }
315a515125bSLeila Ghaffari 
316a515125bSLeila Ghaffari     // -- Density
317a515125bSLeila Ghaffari     // ---- u rho
3182b916ea7SJeremy 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]);
319a515125bSLeila Ghaffari     // -- Momentum
320a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
3212b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3222b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
3232b916ea7SJeremy 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] +
324139613f2SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
3252b916ea7SJeremy L Thompson       }
3262b916ea7SJeremy L Thompson     }
327a515125bSLeila Ghaffari     // -- Total Energy Density
328a515125bSLeila Ghaffari     // ---- (E + P) u
3292b916ea7SJeremy 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]);
330139613f2SLeila Ghaffari 
331139613f2SLeila Ghaffari     // --Stabilization terms
332139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
333139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
334d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
335139613f2SLeila Ghaffari 
336139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
337139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
338493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
339139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
340139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
3412b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
342139613f2SLeila Ghaffari     }
343139613f2SLeila Ghaffari 
344139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
345139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
3462b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3472b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
3482b916ea7SJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
3492b916ea7SJeremy L Thompson       }
3502b916ea7SJeremy L Thompson     }
351139613f2SLeila Ghaffari 
352d8a22b9eSJed Brown     // Stabilization
353d8a22b9eSJed Brown     // -- Tau elements
354d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
355d8a22b9eSJed Brown     CeedScalar       Tau_x[3]    = {0.};
356d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
357139613f2SLeila Ghaffari 
358d8a22b9eSJed Brown     // -- Stabilization method: none or SU
359bb8a0c61SJames Wright     CeedScalar stab[5][3] = {{0.}};
360139613f2SLeila Ghaffari     switch (context->stabilization) {
361139613f2SLeila Ghaffari       case 0:  // Galerkin
362139613f2SLeila Ghaffari         break;
363139613f2SLeila Ghaffari       case 1:  // SU
3642b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
3652b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
3662b916ea7SJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
3672b916ea7SJeremy L Thompson           }
3682b916ea7SJeremy L Thompson         }
369139613f2SLeila Ghaffari 
3702b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
3712b916ea7SJeremy 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]);
3722b916ea7SJeremy L Thompson         }
373139613f2SLeila Ghaffari         break;
374139613f2SLeila Ghaffari       case 2:  // SUPG is not implemented for explicit scheme
375139613f2SLeila Ghaffari         break;
376139613f2SLeila Ghaffari     }
377b193fadcSJames Wright   }
378a515125bSLeila Ghaffari   return 0;
379a515125bSLeila Ghaffari }
380a515125bSLeila Ghaffari 
381a515125bSLeila Ghaffari // *****************************************************************************
38204e40bb6SJeremy L Thompson // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method
383a515125bSLeila Ghaffari // *****************************************************************************
3842b916ea7SJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
3853d65b166SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
3863d65b166SJames Wright   const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
3873d65b166SJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
388ade49511SJames Wright   const CeedScalar(*q_data)            = in[3];
3893d65b166SJames Wright   CeedScalar(*v)[CEED_Q_VLA]           = (CeedScalar(*)[CEED_Q_VLA])out[0];
3903d65b166SJames Wright   CeedScalar(*dv)[5][CEED_Q_VLA]       = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
391a515125bSLeila Ghaffari 
392139613f2SLeila Ghaffari   EulerContext     context = (EulerContext)ctx;
393d8a22b9eSJed Brown   const CeedScalar c_tau   = context->c_tau;
394a515125bSLeila Ghaffari   const CeedScalar gamma   = 1.4;
395a515125bSLeila Ghaffari 
3963d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
397a515125bSLeila Ghaffari     // Setup
398a515125bSLeila Ghaffari     // -- Interp in
399a515125bSLeila Ghaffari     const CeedScalar rho      = q[0][i];
4002b916ea7SJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
401a515125bSLeila Ghaffari     const CeedScalar E        = q[4][i];
4022b916ea7SJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
4032b916ea7SJeremy L Thompson     const CeedScalar dU[3][3] = {
4042b916ea7SJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
4052b916ea7SJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
4062b916ea7SJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
407139613f2SLeila Ghaffari     };
4082b916ea7SJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
409ade49511SJames Wright     CeedScalar       wdetJ, dXdx[3][3];
410ade49511SJames Wright     QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
411139613f2SLeila Ghaffari     // dU/dx
412139613f2SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
413139613f2SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
414139613f2SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
415139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
416493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
417493642f1SJames Wright       for (CeedInt k = 0; k < 3; k++) {
418139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
419139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
420493642f1SJames Wright         for (CeedInt l = 0; l < 3; l++) {
421139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
422139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
423139613f2SLeila Ghaffari         }
424139613f2SLeila Ghaffari       }
425139613f2SLeila Ghaffari     }
4262b916ea7SJeremy 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,
427139613f2SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
428a515125bSLeila Ghaffari 
429a515125bSLeila Ghaffari     // The Physics
430a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
431493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) {
432139613f2SLeila Ghaffari       v[j][i] = 0.;
4332b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
434a515125bSLeila Ghaffari     }
435a515125bSLeila Ghaffari     //-----mass matrix
4362b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i];
437a515125bSLeila Ghaffari 
438a515125bSLeila Ghaffari     // -- Density
439a515125bSLeila Ghaffari     // ---- u rho
4402b916ea7SJeremy 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]);
441a515125bSLeila Ghaffari     // -- Momentum
442a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
4432b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
4442b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
4452b916ea7SJeremy 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] +
446139613f2SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
4472b916ea7SJeremy L Thompson       }
4482b916ea7SJeremy L Thompson     }
449a515125bSLeila Ghaffari     // -- Total Energy Density
450a515125bSLeila Ghaffari     // ---- (E + P) u
4512b916ea7SJeremy 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]);
452139613f2SLeila Ghaffari 
453139613f2SLeila Ghaffari     // -- Stabilization terms
454139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
455139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
456d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
457139613f2SLeila Ghaffari 
458139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
459139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
460493642f1SJames Wright     for (CeedInt j = 0; j < 3; j++) {
461139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
462139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
4632b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
464139613f2SLeila Ghaffari     }
465139613f2SLeila Ghaffari 
466139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
467139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
4682b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
4692b916ea7SJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
4702b916ea7SJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
4712b916ea7SJeremy L Thompson       }
4722b916ea7SJeremy L Thompson     }
473139613f2SLeila Ghaffari 
474139613f2SLeila Ghaffari     // ---- Strong residual
475139613f2SLeila Ghaffari     CeedScalar strong_res[5];
4762b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j];
477139613f2SLeila Ghaffari 
478d8a22b9eSJed Brown     // Stabilization
479d8a22b9eSJed Brown     // -- Tau elements
480d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
481d8a22b9eSJed Brown     CeedScalar       Tau_x[3]    = {0.};
482d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
483139613f2SLeila Ghaffari 
484d8a22b9eSJed Brown     // -- Stabilization method: none, SU, or SUPG
485bb8a0c61SJames Wright     CeedScalar stab[5][3] = {{0.}};
486139613f2SLeila Ghaffari     switch (context->stabilization) {
487139613f2SLeila Ghaffari       case 0:  // Galerkin
488139613f2SLeila Ghaffari         break;
489139613f2SLeila Ghaffari       case 1:  // SU
4902b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
4912b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
4922b916ea7SJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
4932b916ea7SJeremy L Thompson           }
4942b916ea7SJeremy L Thompson         }
495139613f2SLeila Ghaffari 
4962b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
4972b916ea7SJeremy 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]);
4982b916ea7SJeremy L Thompson         }
499139613f2SLeila Ghaffari         break;
500139613f2SLeila Ghaffari       case 2:  // SUPG
5012b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
5022b916ea7SJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
5032b916ea7SJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
5042b916ea7SJeremy L Thompson           }
5052b916ea7SJeremy L Thompson         }
506139613f2SLeila Ghaffari 
5072b916ea7SJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
5082b916ea7SJeremy 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]);
5092b916ea7SJeremy L Thompson         }
510139613f2SLeila Ghaffari         break;
511139613f2SLeila Ghaffari     }
512b193fadcSJames Wright   }
513a515125bSLeila Ghaffari   return 0;
514a515125bSLeila Ghaffari }
515a515125bSLeila Ghaffari // *****************************************************************************
51604e40bb6SJeremy L Thompson // This QFunction sets the inflow boundary conditions for the traveling vortex problem.
517a515125bSLeila Ghaffari //
51804e40bb6SJeremy L Thompson //  Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly.
519a515125bSLeila Ghaffari // *****************************************************************************
5202b916ea7SJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
521ade49511SJames Wright   const CeedScalar(*q_data_sur) = in[2];
522a515125bSLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA]    = (CeedScalar(*)[CEED_Q_VLA])out[0];
523b193fadcSJames Wright 
524a515125bSLeila Ghaffari   EulerContext     context       = (EulerContext)ctx;
525a515125bSLeila Ghaffari   const int        euler_test    = context->euler_test;
526ade49511SJames Wright   const bool       is_implicit   = context->implicit;
527a515125bSLeila Ghaffari   CeedScalar      *mean_velocity = context->mean_velocity;
528a515125bSLeila Ghaffari   const CeedScalar cv            = 2.5;
529a515125bSLeila Ghaffari   const CeedScalar R             = 1.;
530a515125bSLeila Ghaffari   CeedScalar       T_inlet;
531a515125bSLeila Ghaffari   CeedScalar       P_inlet;
532a515125bSLeila Ghaffari 
533a515125bSLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
5342b916ea7SJeremy L Thompson   if (euler_test == 1 || euler_test == 3) {
535a515125bSLeila Ghaffari     for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.;
5362b916ea7SJeremy L Thompson   }
537a515125bSLeila Ghaffari 
538a515125bSLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
539a515125bSLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
540a515125bSLeila Ghaffari   else T_inlet = P_inlet = 1.;
541a515125bSLeila Ghaffari 
5423d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
54378e8b7daSJames Wright     CeedScalar wdetJb, normal[3];
54478e8b7daSJames Wright     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, normal);
545ade49511SJames Wright     wdetJb *= is_implicit ? -1. : 1.;
546a515125bSLeila Ghaffari 
547a515125bSLeila Ghaffari     // face_normal = Normal vector of the face
54878e8b7daSJames Wright     const CeedScalar face_normal = Dot3(normal, mean_velocity);
549a515125bSLeila Ghaffari     // The Physics
550a515125bSLeila Ghaffari     // Zero v so all future terms can safely sum into it
551493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
552a515125bSLeila Ghaffari 
553a515125bSLeila Ghaffari     // Implementing in/outflow BCs
554002797a3SLeila Ghaffari     if (face_normal > 0) {
555a515125bSLeila Ghaffari     } else {  // inflow
556a515125bSLeila Ghaffari       const CeedScalar rho_inlet       = P_inlet / (R * T_inlet);
5572b916ea7SJeremy L Thompson       const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.;
558a515125bSLeila Ghaffari       // incoming total energy
559a515125bSLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
560a515125bSLeila Ghaffari 
561a515125bSLeila Ghaffari       // The Physics
562a515125bSLeila Ghaffari       // -- Density
563a515125bSLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
564a515125bSLeila Ghaffari 
565a515125bSLeila Ghaffari       // -- Momentum
56678e8b7daSJames Wright       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho_inlet * face_normal * mean_velocity[j] + normal[j] * P_inlet);
567a515125bSLeila Ghaffari 
568a515125bSLeila Ghaffari       // -- Total Energy Density
569a515125bSLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
570a515125bSLeila Ghaffari     }
571b193fadcSJames Wright   }
572a515125bSLeila Ghaffari   return 0;
573a515125bSLeila Ghaffari }
574a515125bSLeila Ghaffari 
575a515125bSLeila Ghaffari // *****************************************************************************
57604e40bb6SJeremy L Thompson // This QFunction sets the outflow boundary conditions for the Euler solver.
57768ef3d20SLeila Ghaffari //
57868ef3d20SLeila Ghaffari //  Outflow BCs:
57904e40bb6SJeremy L Thompson //    The validity of the weak form of the governing equations is extended to the outflow.
58068ef3d20SLeila Ghaffari // *****************************************************************************
5812b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
5823d65b166SJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
583ade49511SJames Wright   const CeedScalar(*q_data_sur)    = in[2];
58468ef3d20SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA]       = (CeedScalar(*)[CEED_Q_VLA])out[0];
585b193fadcSJames Wright 
58668ef3d20SLeila Ghaffari   EulerContext context       = (EulerContext)ctx;
587ade49511SJames Wright   const bool   is_implicit   = context->implicit;
58868ef3d20SLeila Ghaffari   CeedScalar  *mean_velocity = context->mean_velocity;
58968ef3d20SLeila Ghaffari 
59068ef3d20SLeila Ghaffari   const CeedScalar gamma = 1.4;
59168ef3d20SLeila Ghaffari 
5923d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
59368ef3d20SLeila Ghaffari     // Setup
59468ef3d20SLeila Ghaffari     // -- Interp in
59568ef3d20SLeila Ghaffari     const CeedScalar rho  = q[0][i];
5962b916ea7SJeremy L Thompson     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
59768ef3d20SLeila Ghaffari     const CeedScalar E    = q[4][i];
59868ef3d20SLeila Ghaffari 
59978e8b7daSJames Wright     CeedScalar wdetJb, normal[3];
60078e8b7daSJames Wright     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, normal);
601ade49511SJames Wright     wdetJb *= is_implicit ? -1. : 1.;
60268ef3d20SLeila Ghaffari 
60368ef3d20SLeila Ghaffari     // face_normal = Normal vector of the face
60478e8b7daSJames Wright     const CeedScalar face_normal = Dot3(normal, mean_velocity);
60568ef3d20SLeila Ghaffari     // The Physics
60668ef3d20SLeila Ghaffari     // Zero v so all future terms can safely sum into it
607493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0;
60868ef3d20SLeila Ghaffari 
60968ef3d20SLeila Ghaffari     // Implementing in/outflow BCs
61068ef3d20SLeila Ghaffari     if (face_normal > 0) {  // outflow
61168ef3d20SLeila Ghaffari       const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.;
61268ef3d20SLeila Ghaffari       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.);  // pressure
61378e8b7daSJames Wright       const CeedScalar u_normal  = Dot3(normal, u);                       // Normal velocity
61468ef3d20SLeila Ghaffari       // The Physics
61568ef3d20SLeila Ghaffari       // -- Density
61668ef3d20SLeila Ghaffari       v[0][i] -= wdetJb * rho * u_normal;
61768ef3d20SLeila Ghaffari 
61868ef3d20SLeila Ghaffari       // -- Momentum
61978e8b7daSJames Wright       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + normal[j] * P);
62068ef3d20SLeila Ghaffari 
62168ef3d20SLeila Ghaffari       // -- Total Energy Density
62268ef3d20SLeila Ghaffari       v[4][i] -= wdetJb * u_normal * (E + P);
62368ef3d20SLeila Ghaffari     }
624b193fadcSJames Wright   }
62568ef3d20SLeila Ghaffari   return 0;
62668ef3d20SLeila Ghaffari }
627