xref: /libCEED/examples/fluids/qfunctions/eulervortex.h (revision ea61e9ac44808524e4667c1525a05976f536c19c)
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:
13*ea61e9acSJeremy 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 // *****************************************************************************
62*ea61e9acSJeremy L Thompson // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling
63*ea61e9acSJeremy 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;
133*ea61e9acSJeremy 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
134*ea61e9acSJeremy 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;
147*ea61e9acSJeremy L Thompson     case 4:  // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant),
148*ea61e9acSJeremy 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 )
215*ea61e9acSJeremy 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 // *****************************************************************************
254*ea61e9acSJeremy L Thompson // This QFunction implements the following formulation of Euler equations with explicit time stepping method
25577841947SLeila Ghaffari //
256*ea61e9acSJeremy 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];
28146603fc5SJames Wright   const CeedScalar(*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])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]};
30577841947SLeila Ghaffari     // -- Interp-to-Interp q_data
30677841947SLeila Ghaffari     const CeedScalar wdetJ = q_data[0][i];
30777841947SLeila Ghaffari     // -- Interp-to-Grad q_data
30877841947SLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
3092b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
3102b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
3112b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
3122b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
31377841947SLeila Ghaffari     };
314e6225c47SLeila Ghaffari     // dU/dx
315e6225c47SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
316e6225c47SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
317e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
318e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
319ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
320ba6664aeSJames Wright       for (CeedInt k = 0; k < 3; k++) {
321e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
322e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
323ba6664aeSJames Wright         for (CeedInt l = 0; l < 3; l++) {
324e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
325e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
326e6225c47SLeila Ghaffari         }
327e6225c47SLeila Ghaffari       }
328e6225c47SLeila Ghaffari     }
329e6225c47SLeila Ghaffari     // Pressure
3302b730f8bSJeremy 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,
331e6225c47SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
33277841947SLeila Ghaffari 
33377841947SLeila Ghaffari     // The Physics
33477841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
335ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) {
336e6225c47SLeila Ghaffari       v[j][i] = 0.;
3372b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
33877841947SLeila Ghaffari     }
33977841947SLeila Ghaffari 
34077841947SLeila Ghaffari     // -- Density
34177841947SLeila Ghaffari     // ---- u rho
3422b730f8bSJeremy 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]);
34377841947SLeila Ghaffari     // -- Momentum
34477841947SLeila Ghaffari     // ---- rho (u x u) + P I3
3452b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3462b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
3472b730f8bSJeremy 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] +
348e6225c47SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
3492b730f8bSJeremy L Thompson       }
3502b730f8bSJeremy L Thompson     }
35177841947SLeila Ghaffari     // -- Total Energy Density
35277841947SLeila Ghaffari     // ---- (E + P) u
3532b730f8bSJeremy 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]);
354e6225c47SLeila Ghaffari 
355e6225c47SLeila Ghaffari     // --Stabilization terms
356e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
357e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
358932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
359e6225c47SLeila Ghaffari 
360e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
361e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
362ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
363e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
364e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
3652b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
366e6225c47SLeila Ghaffari     }
367e6225c47SLeila Ghaffari 
368e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
369e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
3702b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3712b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
3722b730f8bSJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
3732b730f8bSJeremy L Thompson       }
3742b730f8bSJeremy L Thompson     }
375e6225c47SLeila Ghaffari 
376932417b3SJed Brown     // Stabilization
377932417b3SJed Brown     // -- Tau elements
378932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
379932417b3SJed Brown     CeedScalar       Tau_x[3]    = {0.};
380932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
381e6225c47SLeila Ghaffari 
382932417b3SJed Brown     // -- Stabilization method: none or SU
38388626eedSJames Wright     CeedScalar stab[5][3] = {{0.}};
384e6225c47SLeila Ghaffari     switch (context->stabilization) {
385e6225c47SLeila Ghaffari       case 0:  // Galerkin
386e6225c47SLeila Ghaffari         break;
387e6225c47SLeila Ghaffari       case 1:  // SU
3882b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
3892b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
3902b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
3912b730f8bSJeremy L Thompson           }
3922b730f8bSJeremy L Thompson         }
393e6225c47SLeila Ghaffari 
3942b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
3952b730f8bSJeremy 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]);
3962b730f8bSJeremy L Thompson         }
397e6225c47SLeila Ghaffari         break;
398e6225c47SLeila Ghaffari       case 2:  // SUPG is not implemented for explicit scheme
399e6225c47SLeila Ghaffari         break;
400e6225c47SLeila Ghaffari     }
401e6225c47SLeila Ghaffari 
40277841947SLeila Ghaffari   }  // End Quadrature Point Loop
40377841947SLeila Ghaffari 
40477841947SLeila Ghaffari   // Return
40577841947SLeila Ghaffari   return 0;
40677841947SLeila Ghaffari }
40777841947SLeila Ghaffari 
40877841947SLeila Ghaffari // *****************************************************************************
409*ea61e9acSJeremy L Thompson // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method
41077841947SLeila Ghaffari // *****************************************************************************
4112b730f8bSJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
41277841947SLeila Ghaffari   // Inputs
41346603fc5SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]      = (const CeedScalar(*)[CEED_Q_VLA])in[0];
41446603fc5SJames Wright   const CeedScalar(*dq)[5][CEED_Q_VLA]  = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
41546603fc5SJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA]  = (const CeedScalar(*)[CEED_Q_VLA])in[2];
41646603fc5SJames Wright   const CeedScalar(*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
41746603fc5SJames Wright 
41877841947SLeila Ghaffari   // Outputs
41946603fc5SJames Wright   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
42046603fc5SJames Wright   CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
42177841947SLeila Ghaffari 
422e6225c47SLeila Ghaffari   EulerContext     context = (EulerContext)ctx;
423932417b3SJed Brown   const CeedScalar c_tau   = context->c_tau;
42477841947SLeila Ghaffari   const CeedScalar gamma   = 1.4;
42577841947SLeila Ghaffari 
42677841947SLeila Ghaffari   // Quadrature Point Loop
42746603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
42877841947SLeila Ghaffari     // Setup
42977841947SLeila Ghaffari     // -- Interp in
43077841947SLeila Ghaffari     const CeedScalar rho      = q[0][i];
4312b730f8bSJeremy L Thompson     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
43277841947SLeila Ghaffari     const CeedScalar E        = q[4][i];
4332b730f8bSJeremy L Thompson     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
4342b730f8bSJeremy L Thompson     const CeedScalar dU[3][3] = {
4352b730f8bSJeremy L Thompson         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
4362b730f8bSJeremy L Thompson         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
4372b730f8bSJeremy L Thompson         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
438e6225c47SLeila Ghaffari     };
4392b730f8bSJeremy L Thompson     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
44077841947SLeila Ghaffari     // -- Interp-to-Interp q_data
44177841947SLeila Ghaffari     const CeedScalar wdetJ = q_data[0][i];
44277841947SLeila Ghaffari     // -- Interp-to-Grad q_data
44377841947SLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
4442b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
4452b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
4462b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
4472b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
44877841947SLeila Ghaffari     };
449e6225c47SLeila Ghaffari     // dU/dx
450e6225c47SLeila Ghaffari     CeedScalar drhodx[3]       = {0.};
451e6225c47SLeila Ghaffari     CeedScalar dEdx[3]         = {0.};
452e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3]      = {{0.}};
453e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
454ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
455ba6664aeSJames Wright       for (CeedInt k = 0; k < 3; k++) {
456e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
457e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
458ba6664aeSJames Wright         for (CeedInt l = 0; l < 3; l++) {
459e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
460e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
461e6225c47SLeila Ghaffari         }
462e6225c47SLeila Ghaffari       }
463e6225c47SLeila Ghaffari     }
4642b730f8bSJeremy 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,
465e6225c47SLeila Ghaffari                      P = E_internal * (gamma - 1.);  // P = pressure
46677841947SLeila Ghaffari 
46777841947SLeila Ghaffari     // The Physics
46877841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
469ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) {
470e6225c47SLeila Ghaffari       v[j][i] = 0.;
4712b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
47277841947SLeila Ghaffari     }
47377841947SLeila Ghaffari     //-----mass matrix
4742b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i];
47577841947SLeila Ghaffari 
47677841947SLeila Ghaffari     // -- Density
47777841947SLeila Ghaffari     // ---- u rho
4782b730f8bSJeremy 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]);
47977841947SLeila Ghaffari     // -- Momentum
48077841947SLeila Ghaffari     // ---- rho (u x u) + P I3
4812b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
4822b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {
4832b730f8bSJeremy 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] +
484e6225c47SLeila Ghaffari                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
4852b730f8bSJeremy L Thompson       }
4862b730f8bSJeremy L Thompson     }
48777841947SLeila Ghaffari     // -- Total Energy Density
48877841947SLeila Ghaffari     // ---- (E + P) u
4892b730f8bSJeremy 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]);
490e6225c47SLeila Ghaffari 
491e6225c47SLeila Ghaffari     // -- Stabilization terms
492e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
493e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
494932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
495e6225c47SLeila Ghaffari 
496e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
497e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
498ba6664aeSJames Wright     for (CeedInt j = 0; j < 3; j++) {
499e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
500e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
5012b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
502e6225c47SLeila Ghaffari     }
503e6225c47SLeila Ghaffari 
504e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
505e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
5062b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
5072b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < 5; k++) {
5082b730f8bSJeremy L Thompson         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
5092b730f8bSJeremy L Thompson       }
5102b730f8bSJeremy L Thompson     }
511e6225c47SLeila Ghaffari 
512e6225c47SLeila Ghaffari     // ---- Strong residual
513e6225c47SLeila Ghaffari     CeedScalar strong_res[5];
5142b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j];
515e6225c47SLeila Ghaffari 
516932417b3SJed Brown     // Stabilization
517932417b3SJed Brown     // -- Tau elements
518932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
519932417b3SJed Brown     CeedScalar       Tau_x[3]    = {0.};
520932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
521e6225c47SLeila Ghaffari 
522932417b3SJed Brown     // -- Stabilization method: none, SU, or SUPG
52388626eedSJames Wright     CeedScalar stab[5][3] = {{0.}};
524e6225c47SLeila Ghaffari     switch (context->stabilization) {
525e6225c47SLeila Ghaffari       case 0:  // Galerkin
526e6225c47SLeila Ghaffari         break;
527e6225c47SLeila Ghaffari       case 1:  // SU
5282b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
5292b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
5302b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
5312b730f8bSJeremy L Thompson           }
5322b730f8bSJeremy L Thompson         }
533e6225c47SLeila Ghaffari 
5342b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
5352b730f8bSJeremy 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]);
5362b730f8bSJeremy L Thompson         }
537e6225c47SLeila Ghaffari         break;
538e6225c47SLeila Ghaffari       case 2:  // SUPG
5392b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 3; j++) {
5402b730f8bSJeremy L Thompson           for (CeedInt k = 0; k < 5; k++) {
5412b730f8bSJeremy L Thompson             for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
5422b730f8bSJeremy L Thompson           }
5432b730f8bSJeremy L Thompson         }
544e6225c47SLeila Ghaffari 
5452b730f8bSJeremy L Thompson         for (CeedInt j = 0; j < 5; j++) {
5462b730f8bSJeremy 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]);
5472b730f8bSJeremy L Thompson         }
548e6225c47SLeila Ghaffari         break;
549e6225c47SLeila Ghaffari     }
55077841947SLeila Ghaffari   }  // End Quadrature Point Loop
55177841947SLeila Ghaffari 
55277841947SLeila Ghaffari   // Return
55377841947SLeila Ghaffari   return 0;
55477841947SLeila Ghaffari }
55577841947SLeila Ghaffari // *****************************************************************************
556*ea61e9acSJeremy L Thompson // This QFunction sets the inflow boundary conditions for the traveling vortex problem.
55777841947SLeila Ghaffari //
558*ea61e9acSJeremy L Thompson //  Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly.
55977841947SLeila Ghaffari // *****************************************************************************
5602b730f8bSJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
56177841947SLeila Ghaffari   // Inputs
562e8b03feeSJames Wright   const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
56377841947SLeila Ghaffari   // Outputs
56477841947SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
56577841947SLeila Ghaffari   EulerContext     context       = (EulerContext)ctx;
56677841947SLeila Ghaffari   const int        euler_test    = context->euler_test;
56777841947SLeila Ghaffari   const bool       implicit      = context->implicit;
56877841947SLeila Ghaffari   CeedScalar      *mean_velocity = context->mean_velocity;
56977841947SLeila Ghaffari   const CeedScalar cv            = 2.5;
57077841947SLeila Ghaffari   const CeedScalar R             = 1.;
57177841947SLeila Ghaffari   CeedScalar       T_inlet;
57277841947SLeila Ghaffari   CeedScalar       P_inlet;
57377841947SLeila Ghaffari 
57477841947SLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
5752b730f8bSJeremy L Thompson   if (euler_test == 1 || euler_test == 3) {
57677841947SLeila Ghaffari     for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.;
5772b730f8bSJeremy L Thompson   }
57877841947SLeila Ghaffari 
57977841947SLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
58077841947SLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
58177841947SLeila Ghaffari   else T_inlet = P_inlet = 1.;
58277841947SLeila Ghaffari 
58377841947SLeila Ghaffari   // Quadrature Point Loop
58446603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
58577841947SLeila Ghaffari     // Setup
58677841947SLeila Ghaffari     // -- Interp-to-Interp q_data
58777841947SLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
58877841947SLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
58977841947SLeila Ghaffari     // We can effect this by swapping the sign on this weight
59077841947SLeila Ghaffari     const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i];
5912fe7aee7SLeila Ghaffari     // ---- Normal vect
5922b730f8bSJeremy L Thompson     const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
59377841947SLeila Ghaffari 
59477841947SLeila Ghaffari     // face_normal = Normal vector of the face
5952b730f8bSJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
59677841947SLeila Ghaffari     // The Physics
59777841947SLeila Ghaffari     // Zero v so all future terms can safely sum into it
598ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
59977841947SLeila Ghaffari 
60077841947SLeila Ghaffari     // Implementing in/outflow BCs
6012fe7aee7SLeila Ghaffari     if (face_normal > 0) {
60277841947SLeila Ghaffari     } else {  // inflow
60377841947SLeila Ghaffari       const CeedScalar rho_inlet       = P_inlet / (R * T_inlet);
6042b730f8bSJeremy L Thompson       const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.;
60577841947SLeila Ghaffari       // incoming total energy
60677841947SLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
60777841947SLeila Ghaffari 
60877841947SLeila Ghaffari       // The Physics
60977841947SLeila Ghaffari       // -- Density
61077841947SLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
61177841947SLeila Ghaffari 
61277841947SLeila Ghaffari       // -- Momentum
6132b730f8bSJeremy 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);
61477841947SLeila Ghaffari 
61577841947SLeila Ghaffari       // -- Total Energy Density
61677841947SLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
61777841947SLeila Ghaffari     }
61877841947SLeila Ghaffari 
61977841947SLeila Ghaffari   }  // End Quadrature Point Loop
62077841947SLeila Ghaffari   return 0;
62177841947SLeila Ghaffari }
62277841947SLeila Ghaffari 
62377841947SLeila Ghaffari // *****************************************************************************
624*ea61e9acSJeremy L Thompson // This QFunction sets the outflow boundary conditions for the Euler solver.
62555e76554SLeila Ghaffari //
62655e76554SLeila Ghaffari //  Outflow BCs:
627*ea61e9acSJeremy L Thompson //    The validity of the weak form of the governing equations is extended to the outflow.
62855e76554SLeila Ghaffari // *****************************************************************************
6292b730f8bSJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
63055e76554SLeila Ghaffari   // Inputs
63146603fc5SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA])in[0];
63246603fc5SJames Wright   const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
63346603fc5SJames Wright 
63455e76554SLeila Ghaffari   // Outputs
63555e76554SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
63655e76554SLeila Ghaffari   EulerContext context       = (EulerContext)ctx;
63755e76554SLeila Ghaffari   const bool   implicit      = context->implicit;
63855e76554SLeila Ghaffari   CeedScalar  *mean_velocity = context->mean_velocity;
63955e76554SLeila Ghaffari 
64055e76554SLeila Ghaffari   const CeedScalar gamma = 1.4;
64155e76554SLeila Ghaffari 
64255e76554SLeila Ghaffari   // Quadrature Point Loop
64346603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
64455e76554SLeila Ghaffari     // Setup
64555e76554SLeila Ghaffari     // -- Interp in
64655e76554SLeila Ghaffari     const CeedScalar rho  = q[0][i];
6472b730f8bSJeremy L Thompson     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
64855e76554SLeila Ghaffari     const CeedScalar E    = q[4][i];
64955e76554SLeila Ghaffari 
65055e76554SLeila Ghaffari     // -- Interp-to-Interp q_data
65155e76554SLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
65255e76554SLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
65355e76554SLeila Ghaffari     // We can effect this by swapping the sign on this weight
65455e76554SLeila Ghaffari     const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i];
65555e76554SLeila Ghaffari     // ---- Normal vectors
6562b730f8bSJeremy L Thompson     const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
65755e76554SLeila Ghaffari 
65855e76554SLeila Ghaffari     // face_normal = Normal vector of the face
6592b730f8bSJeremy L Thompson     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
66055e76554SLeila Ghaffari     // The Physics
66155e76554SLeila Ghaffari     // Zero v so all future terms can safely sum into it
662ba6664aeSJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0;
66355e76554SLeila Ghaffari 
66455e76554SLeila Ghaffari     // Implementing in/outflow BCs
66555e76554SLeila Ghaffari     if (face_normal > 0) {  // outflow
66655e76554SLeila Ghaffari       const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.;
66755e76554SLeila Ghaffari       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.);              // pressure
6682b730f8bSJeremy L Thompson       const CeedScalar u_normal  = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2];  // Normal velocity
66955e76554SLeila Ghaffari       // The Physics
67055e76554SLeila Ghaffari       // -- Density
67155e76554SLeila Ghaffari       v[0][i] -= wdetJb * rho * u_normal;
67255e76554SLeila Ghaffari 
67355e76554SLeila Ghaffari       // -- Momentum
6742b730f8bSJeremy L Thompson       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P);
67555e76554SLeila Ghaffari 
67655e76554SLeila Ghaffari       // -- Total Energy Density
67755e76554SLeila Ghaffari       v[4][i] -= wdetJb * u_normal * (E + P);
67855e76554SLeila Ghaffari     }
67955e76554SLeila Ghaffari   }  // End Quadrature Point Loop
68055e76554SLeila Ghaffari   return 0;
68155e76554SLeila Ghaffari }
68255e76554SLeila Ghaffari 
68355e76554SLeila Ghaffari // *****************************************************************************
68477841947SLeila Ghaffari 
68577841947SLeila Ghaffari #endif  // eulervortex_h
686