xref: /libCEED/examples/fluids/qfunctions/eulervortex.h (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Euler traveling vortex initial condition and operator for Navier-Stokes
1077841947SLeila Ghaffari /// example using PETSc
1177841947SLeila Ghaffari 
1277841947SLeila Ghaffari // Model from:
13ea61e9acSJeremy L Thompson //   On the Order of Accuracy and Numerical Performance of Two Classes of Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011).
14c0b5abf0SJeremy L Thompson #include <ceed/types.h>
15c0b5abf0SJeremy L Thompson #ifndef CEED_RUNNING_JIT_PASS
16c9c2c079SJeremy L Thompson #include <math.h>
17c0b5abf0SJeremy L Thompson #include <stdbool.h>
18c0b5abf0SJeremy L Thompson #endif
192b730f8bSJeremy L Thompson 
2013fa47b2SJames Wright #include "utils.h"
2177841947SLeila Ghaffari 
2277841947SLeila Ghaffari typedef struct EulerContext_ *EulerContext;
2377841947SLeila Ghaffari struct EulerContext_ {
2477841947SLeila Ghaffari   CeedScalar center[3];
2577841947SLeila Ghaffari   CeedScalar curr_time;
2677841947SLeila Ghaffari   CeedScalar vortex_strength;
27932417b3SJed Brown   CeedScalar c_tau;
2877841947SLeila Ghaffari   CeedScalar mean_velocity[3];
2977841947SLeila Ghaffari   bool       implicit;
30e6225c47SLeila Ghaffari   int        euler_test;
31e6225c47SLeila Ghaffari   int        stabilization;  // See StabilizationType: 0=none, 1=SU, 2=SUPG
3277841947SLeila Ghaffari };
3377841947SLeila Ghaffari 
3477841947SLeila Ghaffari // *****************************************************************************
3577841947SLeila Ghaffari // This function sets the initial conditions
3677841947SLeila Ghaffari //
3777841947SLeila Ghaffari //   Temperature:
3877841947SLeila Ghaffari //     T   = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
3977841947SLeila Ghaffari //   Density:
4077841947SLeila Ghaffari //     rho = (T/S_vortex)^(1 / (gamma - 1))
4177841947SLeila Ghaffari //   Pressure:
4277841947SLeila Ghaffari //     P   = rho * T
4377841947SLeila Ghaffari //   Velocity:
4477841947SLeila Ghaffari //     ui  = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
4577841947SLeila Ghaffari //     r   = sqrt( (x - xc)**2 + (y - yc)**2 )
4677841947SLeila Ghaffari //   Velocity/Momentum Density:
4777841947SLeila Ghaffari //     Ui  = rho ui
4877841947SLeila Ghaffari //   Total Energy:
4977841947SLeila Ghaffari //     E   = P / (gamma - 1) + rho (u u)/2
5077841947SLeila Ghaffari //
5177841947SLeila Ghaffari // Constants:
5277841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
5377841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
5477841947SLeila Ghaffari //   vortex_strength ,  Strength of vortex
5577841947SLeila Ghaffari //   center          ,  Location of bubble center
5677841947SLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
5777841947SLeila Ghaffari //
5877841947SLeila Ghaffari // *****************************************************************************
5977841947SLeila Ghaffari 
6077841947SLeila Ghaffari // *****************************************************************************
61ea61e9acSJeremy L Thompson // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling
62ea61e9acSJeremy L Thompson // vortex
6377841947SLeila Ghaffari // *****************************************************************************
Exact_Euler(CeedInt dim,CeedScalar time,const CeedScalar X[],CeedInt Nf,CeedScalar q[],void * ctx)642b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
6577841947SLeila Ghaffari   // Context
6677841947SLeila Ghaffari   const EulerContext context         = (EulerContext)ctx;
6777841947SLeila Ghaffari   const CeedScalar   vortex_strength = context->vortex_strength;
6877841947SLeila Ghaffari   const CeedScalar  *center          = context->center;  // Center of the domain
6977841947SLeila Ghaffari   const CeedScalar  *mean_velocity   = context->mean_velocity;
7077841947SLeila Ghaffari 
7177841947SLeila Ghaffari   // Setup
7277841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
7377841947SLeila Ghaffari   const CeedScalar cv    = 2.5;
7477841947SLeila Ghaffari   const CeedScalar R     = 1.;
7577841947SLeila Ghaffari   const CeedScalar x = X[0], y = X[1];  // Coordinates
7677841947SLeila Ghaffari   // Vortex center
7777841947SLeila Ghaffari   const CeedScalar xc = center[0] + mean_velocity[0] * time;
7877841947SLeila Ghaffari   const CeedScalar yc = center[1] + mean_velocity[1] * time;
7977841947SLeila Ghaffari 
8077841947SLeila Ghaffari   const CeedScalar x0       = x - xc;
8177841947SLeila Ghaffari   const CeedScalar y0       = y - yc;
8277841947SLeila Ghaffari   const CeedScalar r        = sqrt(x0 * x0 + y0 * y0);
8377841947SLeila Ghaffari   const CeedScalar C        = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI);
842b730f8bSJeremy L Thompson   const CeedScalar delta_T  = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI);
8577841947SLeila Ghaffari   const CeedScalar S_vortex = 1;  // no perturbation in the entropy P / rho^gamma
862b730f8bSJeremy L Thompson   const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI);
8777841947SLeila Ghaffari   CeedScalar       rho, P, T, E, u[3] = {0.};
8877841947SLeila Ghaffari 
8977841947SLeila Ghaffari   // Initial Conditions
9077841947SLeila Ghaffari   switch (context->euler_test) {
9177841947SLeila Ghaffari     case 0:  // Traveling vortex
9277841947SLeila Ghaffari       T = 1 + delta_T;
9377841947SLeila Ghaffari       // P = rho * T
9477841947SLeila Ghaffari       // P = S * rho^gamma
9577841947SLeila Ghaffari       // Solve for rho, then substitute for P
96e6225c47SLeila Ghaffari       rho  = pow(T / S_vortex, 1 / (gamma - 1.));
9777841947SLeila Ghaffari       P    = rho * T;
9877841947SLeila Ghaffari       u[0] = mean_velocity[0] - C * y0;
9977841947SLeila Ghaffari       u[1] = mean_velocity[1] + C * x0;
10077841947SLeila Ghaffari 
10177841947SLeila Ghaffari       // Assign exact solution
10277841947SLeila Ghaffari       q[0] = rho;
10377841947SLeila Ghaffari       q[1] = rho * u[0];
10477841947SLeila Ghaffari       q[2] = rho * u[1];
10577841947SLeila Ghaffari       q[3] = rho * u[2];
10677841947SLeila Ghaffari       q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.;
10777841947SLeila Ghaffari       break;
10877841947SLeila Ghaffari     case 1:  // Constant zero velocity, density constant, total energy constant
10977841947SLeila Ghaffari       rho = 1.;
11077841947SLeila Ghaffari       E   = 2.;
11177841947SLeila Ghaffari 
11277841947SLeila Ghaffari       // Assign exact solution
11377841947SLeila Ghaffari       q[0] = rho;
11477841947SLeila Ghaffari       q[1] = rho * u[0];
11577841947SLeila Ghaffari       q[2] = rho * u[1];
11677841947SLeila Ghaffari       q[3] = rho * u[2];
11777841947SLeila Ghaffari       q[4] = E;
11877841947SLeila Ghaffari       break;
11977841947SLeila Ghaffari     case 2:  // Constant nonzero velocity, density constant, total energy constant
12077841947SLeila Ghaffari       rho  = 1.;
12177841947SLeila Ghaffari       E    = 2.;
12277841947SLeila Ghaffari       u[0] = mean_velocity[0];
12377841947SLeila Ghaffari       u[1] = mean_velocity[1];
12477841947SLeila Ghaffari 
12577841947SLeila Ghaffari       // Assign exact solution
12677841947SLeila Ghaffari       q[0] = rho;
12777841947SLeila Ghaffari       q[1] = rho * u[0];
12877841947SLeila Ghaffari       q[2] = rho * u[1];
12977841947SLeila Ghaffari       q[3] = rho * u[2];
13077841947SLeila Ghaffari       q[4] = E;
13177841947SLeila Ghaffari       break;
132ea61e9acSJeremy 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
133ea61e9acSJeremy L Thompson              // bubble won't diffuse
13477841947SLeila Ghaffari       // (for Euler, where there is no thermal conductivity)
13577841947SLeila Ghaffari       P   = 1.;
13677841947SLeila Ghaffari       T   = 1. - S_bubble * exp(1. - r * r);
13777841947SLeila Ghaffari       rho = P / (R * T);
13877841947SLeila Ghaffari 
13977841947SLeila Ghaffari       // Assign exact solution
14077841947SLeila Ghaffari       q[0] = rho;
14177841947SLeila Ghaffari       q[1] = rho * u[0];
14277841947SLeila Ghaffari       q[2] = rho * u[1];
14377841947SLeila Ghaffari       q[3] = rho * u[2];
14477841947SLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
14577841947SLeila Ghaffari       break;
146ea61e9acSJeremy L Thompson     case 4:  // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant),
147ea61e9acSJeremy L Thompson       // It should be transported across the domain, but velocity stays constant
14877841947SLeila Ghaffari       P    = 1.;
14977841947SLeila Ghaffari       T    = 1. - S_bubble * exp(1. - r * r);
15077841947SLeila Ghaffari       rho  = P / (R * T);
15177841947SLeila Ghaffari       u[0] = mean_velocity[0];
15277841947SLeila Ghaffari       u[1] = mean_velocity[1];
15377841947SLeila Ghaffari 
15477841947SLeila Ghaffari       // Assign exact solution
15577841947SLeila Ghaffari       q[0] = rho;
15677841947SLeila Ghaffari       q[1] = rho * u[0];
15777841947SLeila Ghaffari       q[2] = rho * u[1];
15877841947SLeila Ghaffari       q[3] = rho * u[2];
15977841947SLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
16077841947SLeila Ghaffari       break;
16132f166c6SLeila Ghaffari     case 5:  // non-smooth thermal bubble - cylinder
16232f166c6SLeila Ghaffari       P    = 1.;
16332f166c6SLeila Ghaffari       T    = 1. - (r < 1. ? S_bubble : 0.);
16432f166c6SLeila Ghaffari       rho  = P / (R * T);
16532f166c6SLeila Ghaffari       u[0] = mean_velocity[0];
16632f166c6SLeila Ghaffari       u[1] = mean_velocity[1];
16732f166c6SLeila Ghaffari 
16832f166c6SLeila Ghaffari       // Assign exact solution
16932f166c6SLeila Ghaffari       q[0] = rho;
17032f166c6SLeila Ghaffari       q[1] = rho * u[0];
17132f166c6SLeila Ghaffari       q[2] = rho * u[1];
17232f166c6SLeila Ghaffari       q[3] = rho * u[2];
17332f166c6SLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
17432f166c6SLeila Ghaffari       break;
17577841947SLeila Ghaffari   }
17677841947SLeila Ghaffari   return 0;
17777841947SLeila Ghaffari }
17877841947SLeila Ghaffari 
17977841947SLeila Ghaffari // *****************************************************************************
180e6225c47SLeila Ghaffari // Helper function for computing flux Jacobian
181e6225c47SLeila Ghaffari // *****************************************************************************
ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5],const CeedScalar rho,const CeedScalar u[3],const CeedScalar E,const CeedScalar gamma)1822b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
183e6225c47SLeila Ghaffari                                                         const CeedScalar gamma) {
184e6225c47SLeila Ghaffari   CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2];  // Velocity square
185e6225c47SLeila Ghaffari   for (CeedInt i = 0; i < 3; i++) {                           // Jacobian matrices for 3 directions
186e6225c47SLeila Ghaffari     for (CeedInt j = 0; j < 3; j++) {                         // Rows of each Jacobian matrix
187e6225c47SLeila Ghaffari       dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j];
188e6225c47SLeila Ghaffari       for (CeedInt k = 0; k < 3; k++) {  // Columns of each Jacobian matrix
189e6225c47SLeila Ghaffari         dF[i][0][k + 1]     = ((i == k) ? 1. : 0.);
1902b730f8bSJeremy 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.);
1912b730f8bSJeremy L Thompson         dF[i][4][k + 1]     = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k];
192e6225c47SLeila Ghaffari       }
193e6225c47SLeila Ghaffari       dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.);
194e6225c47SLeila Ghaffari     }
195e6225c47SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho);
196e6225c47SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
197e6225c47SLeila Ghaffari   }
198e6225c47SLeila Ghaffari }
199e6225c47SLeila Ghaffari 
200e6225c47SLeila Ghaffari // *****************************************************************************
201932417b3SJed Brown // Helper function for computing Tau elements (stabilization constant)
202932417b3SJed Brown //   Model from:
203932417b3SJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
204932417b3SJed Brown //
205932417b3SJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
206932417b3SJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
207932417b3SJed Brown //
208932417b3SJed Brown // Where
209932417b3SJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
210932417b3SJed Brown //   h[i]      = 2 length(dxdX[i])
211932417b3SJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
212932417b3SJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
213ea61e9acSJeremy L Thompson //   rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i
214932417b3SJed Brown // *****************************************************************************
Tau_spatial(CeedScalar Tau_x[3],const CeedScalar dXdx[3][3],const CeedScalar u[3],const CeedScalar sound_speed,const CeedScalar c_tau)2152b730f8bSJeremy 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,
2162b730f8bSJeremy L Thompson                                        const CeedScalar c_tau) {
217ba6664aeSJames Wright   for (CeedInt i = 0; i < 3; i++) {
218932417b3SJed Brown     // length of element in direction i
2192b730f8bSJeremy L Thompson     CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]);
220932417b3SJed Brown     // fastest wave in direction i
221932417b3SJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
222932417b3SJed Brown     Tau_x[i]                = c_tau * h / fastest_wave;
223932417b3SJed Brown   }
224932417b3SJed Brown }
225932417b3SJed Brown 
226932417b3SJed Brown // *****************************************************************************
22777841947SLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
22877841947SLeila Ghaffari // *****************************************************************************
ICsEuler(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)2292b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
23077841947SLeila Ghaffari   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
23177841947SLeila Ghaffari   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
232f0b01153SJames Wright 
23377841947SLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
23477841947SLeila Ghaffari 
23546603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
23677841947SLeila Ghaffari     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
237e6225c47SLeila Ghaffari     CeedScalar       q[5] = {0.};
23877841947SLeila Ghaffari 
23977841947SLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
2402b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
241f0b01153SJames Wright   }
24277841947SLeila Ghaffari   return 0;
24377841947SLeila Ghaffari }
24477841947SLeila Ghaffari 
24577841947SLeila Ghaffari // *****************************************************************************
246ea61e9acSJeremy L Thompson // This QFunction implements the following formulation of Euler equations with explicit time stepping method
24777841947SLeila Ghaffari //
248ea61e9acSJeremy L Thompson // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density.
24977841947SLeila Ghaffari //
25077841947SLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
25177841947SLeila Ghaffari //   rho - Mass Density
25277841947SLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
25377841947SLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
25477841947SLeila Ghaffari //
25577841947SLeila Ghaffari // Euler Equations:
25677841947SLeila Ghaffari //   drho/dt + div( U )                   = 0
25777841947SLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
25877841947SLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
25977841947SLeila Ghaffari //
26077841947SLeila Ghaffari // Equation of State:
26177841947SLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
26277841947SLeila Ghaffari //
26377841947SLeila Ghaffari // Constants:
26477841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
26577841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
26677841947SLeila Ghaffari //   g               ,  Gravity
26777841947SLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
26877841947SLeila Ghaffari // *****************************************************************************
Euler(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)2692b730f8bSJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
27046603fc5SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
27146603fc5SJames Wright   const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
272f3e15844SJames Wright   const CeedScalar(*q_data)            = in[2];
27346603fc5SJames Wright   CeedScalar(*v)[CEED_Q_VLA]           = (CeedScalar(*)[CEED_Q_VLA])out[0];
27446603fc5SJames Wright   CeedScalar(*dv)[5][CEED_Q_VLA]       = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
27577841947SLeila Ghaffari 
276e6225c47SLeila Ghaffari   EulerContext     context = (EulerContext)ctx;
277932417b3SJed Brown   const CeedScalar c_tau   = context->c_tau;
27877841947SLeila Ghaffari   const CeedScalar gamma   = 1.4;
27977841947SLeila Ghaffari 
28046603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
28177841947SLeila Ghaffari     // Setup
28277841947SLeila Ghaffari     // -- Interp in
28377841947SLeila Ghaffari     const CeedScalar rho      = q[0][i];
2842b730f8bSJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
28577841947SLeila Ghaffari     const CeedScalar E        = q[4][i];
2862b730f8bSJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
2872b730f8bSJeremy L Thompson     const CeedScalar dU[3][3] = {
2882b730f8bSJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
2892b730f8bSJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
2902b730f8bSJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
291e6225c47SLeila Ghaffari     };
2922b730f8bSJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
293f3e15844SJames Wright     CeedScalar       wdetJ, dXdx[3][3];
294f3e15844SJames Wright     QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
295e6225c47SLeila Ghaffari     // dU/dx
296e6225c47SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
297e6225c47SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
298e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
299e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
300ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
301ba6664aeSJames Wright       for (CeedInt k = 0; k < 3; k++) {
302e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
303e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
304ba6664aeSJames Wright         for (CeedInt l = 0; l < 3; l++) {
305e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
306e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
307e6225c47SLeila Ghaffari         }
308e6225c47SLeila Ghaffari       }
309e6225c47SLeila Ghaffari     }
310e6225c47SLeila Ghaffari     // Pressure
3112b730f8bSJeremy 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,
312e6225c47SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
31377841947SLeila Ghaffari 
31477841947SLeila Ghaffari     // The Physics
31577841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
316ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) {
317e6225c47SLeila Ghaffari       v[j][i] = 0.;
3182b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
31977841947SLeila Ghaffari     }
32077841947SLeila Ghaffari 
32177841947SLeila Ghaffari     // -- Density
32277841947SLeila Ghaffari     // ---- u rho
3232b730f8bSJeremy 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]);
32477841947SLeila Ghaffari     // -- Momentum
32577841947SLeila Ghaffari     // ---- rho (u x u) + P I3
3262b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3272b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
3282b730f8bSJeremy 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] +
329e6225c47SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
3302b730f8bSJeremy L Thompson       }
3312b730f8bSJeremy L Thompson     }
33277841947SLeila Ghaffari     // -- Total Energy Density
33377841947SLeila Ghaffari     // ---- (E + P) u
3342b730f8bSJeremy 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]);
335e6225c47SLeila Ghaffari 
336e6225c47SLeila Ghaffari     // --Stabilization terms
337e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
338e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
339932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
340e6225c47SLeila Ghaffari 
341e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
342e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
343ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
344e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
345e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
3462b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
347e6225c47SLeila Ghaffari     }
348e6225c47SLeila Ghaffari 
349e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
350e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
3512b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3522b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
3532b730f8bSJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
3542b730f8bSJeremy L Thompson       }
3552b730f8bSJeremy L Thompson     }
356e6225c47SLeila Ghaffari 
357932417b3SJed Brown     // Stabilization
358932417b3SJed Brown     // -- Tau elements
359932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
360932417b3SJed Brown     CeedScalar       Tau_x[3]    = {0.};
361932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
362e6225c47SLeila Ghaffari 
363932417b3SJed Brown     // -- Stabilization method: none or SU
36488626eedSJames Wright     CeedScalar stab[5][3] = {{0.}};
365e6225c47SLeila Ghaffari     switch (context->stabilization) {
366e6225c47SLeila Ghaffari       case 0:  // Galerkin
367e6225c47SLeila Ghaffari         break;
368e6225c47SLeila Ghaffari       case 1:  // SU
3692b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
3702b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
3712b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
3722b730f8bSJeremy L Thompson           }
3732b730f8bSJeremy L Thompson         }
374e6225c47SLeila Ghaffari 
3752b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
3762b730f8bSJeremy 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]);
3772b730f8bSJeremy L Thompson         }
378e6225c47SLeila Ghaffari         break;
379e6225c47SLeila Ghaffari       case 2:  // SUPG is not implemented for explicit scheme
380e6225c47SLeila Ghaffari         break;
381e6225c47SLeila Ghaffari     }
382f0b01153SJames Wright   }
38377841947SLeila Ghaffari   return 0;
38477841947SLeila Ghaffari }
38577841947SLeila Ghaffari 
38677841947SLeila Ghaffari // *****************************************************************************
387ea61e9acSJeremy L Thompson // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method
38877841947SLeila Ghaffari // *****************************************************************************
IFunction_Euler(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)3892b730f8bSJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
39046603fc5SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
39146603fc5SJames Wright   const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
39246603fc5SJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
393f3e15844SJames Wright   const CeedScalar(*q_data)            = in[3];
39446603fc5SJames Wright   CeedScalar(*v)[CEED_Q_VLA]           = (CeedScalar(*)[CEED_Q_VLA])out[0];
39546603fc5SJames Wright   CeedScalar(*dv)[5][CEED_Q_VLA]       = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
39629ea4e10SJames Wright   CeedScalar *jac_data                 = out[2];
39777841947SLeila Ghaffari 
398e6225c47SLeila Ghaffari   EulerContext     context   = (EulerContext)ctx;
399932417b3SJed Brown   const CeedScalar c_tau     = context->c_tau;
40077841947SLeila Ghaffari   const CeedScalar gamma     = 1.4;
40129ea4e10SJames Wright   const CeedScalar zeros[14] = {0.};
40277841947SLeila Ghaffari 
40346603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
40477841947SLeila Ghaffari     // Setup
40577841947SLeila Ghaffari     // -- Interp in
40677841947SLeila Ghaffari     const CeedScalar rho      = q[0][i];
4072b730f8bSJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
40877841947SLeila Ghaffari     const CeedScalar E        = q[4][i];
4092b730f8bSJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
4102b730f8bSJeremy L Thompson     const CeedScalar dU[3][3] = {
4112b730f8bSJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
4122b730f8bSJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
4132b730f8bSJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
414e6225c47SLeila Ghaffari     };
4152b730f8bSJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
416f3e15844SJames Wright     CeedScalar       wdetJ, dXdx[3][3];
417f3e15844SJames Wright     QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
418e6225c47SLeila Ghaffari     // dU/dx
419e6225c47SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
420e6225c47SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
421e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
422e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
423ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
424ba6664aeSJames Wright       for (CeedInt k = 0; k < 3; k++) {
425e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
426e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
427ba6664aeSJames Wright         for (CeedInt l = 0; l < 3; l++) {
428e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
429e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
430e6225c47SLeila Ghaffari         }
431e6225c47SLeila Ghaffari       }
432e6225c47SLeila Ghaffari     }
4332b730f8bSJeremy 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,
434e6225c47SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
43577841947SLeila Ghaffari 
43677841947SLeila Ghaffari     // The Physics
43777841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
438ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) {
439e6225c47SLeila Ghaffari       v[j][i] = 0.;
4402b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
44177841947SLeila Ghaffari     }
44277841947SLeila Ghaffari     //-----mass matrix
4432b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i];
44477841947SLeila Ghaffari 
44577841947SLeila Ghaffari     // -- Density
44677841947SLeila Ghaffari     // ---- u rho
4472b730f8bSJeremy 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]);
44877841947SLeila Ghaffari     // -- Momentum
44977841947SLeila Ghaffari     // ---- rho (u x u) + P I3
4502b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
4512b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
4522b730f8bSJeremy 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] +
453e6225c47SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
4542b730f8bSJeremy L Thompson       }
4552b730f8bSJeremy L Thompson     }
45677841947SLeila Ghaffari     // -- Total Energy Density
45777841947SLeila Ghaffari     // ---- (E + P) u
4582b730f8bSJeremy 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]);
459e6225c47SLeila Ghaffari 
460e6225c47SLeila Ghaffari     // -- Stabilization terms
461e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
462e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
463932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
464e6225c47SLeila Ghaffari 
465e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
466e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
467ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
468e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
469e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
4702b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
471e6225c47SLeila Ghaffari     }
472e6225c47SLeila Ghaffari 
473e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
474e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
4752b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
4762b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
4772b730f8bSJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
4782b730f8bSJeremy L Thompson       }
4792b730f8bSJeremy L Thompson     }
480e6225c47SLeila Ghaffari 
481e6225c47SLeila Ghaffari     // ---- Strong residual
482e6225c47SLeila Ghaffari     CeedScalar strong_res[5];
4832b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j];
484e6225c47SLeila Ghaffari 
485932417b3SJed Brown     // Stabilization
486932417b3SJed Brown     // -- Tau elements
487932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
488932417b3SJed Brown     CeedScalar       Tau_x[3]    = {0.};
489932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
490e6225c47SLeila Ghaffari 
491932417b3SJed Brown     // -- Stabilization method: none, SU, or SUPG
49288626eedSJames Wright     CeedScalar stab[5][3] = {{0.}};
493e6225c47SLeila Ghaffari     switch (context->stabilization) {
494e6225c47SLeila Ghaffari       case 0:  // Galerkin
495e6225c47SLeila Ghaffari         break;
496e6225c47SLeila Ghaffari       case 1:  // SU
4972b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
4982b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
4992b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
5002b730f8bSJeremy L Thompson           }
5012b730f8bSJeremy L Thompson         }
502e6225c47SLeila Ghaffari 
5032b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
5042b730f8bSJeremy 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]);
5052b730f8bSJeremy L Thompson         }
506e6225c47SLeila Ghaffari         break;
507e6225c47SLeila Ghaffari       case 2:  // SUPG
5082b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
5092b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
5102b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
5112b730f8bSJeremy L Thompson           }
5122b730f8bSJeremy L Thompson         }
513e6225c47SLeila Ghaffari 
5142b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
5152b730f8bSJeremy 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]);
5162b730f8bSJeremy L Thompson         }
517e6225c47SLeila Ghaffari         break;
518e6225c47SLeila Ghaffari     }
51929ea4e10SJames Wright     StoredValuesPack(Q, i, 0, 14, zeros, jac_data);
520f0b01153SJames Wright   }
52177841947SLeila Ghaffari   return 0;
52277841947SLeila Ghaffari }
52377841947SLeila Ghaffari // *****************************************************************************
524ea61e9acSJeremy L Thompson // This QFunction sets the inflow boundary conditions for the traveling vortex problem.
52577841947SLeila Ghaffari //
526ea61e9acSJeremy L Thompson //  Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly.
52777841947SLeila Ghaffari // *****************************************************************************
TravelingVortex_Inflow(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)5282b730f8bSJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
529f3e15844SJames Wright   const CeedScalar(*q_data_sur) = in[2];
53077841947SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA]    = (CeedScalar(*)[CEED_Q_VLA])out[0];
531f0b01153SJames Wright 
53277841947SLeila Ghaffari   EulerContext     context       = (EulerContext)ctx;
53377841947SLeila Ghaffari   const int        euler_test    = context->euler_test;
534f3e15844SJames Wright   const bool       is_implicit   = context->implicit;
53577841947SLeila Ghaffari   CeedScalar      *mean_velocity = context->mean_velocity;
53677841947SLeila Ghaffari   const CeedScalar cv            = 2.5;
53777841947SLeila Ghaffari   const CeedScalar R             = 1.;
53877841947SLeila Ghaffari   CeedScalar       T_inlet;
53977841947SLeila Ghaffari   CeedScalar       P_inlet;
54077841947SLeila Ghaffari 
54177841947SLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
5422b730f8bSJeremy L Thompson   if (euler_test == 1 || euler_test == 3) {
54377841947SLeila Ghaffari     for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.;
5442b730f8bSJeremy L Thompson   }
54577841947SLeila Ghaffari 
54677841947SLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
54777841947SLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
54877841947SLeila Ghaffari   else T_inlet = P_inlet = 1.;
54977841947SLeila Ghaffari 
55046603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
551f3e15844SJames Wright     CeedScalar wdetJb, norm[3];
552f3e15844SJames Wright     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
553f3e15844SJames Wright     wdetJb *= is_implicit ? -1. : 1.;
55477841947SLeila Ghaffari 
55577841947SLeila Ghaffari     // face_normal = Normal vector of the face
5562b730f8bSJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
55777841947SLeila Ghaffari     // The Physics
55877841947SLeila Ghaffari     // Zero v so all future terms can safely sum into it
559ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
56077841947SLeila Ghaffari 
56177841947SLeila Ghaffari     // Implementing in/outflow BCs
5622fe7aee7SLeila Ghaffari     if (face_normal > 0) {
56377841947SLeila Ghaffari     } else {  // inflow
56477841947SLeila Ghaffari       const CeedScalar rho_inlet       = P_inlet / (R * T_inlet);
5652b730f8bSJeremy L Thompson       const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.;
56677841947SLeila Ghaffari       // incoming total energy
56777841947SLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
56877841947SLeila Ghaffari 
56977841947SLeila Ghaffari       // The Physics
57077841947SLeila Ghaffari       // -- Density
57177841947SLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
57277841947SLeila Ghaffari 
57377841947SLeila Ghaffari       // -- Momentum
5742b730f8bSJeremy 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);
57577841947SLeila Ghaffari 
57677841947SLeila Ghaffari       // -- Total Energy Density
57777841947SLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
57877841947SLeila Ghaffari     }
579f0b01153SJames Wright   }
58077841947SLeila Ghaffari   return 0;
58177841947SLeila Ghaffari }
58277841947SLeila Ghaffari 
58377841947SLeila Ghaffari // *****************************************************************************
584ea61e9acSJeremy L Thompson // This QFunction sets the outflow boundary conditions for the Euler solver.
58555e76554SLeila Ghaffari //
58655e76554SLeila Ghaffari //  Outflow BCs:
587ea61e9acSJeremy L Thompson //    The validity of the weak form of the governing equations is extended to the outflow.
58855e76554SLeila Ghaffari // *****************************************************************************
Euler_Outflow(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)5892b730f8bSJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
59046603fc5SJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
591f3e15844SJames Wright   const CeedScalar(*q_data_sur)    = in[2];
59255e76554SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA]       = (CeedScalar(*)[CEED_Q_VLA])out[0];
593f0b01153SJames Wright 
59455e76554SLeila Ghaffari   EulerContext context       = (EulerContext)ctx;
595f3e15844SJames Wright   const bool   is_implicit   = context->implicit;
59655e76554SLeila Ghaffari   CeedScalar  *mean_velocity = context->mean_velocity;
59755e76554SLeila Ghaffari 
59855e76554SLeila Ghaffari   const CeedScalar gamma = 1.4;
59955e76554SLeila Ghaffari 
60046603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
60155e76554SLeila Ghaffari     // Setup
60255e76554SLeila Ghaffari     // -- Interp in
60355e76554SLeila Ghaffari     const CeedScalar rho  = q[0][i];
6042b730f8bSJeremy L Thompson     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
60555e76554SLeila Ghaffari     const CeedScalar E    = q[4][i];
60655e76554SLeila Ghaffari 
607f3e15844SJames Wright     CeedScalar wdetJb, norm[3];
608f3e15844SJames Wright     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
609f3e15844SJames Wright     wdetJb *= is_implicit ? -1. : 1.;
61055e76554SLeila Ghaffari 
61155e76554SLeila Ghaffari     // face_normal = Normal vector of the face
6122b730f8bSJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
61355e76554SLeila Ghaffari     // The Physics
61455e76554SLeila Ghaffari     // Zero v so all future terms can safely sum into it
615ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0;
61655e76554SLeila Ghaffari 
61755e76554SLeila Ghaffari     // Implementing in/outflow BCs
61855e76554SLeila Ghaffari     if (face_normal > 0) {  // outflow
61955e76554SLeila Ghaffari       const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.;
62055e76554SLeila Ghaffari       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.);              // pressure
6212b730f8bSJeremy L Thompson       const CeedScalar u_normal  = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2];  // Normal velocity
62255e76554SLeila Ghaffari       // The Physics
62355e76554SLeila Ghaffari       // -- Density
62455e76554SLeila Ghaffari       v[0][i] -= wdetJb * rho * u_normal;
62555e76554SLeila Ghaffari 
62655e76554SLeila Ghaffari       // -- Momentum
6272b730f8bSJeremy L Thompson       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P);
62855e76554SLeila Ghaffari 
62955e76554SLeila Ghaffari       // -- Total Energy Density
63055e76554SLeila Ghaffari       v[4][i] -= wdetJb * u_normal * (E + P);
63155e76554SLeila Ghaffari     }
632f0b01153SJames Wright   }
63355e76554SLeila Ghaffari   return 0;
63455e76554SLeila Ghaffari }
635