xref: /libCEED/examples/fluids/qfunctions/eulervortex.h (revision 5aed82e4fa97acf4ba24a7f10a35f5303a6798e0)
1*5aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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).
1477841947SLeila Ghaffari 
1577841947SLeila Ghaffari #ifndef eulervortex_h
1677841947SLeila Ghaffari #define eulervortex_h
1777841947SLeila Ghaffari 
1888b783a1SJames Wright #include <ceed.h>
19c9c2c079SJeremy L Thompson #include <math.h>
202b730f8bSJeremy L Thompson 
2113fa47b2SJames Wright #include "utils.h"
2277841947SLeila Ghaffari 
2377841947SLeila Ghaffari typedef struct EulerContext_ *EulerContext;
2477841947SLeila Ghaffari struct EulerContext_ {
2577841947SLeila Ghaffari   CeedScalar center[3];
2677841947SLeila Ghaffari   CeedScalar curr_time;
2777841947SLeila Ghaffari   CeedScalar vortex_strength;
28932417b3SJed Brown   CeedScalar c_tau;
2977841947SLeila Ghaffari   CeedScalar mean_velocity[3];
3077841947SLeila Ghaffari   bool       implicit;
31e6225c47SLeila Ghaffari   int        euler_test;
32e6225c47SLeila Ghaffari   int        stabilization;  // See StabilizationType: 0=none, 1=SU, 2=SUPG
3377841947SLeila Ghaffari };
3477841947SLeila Ghaffari 
3577841947SLeila Ghaffari // *****************************************************************************
3677841947SLeila Ghaffari // This function sets the initial conditions
3777841947SLeila Ghaffari //
3877841947SLeila Ghaffari //   Temperature:
3977841947SLeila Ghaffari //     T   = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
4077841947SLeila Ghaffari //   Density:
4177841947SLeila Ghaffari //     rho = (T/S_vortex)^(1 / (gamma - 1))
4277841947SLeila Ghaffari //   Pressure:
4377841947SLeila Ghaffari //     P   = rho * T
4477841947SLeila Ghaffari //   Velocity:
4577841947SLeila Ghaffari //     ui  = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
4677841947SLeila Ghaffari //     r   = sqrt( (x - xc)**2 + (y - yc)**2 )
4777841947SLeila Ghaffari //   Velocity/Momentum Density:
4877841947SLeila Ghaffari //     Ui  = rho ui
4977841947SLeila Ghaffari //   Total Energy:
5077841947SLeila Ghaffari //     E   = P / (gamma - 1) + rho (u u)/2
5177841947SLeila Ghaffari //
5277841947SLeila Ghaffari // Constants:
5377841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
5477841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
5577841947SLeila Ghaffari //   vortex_strength ,  Strength of vortex
5677841947SLeila Ghaffari //   center          ,  Location of bubble center
5777841947SLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
5877841947SLeila Ghaffari //
5977841947SLeila Ghaffari // *****************************************************************************
6077841947SLeila Ghaffari 
6177841947SLeila Ghaffari // *****************************************************************************
62ea61e9acSJeremy L Thompson // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling
63ea61e9acSJeremy L Thompson // vortex
6477841947SLeila Ghaffari // *****************************************************************************
652b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
6677841947SLeila Ghaffari   // Context
6777841947SLeila Ghaffari   const EulerContext context         = (EulerContext)ctx;
6877841947SLeila Ghaffari   const CeedScalar   vortex_strength = context->vortex_strength;
6977841947SLeila Ghaffari   const CeedScalar  *center          = context->center;  // Center of the domain
7077841947SLeila Ghaffari   const CeedScalar  *mean_velocity   = context->mean_velocity;
7177841947SLeila Ghaffari 
7277841947SLeila Ghaffari   // Setup
7377841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
7477841947SLeila Ghaffari   const CeedScalar cv    = 2.5;
7577841947SLeila Ghaffari   const CeedScalar R     = 1.;
7677841947SLeila Ghaffari   const CeedScalar x = X[0], y = X[1];  // Coordinates
7777841947SLeila Ghaffari   // Vortex center
7877841947SLeila Ghaffari   const CeedScalar xc = center[0] + mean_velocity[0] * time;
7977841947SLeila Ghaffari   const CeedScalar yc = center[1] + mean_velocity[1] * time;
8077841947SLeila Ghaffari 
8177841947SLeila Ghaffari   const CeedScalar x0       = x - xc;
8277841947SLeila Ghaffari   const CeedScalar y0       = y - yc;
8377841947SLeila Ghaffari   const CeedScalar r        = sqrt(x0 * x0 + y0 * y0);
8477841947SLeila Ghaffari   const CeedScalar C        = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI);
852b730f8bSJeremy L Thompson   const CeedScalar delta_T  = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI);
8677841947SLeila Ghaffari   const CeedScalar S_vortex = 1;  // no perturbation in the entropy P / rho^gamma
872b730f8bSJeremy L Thompson   const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI);
8877841947SLeila Ghaffari   CeedScalar       rho, P, T, E, u[3] = {0.};
8977841947SLeila Ghaffari 
9077841947SLeila Ghaffari   // Initial Conditions
9177841947SLeila Ghaffari   switch (context->euler_test) {
9277841947SLeila Ghaffari     case 0:  // Traveling vortex
9377841947SLeila Ghaffari       T = 1 + delta_T;
9477841947SLeila Ghaffari       // P = rho * T
9577841947SLeila Ghaffari       // P = S * rho^gamma
9677841947SLeila Ghaffari       // Solve for rho, then substitute for P
97e6225c47SLeila Ghaffari       rho  = pow(T / S_vortex, 1 / (gamma - 1.));
9877841947SLeila Ghaffari       P    = rho * T;
9977841947SLeila Ghaffari       u[0] = mean_velocity[0] - C * y0;
10077841947SLeila Ghaffari       u[1] = mean_velocity[1] + C * x0;
10177841947SLeila Ghaffari 
10277841947SLeila Ghaffari       // Assign exact solution
10377841947SLeila Ghaffari       q[0] = rho;
10477841947SLeila Ghaffari       q[1] = rho * u[0];
10577841947SLeila Ghaffari       q[2] = rho * u[1];
10677841947SLeila Ghaffari       q[3] = rho * u[2];
10777841947SLeila Ghaffari       q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.;
10877841947SLeila Ghaffari       break;
10977841947SLeila Ghaffari     case 1:  // Constant zero velocity, density constant, total energy constant
11077841947SLeila Ghaffari       rho = 1.;
11177841947SLeila Ghaffari       E   = 2.;
11277841947SLeila Ghaffari 
11377841947SLeila Ghaffari       // Assign exact solution
11477841947SLeila Ghaffari       q[0] = rho;
11577841947SLeila Ghaffari       q[1] = rho * u[0];
11677841947SLeila Ghaffari       q[2] = rho * u[1];
11777841947SLeila Ghaffari       q[3] = rho * u[2];
11877841947SLeila Ghaffari       q[4] = E;
11977841947SLeila Ghaffari       break;
12077841947SLeila Ghaffari     case 2:  // Constant nonzero velocity, density constant, total energy constant
12177841947SLeila Ghaffari       rho  = 1.;
12277841947SLeila Ghaffari       E    = 2.;
12377841947SLeila Ghaffari       u[0] = mean_velocity[0];
12477841947SLeila Ghaffari       u[1] = mean_velocity[1];
12577841947SLeila Ghaffari 
12677841947SLeila Ghaffari       // Assign exact solution
12777841947SLeila Ghaffari       q[0] = rho;
12877841947SLeila Ghaffari       q[1] = rho * u[0];
12977841947SLeila Ghaffari       q[2] = rho * u[1];
13077841947SLeila Ghaffari       q[3] = rho * u[2];
13177841947SLeila Ghaffari       q[4] = E;
13277841947SLeila Ghaffari       break;
133ea61e9acSJeremy 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
134ea61e9acSJeremy L Thompson              // bubble won't diffuse
13577841947SLeila Ghaffari       // (for Euler, where there is no thermal conductivity)
13677841947SLeila Ghaffari       P   = 1.;
13777841947SLeila Ghaffari       T   = 1. - S_bubble * exp(1. - r * r);
13877841947SLeila Ghaffari       rho = P / (R * T);
13977841947SLeila Ghaffari 
14077841947SLeila Ghaffari       // Assign exact solution
14177841947SLeila Ghaffari       q[0] = rho;
14277841947SLeila Ghaffari       q[1] = rho * u[0];
14377841947SLeila Ghaffari       q[2] = rho * u[1];
14477841947SLeila Ghaffari       q[3] = rho * u[2];
14577841947SLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
14677841947SLeila Ghaffari       break;
147ea61e9acSJeremy L Thompson     case 4:  // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant),
148ea61e9acSJeremy L Thompson       // It should be transported across the domain, but velocity stays constant
14977841947SLeila Ghaffari       P    = 1.;
15077841947SLeila Ghaffari       T    = 1. - S_bubble * exp(1. - r * r);
15177841947SLeila Ghaffari       rho  = P / (R * T);
15277841947SLeila Ghaffari       u[0] = mean_velocity[0];
15377841947SLeila Ghaffari       u[1] = mean_velocity[1];
15477841947SLeila Ghaffari 
15577841947SLeila Ghaffari       // Assign exact solution
15677841947SLeila Ghaffari       q[0] = rho;
15777841947SLeila Ghaffari       q[1] = rho * u[0];
15877841947SLeila Ghaffari       q[2] = rho * u[1];
15977841947SLeila Ghaffari       q[3] = rho * u[2];
16077841947SLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
16177841947SLeila Ghaffari       break;
16232f166c6SLeila Ghaffari     case 5:  // non-smooth thermal bubble - cylinder
16332f166c6SLeila Ghaffari       P    = 1.;
16432f166c6SLeila Ghaffari       T    = 1. - (r < 1. ? S_bubble : 0.);
16532f166c6SLeila Ghaffari       rho  = P / (R * T);
16632f166c6SLeila Ghaffari       u[0] = mean_velocity[0];
16732f166c6SLeila Ghaffari       u[1] = mean_velocity[1];
16832f166c6SLeila Ghaffari 
16932f166c6SLeila Ghaffari       // Assign exact solution
17032f166c6SLeila Ghaffari       q[0] = rho;
17132f166c6SLeila Ghaffari       q[1] = rho * u[0];
17232f166c6SLeila Ghaffari       q[2] = rho * u[1];
17332f166c6SLeila Ghaffari       q[3] = rho * u[2];
17432f166c6SLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
17532f166c6SLeila Ghaffari       break;
17677841947SLeila Ghaffari   }
17777841947SLeila Ghaffari   // Return
17877841947SLeila Ghaffari   return 0;
17977841947SLeila Ghaffari }
18077841947SLeila Ghaffari 
18177841947SLeila Ghaffari // *****************************************************************************
182e6225c47SLeila Ghaffari // Helper function for computing flux Jacobian
183e6225c47SLeila Ghaffari // *****************************************************************************
1842b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
185e6225c47SLeila Ghaffari                                                         const CeedScalar gamma) {
186e6225c47SLeila Ghaffari   CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2];  // Velocity square
187e6225c47SLeila Ghaffari   for (CeedInt i = 0; i < 3; i++) {                           // Jacobian matrices for 3 directions
188e6225c47SLeila Ghaffari     for (CeedInt j = 0; j < 3; j++) {                         // Rows of each Jacobian matrix
189e6225c47SLeila Ghaffari       dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j];
190e6225c47SLeila Ghaffari       for (CeedInt k = 0; k < 3; k++) {  // Columns of each Jacobian matrix
191e6225c47SLeila Ghaffari         dF[i][0][k + 1]     = ((i == k) ? 1. : 0.);
1922b730f8bSJeremy 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.);
1932b730f8bSJeremy L Thompson         dF[i][4][k + 1]     = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k];
194e6225c47SLeila Ghaffari       }
195e6225c47SLeila Ghaffari       dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.);
196e6225c47SLeila Ghaffari     }
197e6225c47SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho);
198e6225c47SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
199e6225c47SLeila Ghaffari   }
200e6225c47SLeila Ghaffari }
201e6225c47SLeila Ghaffari 
202e6225c47SLeila Ghaffari // *****************************************************************************
203932417b3SJed Brown // Helper function for computing Tau elements (stabilization constant)
204932417b3SJed Brown //   Model from:
205932417b3SJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
206932417b3SJed Brown //
207932417b3SJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
208932417b3SJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
209932417b3SJed Brown //
210932417b3SJed Brown // Where
211932417b3SJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
212932417b3SJed Brown //   h[i]      = 2 length(dxdX[i])
213932417b3SJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
214932417b3SJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
215ea61e9acSJeremy L Thompson //   rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i
216932417b3SJed Brown // *****************************************************************************
2172b730f8bSJeremy 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,
2182b730f8bSJeremy L Thompson                                        const CeedScalar c_tau) {
219ba6664aeSJames Wright   for (CeedInt i = 0; i < 3; i++) {
220932417b3SJed Brown     // length of element in direction i
2212b730f8bSJeremy L Thompson     CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]);
222932417b3SJed Brown     // fastest wave in direction i
223932417b3SJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
224932417b3SJed Brown     Tau_x[i]                = c_tau * h / fastest_wave;
225932417b3SJed Brown   }
226932417b3SJed Brown }
227932417b3SJed Brown 
228932417b3SJed Brown // *****************************************************************************
22977841947SLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
23077841947SLeila Ghaffari // *****************************************************************************
2312b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
23277841947SLeila Ghaffari   // Inputs
23377841947SLeila Ghaffari   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
23477841947SLeila Ghaffari 
23577841947SLeila Ghaffari   // Outputs
23677841947SLeila Ghaffari   CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
23777841947SLeila Ghaffari   const EulerContext context  = (EulerContext)ctx;
23877841947SLeila Ghaffari 
23977841947SLeila Ghaffari   // Quadrature Point Loop
24046603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
24177841947SLeila Ghaffari     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
242e6225c47SLeila Ghaffari     CeedScalar       q[5] = {0.};
24377841947SLeila Ghaffari 
24477841947SLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
24577841947SLeila Ghaffari 
2462b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
24777841947SLeila Ghaffari   }  // End of Quadrature Point Loop
24877841947SLeila Ghaffari 
24977841947SLeila Ghaffari   // Return
25077841947SLeila Ghaffari   return 0;
25177841947SLeila Ghaffari }
25277841947SLeila Ghaffari 
25377841947SLeila Ghaffari // *****************************************************************************
254ea61e9acSJeremy L Thompson // This QFunction implements the following formulation of Euler equations with explicit time stepping method
25577841947SLeila Ghaffari //
256ea61e9acSJeremy L Thompson // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density.
25777841947SLeila Ghaffari //
25877841947SLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
25977841947SLeila Ghaffari //   rho - Mass Density
26077841947SLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
26177841947SLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
26277841947SLeila Ghaffari //
26377841947SLeila Ghaffari // Euler Equations:
26477841947SLeila Ghaffari //   drho/dt + div( U )                   = 0
26577841947SLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
26677841947SLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
26777841947SLeila Ghaffari //
26877841947SLeila Ghaffari // Equation of State:
26977841947SLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
27077841947SLeila Ghaffari //
27177841947SLeila Ghaffari // Constants:
27277841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
27377841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
27477841947SLeila Ghaffari //   g               ,  Gravity
27577841947SLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
27677841947SLeila Ghaffari // *****************************************************************************
2772b730f8bSJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
27877841947SLeila Ghaffari   // Inputs
27946603fc5SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
28046603fc5SJames Wright   const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
281f3e15844SJames Wright   const CeedScalar(*q_data)            = in[2];
28246603fc5SJames Wright 
28377841947SLeila Ghaffari   // Outputs
28446603fc5SJames Wright   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
28546603fc5SJames Wright   CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
28677841947SLeila Ghaffari 
287e6225c47SLeila Ghaffari   EulerContext     context = (EulerContext)ctx;
288932417b3SJed Brown   const CeedScalar c_tau   = context->c_tau;
28977841947SLeila Ghaffari   const CeedScalar gamma   = 1.4;
29077841947SLeila Ghaffari 
29177841947SLeila Ghaffari   // Quadrature Point Loop
29246603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
29377841947SLeila Ghaffari     // Setup
29477841947SLeila Ghaffari     // -- Interp in
29577841947SLeila Ghaffari     const CeedScalar rho      = q[0][i];
2962b730f8bSJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
29777841947SLeila Ghaffari     const CeedScalar E        = q[4][i];
2982b730f8bSJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
2992b730f8bSJeremy L Thompson     const CeedScalar dU[3][3] = {
3002b730f8bSJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
3012b730f8bSJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
3022b730f8bSJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
303e6225c47SLeila Ghaffari     };
3042b730f8bSJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
305f3e15844SJames Wright     CeedScalar       wdetJ, dXdx[3][3];
306f3e15844SJames Wright     QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
307e6225c47SLeila Ghaffari     // dU/dx
308e6225c47SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
309e6225c47SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
310e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
311e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
312ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
313ba6664aeSJames Wright       for (CeedInt k = 0; k < 3; k++) {
314e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
315e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
316ba6664aeSJames Wright         for (CeedInt l = 0; l < 3; l++) {
317e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
318e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
319e6225c47SLeila Ghaffari         }
320e6225c47SLeila Ghaffari       }
321e6225c47SLeila Ghaffari     }
322e6225c47SLeila Ghaffari     // Pressure
3232b730f8bSJeremy 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,
324e6225c47SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
32577841947SLeila Ghaffari 
32677841947SLeila Ghaffari     // The Physics
32777841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
328ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) {
329e6225c47SLeila Ghaffari       v[j][i] = 0.;
3302b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
33177841947SLeila Ghaffari     }
33277841947SLeila Ghaffari 
33377841947SLeila Ghaffari     // -- Density
33477841947SLeila Ghaffari     // ---- u rho
3352b730f8bSJeremy 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]);
33677841947SLeila Ghaffari     // -- Momentum
33777841947SLeila Ghaffari     // ---- rho (u x u) + P I3
3382b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3392b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
3402b730f8bSJeremy 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] +
341e6225c47SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
3422b730f8bSJeremy L Thompson       }
3432b730f8bSJeremy L Thompson     }
34477841947SLeila Ghaffari     // -- Total Energy Density
34577841947SLeila Ghaffari     // ---- (E + P) u
3462b730f8bSJeremy 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]);
347e6225c47SLeila Ghaffari 
348e6225c47SLeila Ghaffari     // --Stabilization terms
349e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
350e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
351932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
352e6225c47SLeila Ghaffari 
353e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
354e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
355ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
356e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
357e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
3582b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
359e6225c47SLeila Ghaffari     }
360e6225c47SLeila Ghaffari 
361e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
362e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
3632b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3642b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
3652b730f8bSJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
3662b730f8bSJeremy L Thompson       }
3672b730f8bSJeremy L Thompson     }
368e6225c47SLeila Ghaffari 
369932417b3SJed Brown     // Stabilization
370932417b3SJed Brown     // -- Tau elements
371932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
372932417b3SJed Brown     CeedScalar       Tau_x[3]    = {0.};
373932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
374e6225c47SLeila Ghaffari 
375932417b3SJed Brown     // -- Stabilization method: none or SU
37688626eedSJames Wright     CeedScalar stab[5][3] = {{0.}};
377e6225c47SLeila Ghaffari     switch (context->stabilization) {
378e6225c47SLeila Ghaffari       case 0:  // Galerkin
379e6225c47SLeila Ghaffari         break;
380e6225c47SLeila Ghaffari       case 1:  // SU
3812b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
3822b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
3832b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
3842b730f8bSJeremy L Thompson           }
3852b730f8bSJeremy L Thompson         }
386e6225c47SLeila Ghaffari 
3872b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
3882b730f8bSJeremy 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]);
3892b730f8bSJeremy L Thompson         }
390e6225c47SLeila Ghaffari         break;
391e6225c47SLeila Ghaffari       case 2:  // SUPG is not implemented for explicit scheme
392e6225c47SLeila Ghaffari         break;
393e6225c47SLeila Ghaffari     }
394e6225c47SLeila Ghaffari 
39577841947SLeila Ghaffari   }  // End Quadrature Point Loop
39677841947SLeila Ghaffari 
39777841947SLeila Ghaffari   // Return
39877841947SLeila Ghaffari   return 0;
39977841947SLeila Ghaffari }
40077841947SLeila Ghaffari 
40177841947SLeila Ghaffari // *****************************************************************************
402ea61e9acSJeremy L Thompson // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method
40377841947SLeila Ghaffari // *****************************************************************************
4042b730f8bSJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
40577841947SLeila Ghaffari   // Inputs
40646603fc5SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
40746603fc5SJames Wright   const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
40846603fc5SJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
409f3e15844SJames Wright   const CeedScalar(*q_data)            = in[3];
41046603fc5SJames Wright 
41177841947SLeila Ghaffari   // Outputs
41246603fc5SJames Wright   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
41346603fc5SJames Wright   CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
41429ea4e10SJames Wright   CeedScalar *jac_data           = out[2];
41577841947SLeila Ghaffari 
416e6225c47SLeila Ghaffari   EulerContext     context   = (EulerContext)ctx;
417932417b3SJed Brown   const CeedScalar c_tau     = context->c_tau;
41877841947SLeila Ghaffari   const CeedScalar gamma     = 1.4;
41929ea4e10SJames Wright   const CeedScalar zeros[14] = {0.};
42077841947SLeila Ghaffari 
42177841947SLeila Ghaffari   // Quadrature Point Loop
42246603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
42377841947SLeila Ghaffari     // Setup
42477841947SLeila Ghaffari     // -- Interp in
42577841947SLeila Ghaffari     const CeedScalar rho      = q[0][i];
4262b730f8bSJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
42777841947SLeila Ghaffari     const CeedScalar E        = q[4][i];
4282b730f8bSJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
4292b730f8bSJeremy L Thompson     const CeedScalar dU[3][3] = {
4302b730f8bSJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
4312b730f8bSJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
4322b730f8bSJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
433e6225c47SLeila Ghaffari     };
4342b730f8bSJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
435f3e15844SJames Wright     CeedScalar       wdetJ, dXdx[3][3];
436f3e15844SJames Wright     QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
437e6225c47SLeila Ghaffari     // dU/dx
438e6225c47SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
439e6225c47SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
440e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
441e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
442ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
443ba6664aeSJames Wright       for (CeedInt k = 0; k < 3; k++) {
444e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
445e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
446ba6664aeSJames Wright         for (CeedInt l = 0; l < 3; l++) {
447e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
448e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
449e6225c47SLeila Ghaffari         }
450e6225c47SLeila Ghaffari       }
451e6225c47SLeila Ghaffari     }
4522b730f8bSJeremy 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,
453e6225c47SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
45477841947SLeila Ghaffari 
45577841947SLeila Ghaffari     // The Physics
45677841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
457ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) {
458e6225c47SLeila Ghaffari       v[j][i] = 0.;
4592b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
46077841947SLeila Ghaffari     }
46177841947SLeila Ghaffari     //-----mass matrix
4622b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i];
46377841947SLeila Ghaffari 
46477841947SLeila Ghaffari     // -- Density
46577841947SLeila Ghaffari     // ---- u rho
4662b730f8bSJeremy 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]);
46777841947SLeila Ghaffari     // -- Momentum
46877841947SLeila Ghaffari     // ---- rho (u x u) + P I3
4692b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
4702b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
4712b730f8bSJeremy 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] +
472e6225c47SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
4732b730f8bSJeremy L Thompson       }
4742b730f8bSJeremy L Thompson     }
47577841947SLeila Ghaffari     // -- Total Energy Density
47677841947SLeila Ghaffari     // ---- (E + P) u
4772b730f8bSJeremy 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]);
478e6225c47SLeila Ghaffari 
479e6225c47SLeila Ghaffari     // -- Stabilization terms
480e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
481e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
482932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
483e6225c47SLeila Ghaffari 
484e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
485e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
486ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
487e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
488e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
4892b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
490e6225c47SLeila Ghaffari     }
491e6225c47SLeila Ghaffari 
492e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
493e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
4942b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
4952b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
4962b730f8bSJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
4972b730f8bSJeremy L Thompson       }
4982b730f8bSJeremy L Thompson     }
499e6225c47SLeila Ghaffari 
500e6225c47SLeila Ghaffari     // ---- Strong residual
501e6225c47SLeila Ghaffari     CeedScalar strong_res[5];
5022b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j];
503e6225c47SLeila Ghaffari 
504932417b3SJed Brown     // Stabilization
505932417b3SJed Brown     // -- Tau elements
506932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
507932417b3SJed Brown     CeedScalar       Tau_x[3]    = {0.};
508932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
509e6225c47SLeila Ghaffari 
510932417b3SJed Brown     // -- Stabilization method: none, SU, or SUPG
51188626eedSJames Wright     CeedScalar stab[5][3] = {{0.}};
512e6225c47SLeila Ghaffari     switch (context->stabilization) {
513e6225c47SLeila Ghaffari       case 0:  // Galerkin
514e6225c47SLeila Ghaffari         break;
515e6225c47SLeila Ghaffari       case 1:  // SU
5162b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
5172b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
5182b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
5192b730f8bSJeremy L Thompson           }
5202b730f8bSJeremy L Thompson         }
521e6225c47SLeila Ghaffari 
5222b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
5232b730f8bSJeremy 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]);
5242b730f8bSJeremy L Thompson         }
525e6225c47SLeila Ghaffari         break;
526e6225c47SLeila Ghaffari       case 2:  // SUPG
5272b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
5282b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
5292b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
5302b730f8bSJeremy L Thompson           }
5312b730f8bSJeremy L Thompson         }
532e6225c47SLeila Ghaffari 
5332b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
5342b730f8bSJeremy 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]);
5352b730f8bSJeremy L Thompson         }
536e6225c47SLeila Ghaffari         break;
537e6225c47SLeila Ghaffari     }
53829ea4e10SJames Wright     StoredValuesPack(Q, i, 0, 14, zeros, jac_data);
53977841947SLeila Ghaffari   }  // End Quadrature Point Loop
54077841947SLeila Ghaffari 
54177841947SLeila Ghaffari   // Return
54277841947SLeila Ghaffari   return 0;
54377841947SLeila Ghaffari }
54477841947SLeila Ghaffari // *****************************************************************************
545ea61e9acSJeremy L Thompson // This QFunction sets the inflow boundary conditions for the traveling vortex problem.
54677841947SLeila Ghaffari //
547ea61e9acSJeremy L Thompson //  Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly.
54877841947SLeila Ghaffari // *****************************************************************************
5492b730f8bSJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
55077841947SLeila Ghaffari   // Inputs
551f3e15844SJames Wright   const CeedScalar(*q_data_sur) = in[2];
55277841947SLeila Ghaffari   // Outputs
55377841947SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
55477841947SLeila Ghaffari   EulerContext     context       = (EulerContext)ctx;
55577841947SLeila Ghaffari   const int        euler_test    = context->euler_test;
556f3e15844SJames Wright   const bool       is_implicit   = context->implicit;
55777841947SLeila Ghaffari   CeedScalar      *mean_velocity = context->mean_velocity;
55877841947SLeila Ghaffari   const CeedScalar cv            = 2.5;
55977841947SLeila Ghaffari   const CeedScalar R             = 1.;
56077841947SLeila Ghaffari   CeedScalar       T_inlet;
56177841947SLeila Ghaffari   CeedScalar       P_inlet;
56277841947SLeila Ghaffari 
56377841947SLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
5642b730f8bSJeremy L Thompson   if (euler_test == 1 || euler_test == 3) {
56577841947SLeila Ghaffari     for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.;
5662b730f8bSJeremy L Thompson   }
56777841947SLeila Ghaffari 
56877841947SLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
56977841947SLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
57077841947SLeila Ghaffari   else T_inlet = P_inlet = 1.;
57177841947SLeila Ghaffari 
57277841947SLeila Ghaffari   // Quadrature Point Loop
57346603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
574f3e15844SJames Wright     CeedScalar wdetJb, norm[3];
575f3e15844SJames Wright     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
576f3e15844SJames Wright     wdetJb *= is_implicit ? -1. : 1.;
57777841947SLeila Ghaffari 
57877841947SLeila Ghaffari     // face_normal = Normal vector of the face
5792b730f8bSJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
58077841947SLeila Ghaffari     // The Physics
58177841947SLeila Ghaffari     // Zero v so all future terms can safely sum into it
582ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
58377841947SLeila Ghaffari 
58477841947SLeila Ghaffari     // Implementing in/outflow BCs
5852fe7aee7SLeila Ghaffari     if (face_normal > 0) {
58677841947SLeila Ghaffari     } else {  // inflow
58777841947SLeila Ghaffari       const CeedScalar rho_inlet       = P_inlet / (R * T_inlet);
5882b730f8bSJeremy L Thompson       const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.;
58977841947SLeila Ghaffari       // incoming total energy
59077841947SLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
59177841947SLeila Ghaffari 
59277841947SLeila Ghaffari       // The Physics
59377841947SLeila Ghaffari       // -- Density
59477841947SLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
59577841947SLeila Ghaffari 
59677841947SLeila Ghaffari       // -- Momentum
5972b730f8bSJeremy 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);
59877841947SLeila Ghaffari 
59977841947SLeila Ghaffari       // -- Total Energy Density
60077841947SLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
60177841947SLeila Ghaffari     }
60277841947SLeila Ghaffari 
60377841947SLeila Ghaffari   }  // End Quadrature Point Loop
60477841947SLeila Ghaffari   return 0;
60577841947SLeila Ghaffari }
60677841947SLeila Ghaffari 
60777841947SLeila Ghaffari // *****************************************************************************
608ea61e9acSJeremy L Thompson // This QFunction sets the outflow boundary conditions for the Euler solver.
60955e76554SLeila Ghaffari //
61055e76554SLeila Ghaffari //  Outflow BCs:
611ea61e9acSJeremy L Thompson //    The validity of the weak form of the governing equations is extended to the outflow.
61255e76554SLeila Ghaffari // *****************************************************************************
6132b730f8bSJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
61455e76554SLeila Ghaffari   // Inputs
61546603fc5SJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
616f3e15844SJames Wright   const CeedScalar(*q_data_sur)    = in[2];
61746603fc5SJames Wright 
61855e76554SLeila Ghaffari   // Outputs
61955e76554SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
62055e76554SLeila Ghaffari   EulerContext context       = (EulerContext)ctx;
621f3e15844SJames Wright   const bool   is_implicit   = context->implicit;
62255e76554SLeila Ghaffari   CeedScalar  *mean_velocity = context->mean_velocity;
62355e76554SLeila Ghaffari 
62455e76554SLeila Ghaffari   const CeedScalar gamma = 1.4;
62555e76554SLeila Ghaffari 
62655e76554SLeila Ghaffari   // Quadrature Point Loop
62746603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
62855e76554SLeila Ghaffari     // Setup
62955e76554SLeila Ghaffari     // -- Interp in
63055e76554SLeila Ghaffari     const CeedScalar rho  = q[0][i];
6312b730f8bSJeremy L Thompson     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
63255e76554SLeila Ghaffari     const CeedScalar E    = q[4][i];
63355e76554SLeila Ghaffari 
634f3e15844SJames Wright     CeedScalar wdetJb, norm[3];
635f3e15844SJames Wright     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
636f3e15844SJames Wright     wdetJb *= is_implicit ? -1. : 1.;
63755e76554SLeila Ghaffari 
63855e76554SLeila Ghaffari     // face_normal = Normal vector of the face
6392b730f8bSJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
64055e76554SLeila Ghaffari     // The Physics
64155e76554SLeila Ghaffari     // Zero v so all future terms can safely sum into it
642ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0;
64355e76554SLeila Ghaffari 
64455e76554SLeila Ghaffari     // Implementing in/outflow BCs
64555e76554SLeila Ghaffari     if (face_normal > 0) {  // outflow
64655e76554SLeila Ghaffari       const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.;
64755e76554SLeila Ghaffari       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.);              // pressure
6482b730f8bSJeremy L Thompson       const CeedScalar u_normal  = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2];  // Normal velocity
64955e76554SLeila Ghaffari       // The Physics
65055e76554SLeila Ghaffari       // -- Density
65155e76554SLeila Ghaffari       v[0][i] -= wdetJb * rho * u_normal;
65255e76554SLeila Ghaffari 
65355e76554SLeila Ghaffari       // -- Momentum
6542b730f8bSJeremy L Thompson       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P);
65555e76554SLeila Ghaffari 
65655e76554SLeila Ghaffari       // -- Total Energy Density
65755e76554SLeila Ghaffari       v[4][i] -= wdetJb * u_normal * (E + P);
65855e76554SLeila Ghaffari     }
65955e76554SLeila Ghaffari   }  // End Quadrature Point Loop
66055e76554SLeila Ghaffari   return 0;
66155e76554SLeila Ghaffari }
66255e76554SLeila Ghaffari 
66355e76554SLeila Ghaffari // *****************************************************************************
66477841947SLeila Ghaffari 
66577841947SLeila Ghaffari #endif  // eulervortex_h
666