xref: /libCEED/examples/fluids/qfunctions/eulervortex.h (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, 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:
1377841947SLeila Ghaffari //   On the Order of Accuracy and Numerical Performance of Two Classes of
1477841947SLeila Ghaffari //   Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011).
1577841947SLeila Ghaffari 
1677841947SLeila Ghaffari #ifndef eulervortex_h
1777841947SLeila Ghaffari #define eulervortex_h
1877841947SLeila Ghaffari 
1988b783a1SJames Wright #include <ceed.h>
20c9c2c079SJeremy L Thompson #include <math.h>
21*2b730f8bSJeremy L Thompson 
2213fa47b2SJames Wright #include "utils.h"
2377841947SLeila Ghaffari 
2477841947SLeila Ghaffari typedef struct EulerContext_ *EulerContext;
2577841947SLeila Ghaffari struct EulerContext_ {
2677841947SLeila Ghaffari   CeedScalar center[3];
2777841947SLeila Ghaffari   CeedScalar curr_time;
2877841947SLeila Ghaffari   CeedScalar vortex_strength;
29932417b3SJed Brown   CeedScalar c_tau;
3077841947SLeila Ghaffari   CeedScalar mean_velocity[3];
3177841947SLeila Ghaffari   bool       implicit;
32e6225c47SLeila Ghaffari   int        euler_test;
33e6225c47SLeila Ghaffari   int        stabilization;  // See StabilizationType: 0=none, 1=SU, 2=SUPG
3477841947SLeila Ghaffari };
3577841947SLeila Ghaffari 
3677841947SLeila Ghaffari // *****************************************************************************
3777841947SLeila Ghaffari // This function sets the initial conditions
3877841947SLeila Ghaffari //
3977841947SLeila Ghaffari //   Temperature:
4077841947SLeila Ghaffari //     T   = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
4177841947SLeila Ghaffari //   Density:
4277841947SLeila Ghaffari //     rho = (T/S_vortex)^(1 / (gamma - 1))
4377841947SLeila Ghaffari //   Pressure:
4477841947SLeila Ghaffari //     P   = rho * T
4577841947SLeila Ghaffari //   Velocity:
4677841947SLeila Ghaffari //     ui  = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
4777841947SLeila Ghaffari //     r   = sqrt( (x - xc)**2 + (y - yc)**2 )
4877841947SLeila Ghaffari //   Velocity/Momentum Density:
4977841947SLeila Ghaffari //     Ui  = rho ui
5077841947SLeila Ghaffari //   Total Energy:
5177841947SLeila Ghaffari //     E   = P / (gamma - 1) + rho (u u)/2
5277841947SLeila Ghaffari //
5377841947SLeila Ghaffari // Constants:
5477841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
5577841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
5677841947SLeila Ghaffari //   vortex_strength ,  Strength of vortex
5777841947SLeila Ghaffari //   center          ,  Location of bubble center
5877841947SLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
5977841947SLeila Ghaffari //
6077841947SLeila Ghaffari // *****************************************************************************
6177841947SLeila Ghaffari 
6277841947SLeila Ghaffari // *****************************************************************************
6377841947SLeila Ghaffari // This helper function provides support for the exact, time-dependent solution
6477841947SLeila Ghaffari //   (currently not implemented) and IC formulation for Euler traveling vortex
6577841947SLeila Ghaffari // *****************************************************************************
66*2b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
6777841947SLeila Ghaffari   // Context
6877841947SLeila Ghaffari   const EulerContext context         = (EulerContext)ctx;
6977841947SLeila Ghaffari   const CeedScalar   vortex_strength = context->vortex_strength;
7077841947SLeila Ghaffari   const CeedScalar  *center          = context->center;  // Center of the domain
7177841947SLeila Ghaffari   const CeedScalar  *mean_velocity   = context->mean_velocity;
7277841947SLeila Ghaffari 
7377841947SLeila Ghaffari   // Setup
7477841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
7577841947SLeila Ghaffari   const CeedScalar cv    = 2.5;
7677841947SLeila Ghaffari   const CeedScalar R     = 1.;
7777841947SLeila Ghaffari   const CeedScalar x = X[0], y = X[1];  // Coordinates
7877841947SLeila Ghaffari   // Vortex center
7977841947SLeila Ghaffari   const CeedScalar xc = center[0] + mean_velocity[0] * time;
8077841947SLeila Ghaffari   const CeedScalar yc = center[1] + mean_velocity[1] * time;
8177841947SLeila Ghaffari 
8277841947SLeila Ghaffari   const CeedScalar x0       = x - xc;
8377841947SLeila Ghaffari   const CeedScalar y0       = y - yc;
8477841947SLeila Ghaffari   const CeedScalar r        = sqrt(x0 * x0 + y0 * y0);
8577841947SLeila Ghaffari   const CeedScalar C        = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI);
86*2b730f8bSJeremy L Thompson   const CeedScalar delta_T  = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI);
8777841947SLeila Ghaffari   const CeedScalar S_vortex = 1;  // no perturbation in the entropy P / rho^gamma
88*2b730f8bSJeremy L Thompson   const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI);
8977841947SLeila Ghaffari   CeedScalar       rho, P, T, E, u[3] = {0.};
9077841947SLeila Ghaffari 
9177841947SLeila Ghaffari   // Initial Conditions
9277841947SLeila Ghaffari   switch (context->euler_test) {
9377841947SLeila Ghaffari     case 0:  // Traveling vortex
9477841947SLeila Ghaffari       T = 1 + delta_T;
9577841947SLeila Ghaffari       // P = rho * T
9677841947SLeila Ghaffari       // P = S * rho^gamma
9777841947SLeila Ghaffari       // Solve for rho, then substitute for P
98e6225c47SLeila Ghaffari       rho  = pow(T / S_vortex, 1 / (gamma - 1.));
9977841947SLeila Ghaffari       P    = rho * T;
10077841947SLeila Ghaffari       u[0] = mean_velocity[0] - C * y0;
10177841947SLeila Ghaffari       u[1] = mean_velocity[1] + C * x0;
10277841947SLeila Ghaffari 
10377841947SLeila Ghaffari       // Assign exact solution
10477841947SLeila Ghaffari       q[0] = rho;
10577841947SLeila Ghaffari       q[1] = rho * u[0];
10677841947SLeila Ghaffari       q[2] = rho * u[1];
10777841947SLeila Ghaffari       q[3] = rho * u[2];
10877841947SLeila Ghaffari       q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.;
10977841947SLeila Ghaffari       break;
11077841947SLeila Ghaffari     case 1:  // Constant zero velocity, density constant, total energy constant
11177841947SLeila Ghaffari       rho = 1.;
11277841947SLeila Ghaffari       E   = 2.;
11377841947SLeila Ghaffari 
11477841947SLeila Ghaffari       // Assign exact solution
11577841947SLeila Ghaffari       q[0] = rho;
11677841947SLeila Ghaffari       q[1] = rho * u[0];
11777841947SLeila Ghaffari       q[2] = rho * u[1];
11877841947SLeila Ghaffari       q[3] = rho * u[2];
11977841947SLeila Ghaffari       q[4] = E;
12077841947SLeila Ghaffari       break;
12177841947SLeila Ghaffari     case 2:  // Constant nonzero velocity, density constant, total energy constant
12277841947SLeila Ghaffari       rho  = 1.;
12377841947SLeila Ghaffari       E    = 2.;
12477841947SLeila Ghaffari       u[0] = mean_velocity[0];
12577841947SLeila Ghaffari       u[1] = mean_velocity[1];
12677841947SLeila Ghaffari 
12777841947SLeila Ghaffari       // Assign exact solution
12877841947SLeila Ghaffari       q[0] = rho;
12977841947SLeila Ghaffari       q[1] = rho * u[0];
13077841947SLeila Ghaffari       q[2] = rho * u[1];
13177841947SLeila Ghaffari       q[3] = rho * u[2];
13277841947SLeila Ghaffari       q[4] = E;
13377841947SLeila Ghaffari       break;
13477841947SLeila Ghaffari     case 3:  // Velocity zero, pressure constant
13577841947SLeila Ghaffari       // (so density and internal energy will be non-constant),
13677841947SLeila Ghaffari       // but the velocity should stay zero and the bubble won't diffuse
13777841947SLeila Ghaffari       // (for Euler, where there is no thermal conductivity)
13877841947SLeila Ghaffari       P   = 1.;
13977841947SLeila Ghaffari       T   = 1. - S_bubble * exp(1. - r * r);
14077841947SLeila Ghaffari       rho = P / (R * T);
14177841947SLeila Ghaffari 
14277841947SLeila Ghaffari       // Assign exact solution
14377841947SLeila Ghaffari       q[0] = rho;
14477841947SLeila Ghaffari       q[1] = rho * u[0];
14577841947SLeila Ghaffari       q[2] = rho * u[1];
14677841947SLeila Ghaffari       q[3] = rho * u[2];
14777841947SLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
14877841947SLeila Ghaffari       break;
14977841947SLeila Ghaffari     case 4:  // Constant nonzero velocity, pressure constant
15077841947SLeila Ghaffari       // (so density and internal energy will be non-constant),
15177841947SLeila Ghaffari       // it should be transported across the domain, but velocity stays constant
15277841947SLeila Ghaffari       P    = 1.;
15377841947SLeila Ghaffari       T    = 1. - S_bubble * exp(1. - r * r);
15477841947SLeila Ghaffari       rho  = P / (R * T);
15577841947SLeila Ghaffari       u[0] = mean_velocity[0];
15677841947SLeila Ghaffari       u[1] = mean_velocity[1];
15777841947SLeila Ghaffari 
15877841947SLeila Ghaffari       // Assign exact solution
15977841947SLeila Ghaffari       q[0] = rho;
16077841947SLeila Ghaffari       q[1] = rho * u[0];
16177841947SLeila Ghaffari       q[2] = rho * u[1];
16277841947SLeila Ghaffari       q[3] = rho * u[2];
16377841947SLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
16477841947SLeila Ghaffari       break;
16532f166c6SLeila Ghaffari     case 5:  // non-smooth thermal bubble - cylinder
16632f166c6SLeila Ghaffari       P    = 1.;
16732f166c6SLeila Ghaffari       T    = 1. - (r < 1. ? S_bubble : 0.);
16832f166c6SLeila Ghaffari       rho  = P / (R * T);
16932f166c6SLeila Ghaffari       u[0] = mean_velocity[0];
17032f166c6SLeila Ghaffari       u[1] = mean_velocity[1];
17132f166c6SLeila Ghaffari 
17232f166c6SLeila Ghaffari       // Assign exact solution
17332f166c6SLeila Ghaffari       q[0] = rho;
17432f166c6SLeila Ghaffari       q[1] = rho * u[0];
17532f166c6SLeila Ghaffari       q[2] = rho * u[1];
17632f166c6SLeila Ghaffari       q[3] = rho * u[2];
17732f166c6SLeila Ghaffari       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
17832f166c6SLeila Ghaffari       break;
17977841947SLeila Ghaffari   }
18077841947SLeila Ghaffari   // Return
18177841947SLeila Ghaffari   return 0;
18277841947SLeila Ghaffari }
18377841947SLeila Ghaffari 
18477841947SLeila Ghaffari // *****************************************************************************
185e6225c47SLeila Ghaffari // Helper function for computing flux Jacobian
186e6225c47SLeila Ghaffari // *****************************************************************************
187*2b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
188e6225c47SLeila Ghaffari                                                         const CeedScalar gamma) {
189e6225c47SLeila Ghaffari   CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2];  // Velocity square
190e6225c47SLeila Ghaffari   for (CeedInt i = 0; i < 3; i++) {                           // Jacobian matrices for 3 directions
191e6225c47SLeila Ghaffari     for (CeedInt j = 0; j < 3; j++) {                         // Rows of each Jacobian matrix
192e6225c47SLeila Ghaffari       dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j];
193e6225c47SLeila Ghaffari       for (CeedInt k = 0; k < 3; k++) {  // Columns of each Jacobian matrix
194e6225c47SLeila Ghaffari         dF[i][0][k + 1]     = ((i == k) ? 1. : 0.);
195*2b730f8bSJeremy 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.);
196*2b730f8bSJeremy L Thompson         dF[i][4][k + 1]     = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k];
197e6225c47SLeila Ghaffari       }
198e6225c47SLeila Ghaffari       dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.);
199e6225c47SLeila Ghaffari     }
200e6225c47SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho);
201e6225c47SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
202e6225c47SLeila Ghaffari   }
203e6225c47SLeila Ghaffari }
204e6225c47SLeila Ghaffari 
205e6225c47SLeila Ghaffari // *****************************************************************************
206932417b3SJed Brown // Helper function for computing Tau elements (stabilization constant)
207932417b3SJed Brown //   Model from:
208932417b3SJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
209932417b3SJed Brown //
210932417b3SJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
211932417b3SJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
212932417b3SJed Brown //
213932417b3SJed Brown // Where
214932417b3SJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
215932417b3SJed Brown //   h[i]      = 2 length(dxdX[i])
216932417b3SJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
217932417b3SJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
218932417b3SJed Brown //   rho(A[i]) = spectral radius of the convective flux Jacobian i,
219932417b3SJed Brown //               wave speed in direction i
220932417b3SJed Brown // *****************************************************************************
221*2b730f8bSJeremy 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,
222*2b730f8bSJeremy L Thompson                                        const CeedScalar c_tau) {
223ba6664aeSJames Wright   for (CeedInt i = 0; i < 3; i++) {
224932417b3SJed Brown     // length of element in direction i
225*2b730f8bSJeremy L Thompson     CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]);
226932417b3SJed Brown     // fastest wave in direction i
227932417b3SJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
228932417b3SJed Brown     Tau_x[i]                = c_tau * h / fastest_wave;
229932417b3SJed Brown   }
230932417b3SJed Brown }
231932417b3SJed Brown 
232932417b3SJed Brown // *****************************************************************************
23377841947SLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
23477841947SLeila Ghaffari // *****************************************************************************
235*2b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
23677841947SLeila Ghaffari   // Inputs
23777841947SLeila Ghaffari   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
23877841947SLeila Ghaffari 
23977841947SLeila Ghaffari   // Outputs
24077841947SLeila Ghaffari   CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
24177841947SLeila Ghaffari   const EulerContext context  = (EulerContext)ctx;
24277841947SLeila Ghaffari 
24377841947SLeila Ghaffari   CeedPragmaSIMD
24477841947SLeila Ghaffari       // Quadrature Point Loop
24577841947SLeila Ghaffari       for (CeedInt i = 0; i < Q; i++) {
24677841947SLeila Ghaffari     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
247e6225c47SLeila Ghaffari     CeedScalar       q[5] = {0.};
24877841947SLeila Ghaffari 
24977841947SLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
25077841947SLeila Ghaffari 
251*2b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
25277841947SLeila Ghaffari   }  // End of Quadrature Point Loop
25377841947SLeila Ghaffari 
25477841947SLeila Ghaffari   // Return
25577841947SLeila Ghaffari   return 0;
25677841947SLeila Ghaffari }
25777841947SLeila Ghaffari 
25877841947SLeila Ghaffari // *****************************************************************************
25977841947SLeila Ghaffari // This QFunction implements the following formulation of Euler equations
26077841947SLeila Ghaffari //   with explicit time stepping method
26177841947SLeila Ghaffari //
26277841947SLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation
26377841947SLeila Ghaffari //   form with state variables of density, momentum density, and total
26477841947SLeila Ghaffari //   energy density.
26577841947SLeila Ghaffari //
26677841947SLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
26777841947SLeila Ghaffari //   rho - Mass Density
26877841947SLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
26977841947SLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
27077841947SLeila Ghaffari //
27177841947SLeila Ghaffari // Euler Equations:
27277841947SLeila Ghaffari //   drho/dt + div( U )                   = 0
27377841947SLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
27477841947SLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
27577841947SLeila Ghaffari //
27677841947SLeila Ghaffari // Equation of State:
27777841947SLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
27877841947SLeila Ghaffari //
27977841947SLeila Ghaffari // Constants:
28077841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
28177841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
28277841947SLeila Ghaffari //   g               ,  Gravity
28377841947SLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
28477841947SLeila Ghaffari // *****************************************************************************
285*2b730f8bSJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
28677841947SLeila Ghaffari   // *INDENT-OFF*
28777841947SLeila Ghaffari   // Inputs
288*2b730f8bSJeremy L Thompson   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
28977841947SLeila Ghaffari         (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
29077841947SLeila Ghaffari   // Outputs
291*2b730f8bSJeremy L Thompson   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
29277841947SLeila Ghaffari 
293e6225c47SLeila Ghaffari   EulerContext     context = (EulerContext)ctx;
294932417b3SJed Brown   const CeedScalar c_tau   = context->c_tau;
29577841947SLeila Ghaffari   const CeedScalar gamma   = 1.4;
29677841947SLeila Ghaffari 
29777841947SLeila Ghaffari   CeedPragmaSIMD
29877841947SLeila Ghaffari       // Quadrature Point Loop
29977841947SLeila Ghaffari       for (CeedInt i = 0; i < Q; i++) {
30077841947SLeila Ghaffari     // *INDENT-OFF*
30177841947SLeila Ghaffari     // Setup
30277841947SLeila Ghaffari     // -- Interp in
30377841947SLeila Ghaffari     const CeedScalar rho      = q[0][i];
304*2b730f8bSJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
30577841947SLeila Ghaffari     const CeedScalar E        = q[4][i];
306*2b730f8bSJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
307*2b730f8bSJeremy L Thompson     const CeedScalar dU[3][3] = {
308*2b730f8bSJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
309*2b730f8bSJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
310*2b730f8bSJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
311e6225c47SLeila Ghaffari     };
312*2b730f8bSJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
31377841947SLeila Ghaffari     // -- Interp-to-Interp q_data
31477841947SLeila Ghaffari     const CeedScalar wdetJ = q_data[0][i];
31577841947SLeila Ghaffari     // -- Interp-to-Grad q_data
31677841947SLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
31777841947SLeila Ghaffari     // *INDENT-OFF*
318*2b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
319*2b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
320*2b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
321*2b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
32277841947SLeila Ghaffari     };
32377841947SLeila Ghaffari     // *INDENT-ON*
324e6225c47SLeila Ghaffari     // dU/dx
325e6225c47SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
326e6225c47SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
327e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
328e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
329ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
330ba6664aeSJames Wright       for (CeedInt k = 0; k < 3; k++) {
331e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
332e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
333ba6664aeSJames Wright         for (CeedInt l = 0; l < 3; l++) {
334e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
335e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
336e6225c47SLeila Ghaffari         }
337e6225c47SLeila Ghaffari       }
338e6225c47SLeila Ghaffari     }
339e6225c47SLeila Ghaffari     // Pressure
340*2b730f8bSJeremy 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,
341e6225c47SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
34277841947SLeila Ghaffari 
34377841947SLeila Ghaffari     // The Physics
34477841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
345ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) {
346e6225c47SLeila Ghaffari       v[j][i] = 0.;
347*2b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
34877841947SLeila Ghaffari     }
34977841947SLeila Ghaffari 
35077841947SLeila Ghaffari     // -- Density
35177841947SLeila Ghaffari     // ---- u rho
352*2b730f8bSJeremy 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]);
35377841947SLeila Ghaffari     // -- Momentum
35477841947SLeila Ghaffari     // ---- rho (u x u) + P I3
355*2b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
356*2b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
357*2b730f8bSJeremy 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] +
358e6225c47SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
359*2b730f8bSJeremy L Thompson       }
360*2b730f8bSJeremy L Thompson     }
36177841947SLeila Ghaffari     // -- Total Energy Density
36277841947SLeila Ghaffari     // ---- (E + P) u
363*2b730f8bSJeremy 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]);
364e6225c47SLeila Ghaffari 
365e6225c47SLeila Ghaffari     // --Stabilization terms
366e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
367e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
368932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
369e6225c47SLeila Ghaffari 
370e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
371e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
372ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
373e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
374e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
375*2b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
376e6225c47SLeila Ghaffari     }
377e6225c47SLeila Ghaffari 
378e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
379e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
380*2b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
381*2b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
382*2b730f8bSJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
383*2b730f8bSJeremy L Thompson       }
384*2b730f8bSJeremy L Thompson     }
385e6225c47SLeila Ghaffari 
386932417b3SJed Brown     // Stabilization
387932417b3SJed Brown     // -- Tau elements
388932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
389932417b3SJed Brown     CeedScalar       Tau_x[3]    = {0.};
390932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
391e6225c47SLeila Ghaffari 
392932417b3SJed Brown     // -- Stabilization method: none or SU
39388626eedSJames Wright     CeedScalar stab[5][3] = {{0.}};
394e6225c47SLeila Ghaffari     switch (context->stabilization) {
395e6225c47SLeila Ghaffari       case 0:  // Galerkin
396e6225c47SLeila Ghaffari         break;
397e6225c47SLeila Ghaffari       case 1:  // SU
398*2b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
399*2b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
400*2b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
401*2b730f8bSJeremy L Thompson           }
402*2b730f8bSJeremy L Thompson         }
403e6225c47SLeila Ghaffari 
404*2b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
405*2b730f8bSJeremy 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]);
406*2b730f8bSJeremy L Thompson         }
407e6225c47SLeila Ghaffari         break;
408e6225c47SLeila Ghaffari       case 2:  // SUPG is not implemented for explicit scheme
409e6225c47SLeila Ghaffari         break;
410e6225c47SLeila Ghaffari     }
411e6225c47SLeila Ghaffari 
41277841947SLeila Ghaffari   }  // End Quadrature Point Loop
41377841947SLeila Ghaffari 
41477841947SLeila Ghaffari   // Return
41577841947SLeila Ghaffari   return 0;
41677841947SLeila Ghaffari }
41777841947SLeila Ghaffari 
41877841947SLeila Ghaffari // *****************************************************************************
41977841947SLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above)
42077841947SLeila Ghaffari //   with implicit time stepping method
42177841947SLeila Ghaffari //
42277841947SLeila Ghaffari // *****************************************************************************
423*2b730f8bSJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
42477841947SLeila Ghaffari   // *INDENT-OFF*
42577841947SLeila Ghaffari   // Inputs
426*2b730f8bSJeremy L Thompson   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
427*2b730f8bSJeremy L Thompson         (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
42877841947SLeila Ghaffari   // Outputs
429*2b730f8bSJeremy L Thompson   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
43077841947SLeila Ghaffari 
431e6225c47SLeila Ghaffari   EulerContext     context = (EulerContext)ctx;
432932417b3SJed Brown   const CeedScalar c_tau   = context->c_tau;
43377841947SLeila Ghaffari   const CeedScalar gamma   = 1.4;
43477841947SLeila Ghaffari 
43577841947SLeila Ghaffari   CeedPragmaSIMD
43677841947SLeila Ghaffari       // Quadrature Point Loop
43777841947SLeila Ghaffari       for (CeedInt i = 0; i < Q; i++) {
43877841947SLeila Ghaffari     // *INDENT-OFF*
43977841947SLeila Ghaffari     // Setup
44077841947SLeila Ghaffari     // -- Interp in
44177841947SLeila Ghaffari     const CeedScalar rho      = q[0][i];
442*2b730f8bSJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
44377841947SLeila Ghaffari     const CeedScalar E        = q[4][i];
444*2b730f8bSJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
445*2b730f8bSJeremy L Thompson     const CeedScalar dU[3][3] = {
446*2b730f8bSJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
447*2b730f8bSJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
448*2b730f8bSJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
449e6225c47SLeila Ghaffari     };
450*2b730f8bSJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
45177841947SLeila Ghaffari     // -- Interp-to-Interp q_data
45277841947SLeila Ghaffari     const CeedScalar wdetJ = q_data[0][i];
45377841947SLeila Ghaffari     // -- Interp-to-Grad q_data
45477841947SLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
45577841947SLeila Ghaffari     // *INDENT-OFF*
456*2b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
457*2b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
458*2b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
459*2b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
46077841947SLeila Ghaffari     };
46177841947SLeila Ghaffari     // *INDENT-ON*
462e6225c47SLeila Ghaffari     // dU/dx
463e6225c47SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
464e6225c47SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
465e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
466e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
467ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
468ba6664aeSJames Wright       for (CeedInt k = 0; k < 3; k++) {
469e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
470e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
471ba6664aeSJames Wright         for (CeedInt l = 0; l < 3; l++) {
472e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
473e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
474e6225c47SLeila Ghaffari         }
475e6225c47SLeila Ghaffari       }
476e6225c47SLeila Ghaffari     }
477*2b730f8bSJeremy 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,
478e6225c47SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
47977841947SLeila Ghaffari 
48077841947SLeila Ghaffari     // The Physics
48177841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
482ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) {
483e6225c47SLeila Ghaffari       v[j][i] = 0.;
484*2b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
48577841947SLeila Ghaffari     }
48677841947SLeila Ghaffari     //-----mass matrix
487*2b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i];
48877841947SLeila Ghaffari 
48977841947SLeila Ghaffari     // -- Density
49077841947SLeila Ghaffari     // ---- u rho
491*2b730f8bSJeremy 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]);
49277841947SLeila Ghaffari     // -- Momentum
49377841947SLeila Ghaffari     // ---- rho (u x u) + P I3
494*2b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
495*2b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
496*2b730f8bSJeremy 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] +
497e6225c47SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
498*2b730f8bSJeremy L Thompson       }
499*2b730f8bSJeremy L Thompson     }
50077841947SLeila Ghaffari     // -- Total Energy Density
50177841947SLeila Ghaffari     // ---- (E + P) u
502*2b730f8bSJeremy 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]);
503e6225c47SLeila Ghaffari 
504e6225c47SLeila Ghaffari     // -- Stabilization terms
505e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
506e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
507932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
508e6225c47SLeila Ghaffari 
509e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
510e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
511ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
512e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
513e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
514*2b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
515e6225c47SLeila Ghaffari     }
516e6225c47SLeila Ghaffari 
517e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
518e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
519*2b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
520*2b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
521*2b730f8bSJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
522*2b730f8bSJeremy L Thompson       }
523*2b730f8bSJeremy L Thompson     }
524e6225c47SLeila Ghaffari 
525e6225c47SLeila Ghaffari     // ---- Strong residual
526e6225c47SLeila Ghaffari     CeedScalar strong_res[5];
527*2b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j];
528e6225c47SLeila Ghaffari 
529932417b3SJed Brown     // Stabilization
530932417b3SJed Brown     // -- Tau elements
531932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
532932417b3SJed Brown     CeedScalar       Tau_x[3]    = {0.};
533932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
534e6225c47SLeila Ghaffari 
535932417b3SJed Brown     // -- Stabilization method: none, SU, or SUPG
53688626eedSJames Wright     CeedScalar stab[5][3] = {{0.}};
537e6225c47SLeila Ghaffari     switch (context->stabilization) {
538e6225c47SLeila Ghaffari       case 0:  // Galerkin
539e6225c47SLeila Ghaffari         break;
540e6225c47SLeila Ghaffari       case 1:  // SU
541*2b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
542*2b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
543*2b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
544*2b730f8bSJeremy L Thompson           }
545*2b730f8bSJeremy L Thompson         }
546e6225c47SLeila Ghaffari 
547*2b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
548*2b730f8bSJeremy 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]);
549*2b730f8bSJeremy L Thompson         }
550e6225c47SLeila Ghaffari         break;
551e6225c47SLeila Ghaffari       case 2:  // SUPG
552*2b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
553*2b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
554*2b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
555*2b730f8bSJeremy L Thompson           }
556*2b730f8bSJeremy L Thompson         }
557e6225c47SLeila Ghaffari 
558*2b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
559*2b730f8bSJeremy 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]);
560*2b730f8bSJeremy L Thompson         }
561e6225c47SLeila Ghaffari         break;
562e6225c47SLeila Ghaffari     }
56377841947SLeila Ghaffari   }  // End Quadrature Point Loop
56477841947SLeila Ghaffari 
56577841947SLeila Ghaffari   // Return
56677841947SLeila Ghaffari   return 0;
56777841947SLeila Ghaffari }
56877841947SLeila Ghaffari // *****************************************************************************
5692fe7aee7SLeila Ghaffari // This QFunction sets the inflow boundary conditions for
5702fe7aee7SLeila Ghaffari //   the traveling vortex problem.
57177841947SLeila Ghaffari //
57277841947SLeila Ghaffari //  Prescribed T_inlet and P_inlet are converted to conservative variables
57377841947SLeila Ghaffari //      and applied weakly.
57477841947SLeila Ghaffari //
57577841947SLeila Ghaffari // *****************************************************************************
576*2b730f8bSJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
57777841947SLeila Ghaffari   // *INDENT-OFF*
57877841947SLeila Ghaffari   // Inputs
579e8b03feeSJames Wright   const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
58077841947SLeila Ghaffari   // Outputs
58177841947SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
58277841947SLeila Ghaffari   // *INDENT-ON*
58377841947SLeila Ghaffari   EulerContext     context       = (EulerContext)ctx;
58477841947SLeila Ghaffari   const int        euler_test    = context->euler_test;
58577841947SLeila Ghaffari   const bool       implicit      = context->implicit;
58677841947SLeila Ghaffari   CeedScalar      *mean_velocity = context->mean_velocity;
58777841947SLeila Ghaffari   const CeedScalar cv            = 2.5;
58877841947SLeila Ghaffari   const CeedScalar R             = 1.;
58977841947SLeila Ghaffari   CeedScalar       T_inlet;
59077841947SLeila Ghaffari   CeedScalar       P_inlet;
59177841947SLeila Ghaffari 
59277841947SLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
593*2b730f8bSJeremy L Thompson   if (euler_test == 1 || euler_test == 3) {
59477841947SLeila Ghaffari     for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.;
595*2b730f8bSJeremy L Thompson   }
59677841947SLeila Ghaffari 
59777841947SLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
59877841947SLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
59977841947SLeila Ghaffari   else T_inlet = P_inlet = 1.;
60077841947SLeila Ghaffari 
60177841947SLeila Ghaffari   CeedPragmaSIMD
60277841947SLeila Ghaffari       // Quadrature Point Loop
60377841947SLeila Ghaffari       for (CeedInt i = 0; i < Q; i++) {
60477841947SLeila Ghaffari     // Setup
60577841947SLeila Ghaffari     // -- Interp-to-Interp q_data
60677841947SLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
60777841947SLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
60877841947SLeila Ghaffari     // We can effect this by swapping the sign on this weight
60977841947SLeila Ghaffari     const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i];
6102fe7aee7SLeila Ghaffari     // ---- Normal vect
611*2b730f8bSJeremy L Thompson     const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
61277841947SLeila Ghaffari 
61377841947SLeila Ghaffari     // face_normal = Normal vector of the face
614*2b730f8bSJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
61577841947SLeila Ghaffari     // The Physics
61677841947SLeila Ghaffari     // Zero v so all future terms can safely sum into it
617ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
61877841947SLeila Ghaffari 
61977841947SLeila Ghaffari     // Implementing in/outflow BCs
6202fe7aee7SLeila Ghaffari     if (face_normal > 0) {
62177841947SLeila Ghaffari     } else {  // inflow
62277841947SLeila Ghaffari       const CeedScalar rho_inlet       = P_inlet / (R * T_inlet);
623*2b730f8bSJeremy L Thompson       const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.;
62477841947SLeila Ghaffari       // incoming total energy
62577841947SLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
62677841947SLeila Ghaffari 
62777841947SLeila Ghaffari       // The Physics
62877841947SLeila Ghaffari       // -- Density
62977841947SLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
63077841947SLeila Ghaffari 
63177841947SLeila Ghaffari       // -- Momentum
632*2b730f8bSJeremy 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);
63377841947SLeila Ghaffari 
63477841947SLeila Ghaffari       // -- Total Energy Density
63577841947SLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
63677841947SLeila Ghaffari     }
63777841947SLeila Ghaffari 
63877841947SLeila Ghaffari   }  // End Quadrature Point Loop
63977841947SLeila Ghaffari   return 0;
64077841947SLeila Ghaffari }
64177841947SLeila Ghaffari 
64277841947SLeila Ghaffari // *****************************************************************************
64355e76554SLeila Ghaffari // This QFunction sets the outflow boundary conditions for
64455e76554SLeila Ghaffari //   the Euler solver.
64555e76554SLeila Ghaffari //
64655e76554SLeila Ghaffari //  Outflow BCs:
64755e76554SLeila Ghaffari //    The validity of the weak form of the governing equations is
64855e76554SLeila Ghaffari //      extended to the outflow.
64955e76554SLeila Ghaffari //
65055e76554SLeila Ghaffari // *****************************************************************************
651*2b730f8bSJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
65255e76554SLeila Ghaffari   // *INDENT-OFF*
65355e76554SLeila Ghaffari   // Inputs
654*2b730f8bSJeremy L Thompson   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
65555e76554SLeila Ghaffari   // Outputs
65655e76554SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
65755e76554SLeila Ghaffari   // *INDENT-ON*
65855e76554SLeila Ghaffari   EulerContext context       = (EulerContext)ctx;
65955e76554SLeila Ghaffari   const bool   implicit      = context->implicit;
66055e76554SLeila Ghaffari   CeedScalar  *mean_velocity = context->mean_velocity;
66155e76554SLeila Ghaffari 
66255e76554SLeila Ghaffari   const CeedScalar gamma = 1.4;
66355e76554SLeila Ghaffari 
66455e76554SLeila Ghaffari   CeedPragmaSIMD
66555e76554SLeila Ghaffari       // Quadrature Point Loop
66655e76554SLeila Ghaffari       for (CeedInt i = 0; i < Q; i++) {
66755e76554SLeila Ghaffari     // Setup
66855e76554SLeila Ghaffari     // -- Interp in
66955e76554SLeila Ghaffari     const CeedScalar rho  = q[0][i];
670*2b730f8bSJeremy L Thompson     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
67155e76554SLeila Ghaffari     const CeedScalar E    = q[4][i];
67255e76554SLeila Ghaffari 
67355e76554SLeila Ghaffari     // -- Interp-to-Interp q_data
67455e76554SLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
67555e76554SLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
67655e76554SLeila Ghaffari     // We can effect this by swapping the sign on this weight
67755e76554SLeila Ghaffari     const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i];
67855e76554SLeila Ghaffari     // ---- Normal vectors
679*2b730f8bSJeremy L Thompson     const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
68055e76554SLeila Ghaffari 
68155e76554SLeila Ghaffari     // face_normal = Normal vector of the face
682*2b730f8bSJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
68355e76554SLeila Ghaffari     // The Physics
68455e76554SLeila Ghaffari     // Zero v so all future terms can safely sum into it
685ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0;
68655e76554SLeila Ghaffari 
68755e76554SLeila Ghaffari     // Implementing in/outflow BCs
68855e76554SLeila Ghaffari     if (face_normal > 0) {  // outflow
68955e76554SLeila Ghaffari       const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.;
69055e76554SLeila Ghaffari       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.);              // pressure
691*2b730f8bSJeremy L Thompson       const CeedScalar u_normal  = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2];  // Normal velocity
69255e76554SLeila Ghaffari       // The Physics
69355e76554SLeila Ghaffari       // -- Density
69455e76554SLeila Ghaffari       v[0][i] -= wdetJb * rho * u_normal;
69555e76554SLeila Ghaffari 
69655e76554SLeila Ghaffari       // -- Momentum
697*2b730f8bSJeremy L Thompson       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P);
69855e76554SLeila Ghaffari 
69955e76554SLeila Ghaffari       // -- Total Energy Density
70055e76554SLeila Ghaffari       v[4][i] -= wdetJb * u_normal * (E + P);
70155e76554SLeila Ghaffari     }
70255e76554SLeila Ghaffari   }  // End Quadrature Point Loop
70355e76554SLeila Ghaffari   return 0;
70455e76554SLeila Ghaffari }
70555e76554SLeila Ghaffari 
70655e76554SLeila Ghaffari // *****************************************************************************
70777841947SLeila Ghaffari 
70877841947SLeila Ghaffari #endif  // eulervortex_h
709