xref: /libCEED/examples/fluids/qfunctions/eulervortex.h (revision 3d8e882215d238700cdceb37404f76ca7fa24eaa)
1*3d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*3d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
4*3d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
6*3d8e8822SJeremy 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 
1977841947SLeila Ghaffari #include <math.h>
2088b783a1SJames Wright #include <ceed.h>
2177841947SLeila Ghaffari 
2277841947SLeila Ghaffari #ifndef M_PI
2377841947SLeila Ghaffari #define M_PI    3.14159265358979323846
2477841947SLeila Ghaffari #endif
2577841947SLeila Ghaffari 
2677841947SLeila Ghaffari #ifndef euler_context_struct
2777841947SLeila Ghaffari #define euler_context_struct
2877841947SLeila Ghaffari typedef struct EulerContext_ *EulerContext;
2977841947SLeila Ghaffari struct EulerContext_ {
3077841947SLeila Ghaffari   CeedScalar center[3];
3177841947SLeila Ghaffari   CeedScalar curr_time;
3277841947SLeila Ghaffari   CeedScalar vortex_strength;
33932417b3SJed Brown   CeedScalar c_tau;
3477841947SLeila Ghaffari   CeedScalar mean_velocity[3];
3577841947SLeila Ghaffari   bool implicit;
36e6225c47SLeila Ghaffari   int euler_test;
37e6225c47SLeila Ghaffari   int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG
3877841947SLeila Ghaffari };
3977841947SLeila Ghaffari #endif
4077841947SLeila Ghaffari 
4177841947SLeila Ghaffari // *****************************************************************************
4277841947SLeila Ghaffari // This function sets the initial conditions
4377841947SLeila Ghaffari //
4477841947SLeila Ghaffari //   Temperature:
4577841947SLeila Ghaffari //     T   = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
4677841947SLeila Ghaffari //   Density:
4777841947SLeila Ghaffari //     rho = (T/S_vortex)^(1 / (gamma - 1))
4877841947SLeila Ghaffari //   Pressure:
4977841947SLeila Ghaffari //     P   = rho * T
5077841947SLeila Ghaffari //   Velocity:
5177841947SLeila Ghaffari //     ui  = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
5277841947SLeila Ghaffari //     r   = sqrt( (x - xc)**2 + (y - yc)**2 )
5377841947SLeila Ghaffari //   Velocity/Momentum Density:
5477841947SLeila Ghaffari //     Ui  = rho ui
5577841947SLeila Ghaffari //   Total Energy:
5677841947SLeila Ghaffari //     E   = P / (gamma - 1) + rho (u u)/2
5777841947SLeila Ghaffari //
5877841947SLeila Ghaffari // Constants:
5977841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
6077841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
6177841947SLeila Ghaffari //   vortex_strength ,  Strength of vortex
6277841947SLeila Ghaffari //   center          ,  Location of bubble center
6377841947SLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
6477841947SLeila Ghaffari //
6577841947SLeila Ghaffari // *****************************************************************************
6677841947SLeila Ghaffari 
6777841947SLeila Ghaffari // *****************************************************************************
6877841947SLeila Ghaffari // This helper function provides support for the exact, time-dependent solution
6977841947SLeila Ghaffari //   (currently not implemented) and IC formulation for Euler traveling vortex
7077841947SLeila Ghaffari // *****************************************************************************
7177841947SLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time,
7277841947SLeila Ghaffari                                       const CeedScalar X[], CeedInt Nf, CeedScalar q[],
7377841947SLeila Ghaffari                                       void *ctx) {
7477841947SLeila Ghaffari   // Context
7577841947SLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
7677841947SLeila Ghaffari   const CeedScalar vortex_strength    = context->vortex_strength;
7777841947SLeila Ghaffari   const CeedScalar *center            = context->center; // Center of the domain
7877841947SLeila Ghaffari   const CeedScalar *mean_velocity = context->mean_velocity;
7977841947SLeila Ghaffari 
8077841947SLeila Ghaffari   // Setup
8177841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
8277841947SLeila Ghaffari   const CeedScalar cv    = 2.5;
8377841947SLeila Ghaffari   const CeedScalar R     = 1.;
8477841947SLeila Ghaffari   const CeedScalar x     = X[0], y = X[1]; // Coordinates
8577841947SLeila Ghaffari   // Vortex center
8677841947SLeila Ghaffari   const CeedScalar xc = center[0] + mean_velocity[0] * time;
8777841947SLeila Ghaffari   const CeedScalar yc = center[1] + mean_velocity[1] * time;
8877841947SLeila Ghaffari 
8977841947SLeila Ghaffari   const CeedScalar x0       = x - xc;
9077841947SLeila Ghaffari   const CeedScalar y0       = y - yc;
9177841947SLeila Ghaffari   const CeedScalar r        = sqrt( x0*x0 + y0*y0 );
9277841947SLeila Ghaffari   const CeedScalar C        = vortex_strength * exp((1. - r*r)/2.) / (2. * M_PI);
93e6225c47SLeila Ghaffari   const CeedScalar delta_T  = - (gamma - 1.) * vortex_strength * vortex_strength *
94e6225c47SLeila Ghaffari                               exp(1 - r*r) / (8. * gamma * M_PI * M_PI);
9577841947SLeila Ghaffari   const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma
9677841947SLeila Ghaffari   const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength /
9777841947SLeila Ghaffari                               (8.*gamma*M_PI*M_PI);
9877841947SLeila Ghaffari   CeedScalar rho, P, T, E, u[3] = {0.};
9977841947SLeila Ghaffari 
10077841947SLeila Ghaffari   // Initial Conditions
10177841947SLeila Ghaffari   switch (context->euler_test) {
10277841947SLeila Ghaffari   case 0: // Traveling vortex
10377841947SLeila Ghaffari     T = 1 + delta_T;
10477841947SLeila Ghaffari     // P = rho * T
10577841947SLeila Ghaffari     // P = S * rho^gamma
10677841947SLeila Ghaffari     // Solve for rho, then substitute for P
107e6225c47SLeila Ghaffari     rho  = pow(T/S_vortex, 1 / (gamma - 1.));
10877841947SLeila Ghaffari     P    = rho * T;
10977841947SLeila Ghaffari     u[0] = mean_velocity[0] - C*y0;
11077841947SLeila Ghaffari     u[1] = mean_velocity[1] + C*x0;
11177841947SLeila Ghaffari 
11277841947SLeila Ghaffari     // Assign exact solution
11377841947SLeila Ghaffari     q[0] = rho;
11477841947SLeila Ghaffari     q[1] = rho * u[0];
11577841947SLeila Ghaffari     q[2] = rho * u[1];
11677841947SLeila Ghaffari     q[3] = rho * u[2];
11777841947SLeila Ghaffari     q[4] = P / (gamma - 1.) + rho * (u[0]*u[0] + u[1]*u[1]) / 2.;
11877841947SLeila Ghaffari     break;
11977841947SLeila Ghaffari   case 1: // Constant zero velocity, density constant, total energy constant
12077841947SLeila Ghaffari     rho  = 1.;
12177841947SLeila Ghaffari     E    = 2.;
12277841947SLeila Ghaffari 
12377841947SLeila Ghaffari     // Assign exact solution
12477841947SLeila Ghaffari     q[0] = rho;
12577841947SLeila Ghaffari     q[1] = rho * u[0];
12677841947SLeila Ghaffari     q[2] = rho * u[1];
12777841947SLeila Ghaffari     q[3] = rho * u[2];
12877841947SLeila Ghaffari     q[4] = E;
12977841947SLeila Ghaffari     break;
13077841947SLeila Ghaffari   case 2: // Constant nonzero velocity, density constant, total energy constant
13177841947SLeila Ghaffari     rho  = 1.;
13277841947SLeila Ghaffari     E    = 2.;
13377841947SLeila Ghaffari     u[0] = mean_velocity[0];
13477841947SLeila Ghaffari     u[1] = mean_velocity[1];
13577841947SLeila Ghaffari 
13677841947SLeila Ghaffari     // Assign exact solution
13777841947SLeila Ghaffari     q[0] = rho;
13877841947SLeila Ghaffari     q[1] = rho * u[0];
13977841947SLeila Ghaffari     q[2] = rho * u[1];
14077841947SLeila Ghaffari     q[3] = rho * u[2];
14177841947SLeila Ghaffari     q[4] = E;
14277841947SLeila Ghaffari     break;
14377841947SLeila Ghaffari   case 3: // Velocity zero, pressure constant
14477841947SLeila Ghaffari     // (so density and internal energy will be non-constant),
14577841947SLeila Ghaffari     // but the velocity should stay zero and the bubble won't diffuse
14677841947SLeila Ghaffari     // (for Euler, where there is no thermal conductivity)
14777841947SLeila Ghaffari     P    = 1.;
14877841947SLeila Ghaffari     T    = 1. - S_bubble * exp(1. - r*r);
14977841947SLeila Ghaffari     rho  = P / (R*T);
15077841947SLeila Ghaffari 
15177841947SLeila Ghaffari     // Assign exact solution
15277841947SLeila Ghaffari     q[0] = rho;
15377841947SLeila Ghaffari     q[1] = rho * u[0];
15477841947SLeila Ghaffari     q[2] = rho * u[1];
15577841947SLeila Ghaffari     q[3] = rho * u[2];
15677841947SLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
15777841947SLeila Ghaffari     break;
15877841947SLeila Ghaffari   case 4: // Constant nonzero velocity, pressure constant
15977841947SLeila Ghaffari     // (so density and internal energy will be non-constant),
16077841947SLeila Ghaffari     // it should be transported across the domain, but velocity stays constant
16177841947SLeila Ghaffari     P    = 1.;
16277841947SLeila Ghaffari     T    = 1. - S_bubble * exp(1. - r*r);
16377841947SLeila Ghaffari     rho  = P / (R*T);
16477841947SLeila Ghaffari     u[0] = mean_velocity[0];
16577841947SLeila Ghaffari     u[1] = mean_velocity[1];
16677841947SLeila Ghaffari 
16777841947SLeila Ghaffari     // Assign exact solution
16877841947SLeila Ghaffari     q[0] = rho;
16977841947SLeila Ghaffari     q[1] = rho * u[0];
17077841947SLeila Ghaffari     q[2] = rho * u[1];
17177841947SLeila Ghaffari     q[3] = rho * u[2];
17277841947SLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
17377841947SLeila Ghaffari     break;
17432f166c6SLeila Ghaffari   case 5: // non-smooth thermal bubble - cylinder
17532f166c6SLeila Ghaffari     P    = 1.;
17632f166c6SLeila Ghaffari     T = 1. - (r < 1. ? S_bubble : 0.);
17732f166c6SLeila Ghaffari     rho  = P / (R*T);
17832f166c6SLeila Ghaffari     u[0] = mean_velocity[0];
17932f166c6SLeila Ghaffari     u[1] = mean_velocity[1];
18032f166c6SLeila Ghaffari 
18132f166c6SLeila Ghaffari     // Assign exact solution
18232f166c6SLeila Ghaffari     q[0] = rho;
18332f166c6SLeila Ghaffari     q[1] = rho * u[0];
18432f166c6SLeila Ghaffari     q[2] = rho * u[1];
18532f166c6SLeila Ghaffari     q[3] = rho * u[2];
18632f166c6SLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
18732f166c6SLeila Ghaffari     break;
18877841947SLeila Ghaffari   }
18977841947SLeila Ghaffari   // Return
19077841947SLeila Ghaffari   return 0;
19177841947SLeila Ghaffari }
19277841947SLeila Ghaffari 
19377841947SLeila Ghaffari // *****************************************************************************
194e6225c47SLeila Ghaffari // Helper function for computing flux Jacobian
195e6225c47SLeila Ghaffari // *****************************************************************************
196932417b3SJed Brown CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5],
197e6225c47SLeila Ghaffari     const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
198e6225c47SLeila Ghaffari     const CeedScalar gamma) {
199e6225c47SLeila Ghaffari   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square
200e6225c47SLeila Ghaffari   for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions
201e6225c47SLeila Ghaffari     for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix
202e6225c47SLeila Ghaffari       dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j];
203e6225c47SLeila Ghaffari       for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix
204e6225c47SLeila Ghaffari         dF[i][0][k+1]   = ((i==k) ? 1. : 0.);
205e6225c47SLeila Ghaffari         dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) +
206e6225c47SLeila Ghaffari                           ((i==k) ? u[j] : 0.) -
207e6225c47SLeila Ghaffari                           ((i==j) ? u[k] : 0.) * (gamma-1.);
208e6225c47SLeila Ghaffari         dF[i][4][k+1]   = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) -
209e6225c47SLeila Ghaffari                           (gamma-1.)*u[i]*u[k];
210e6225c47SLeila Ghaffari       }
211e6225c47SLeila Ghaffari       dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.);
212e6225c47SLeila Ghaffari     }
213e6225c47SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho);
214e6225c47SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
215e6225c47SLeila Ghaffari   }
216e6225c47SLeila Ghaffari }
217e6225c47SLeila Ghaffari 
218e6225c47SLeila Ghaffari // *****************************************************************************
219932417b3SJed Brown // Helper function for computing Tau elements (stabilization constant)
220932417b3SJed Brown //   Model from:
221932417b3SJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
222932417b3SJed Brown //
223932417b3SJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
224932417b3SJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
225932417b3SJed Brown //
226932417b3SJed Brown // Where
227932417b3SJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
228932417b3SJed Brown //   h[i]      = 2 length(dxdX[i])
229932417b3SJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
230932417b3SJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
231932417b3SJed Brown //   rho(A[i]) = spectral radius of the convective flux Jacobian i,
232932417b3SJed Brown //               wave speed in direction i
233932417b3SJed Brown // *****************************************************************************
234932417b3SJed Brown CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3],
235932417b3SJed Brown                                        const CeedScalar dXdx[3][3], const CeedScalar u[3],
236932417b3SJed Brown                                        const CeedScalar sound_speed, const CeedScalar c_tau) {
237932417b3SJed Brown   for (int i=0; i<3; i++) {
238932417b3SJed Brown     // length of element in direction i
239932417b3SJed Brown     CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] +
240932417b3SJed Brown                             dXdx[2][i]*dXdx[2][i]);
241932417b3SJed Brown     // fastest wave in direction i
242932417b3SJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
243932417b3SJed Brown     Tau_x[i] = c_tau * h / fastest_wave;
244932417b3SJed Brown   }
245932417b3SJed Brown }
246932417b3SJed Brown 
247932417b3SJed Brown // *****************************************************************************
24877841947SLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
24977841947SLeila Ghaffari // *****************************************************************************
25077841947SLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q,
25177841947SLeila Ghaffari                          const CeedScalar *const *in, CeedScalar *const *out) {
25277841947SLeila Ghaffari   // Inputs
25377841947SLeila Ghaffari   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
25477841947SLeila Ghaffari 
25577841947SLeila Ghaffari   // Outputs
25677841947SLeila Ghaffari   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
25777841947SLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
25877841947SLeila Ghaffari 
25977841947SLeila Ghaffari   CeedPragmaSIMD
26077841947SLeila Ghaffari   // Quadrature Point Loop
26177841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
26277841947SLeila Ghaffari     const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
263e6225c47SLeila Ghaffari     CeedScalar q[5] = {0.};
26477841947SLeila Ghaffari 
26577841947SLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
26677841947SLeila Ghaffari 
26777841947SLeila Ghaffari     for (CeedInt j=0; j<5; j++)
26877841947SLeila Ghaffari       q0[j][i] = q[j];
26977841947SLeila Ghaffari   } // End of Quadrature Point Loop
27077841947SLeila Ghaffari 
27177841947SLeila Ghaffari   // Return
27277841947SLeila Ghaffari   return 0;
27377841947SLeila Ghaffari }
27477841947SLeila Ghaffari 
27577841947SLeila Ghaffari // *****************************************************************************
27677841947SLeila Ghaffari // This QFunction implements the following formulation of Euler equations
27777841947SLeila Ghaffari //   with explicit time stepping method
27877841947SLeila Ghaffari //
27977841947SLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation
28077841947SLeila Ghaffari //   form with state variables of density, momentum density, and total
28177841947SLeila Ghaffari //   energy density.
28277841947SLeila Ghaffari //
28377841947SLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
28477841947SLeila Ghaffari //   rho - Mass Density
28577841947SLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
28677841947SLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
28777841947SLeila Ghaffari //
28877841947SLeila Ghaffari // Euler Equations:
28977841947SLeila Ghaffari //   drho/dt + div( U )                   = 0
29077841947SLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
29177841947SLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
29277841947SLeila Ghaffari //
29377841947SLeila Ghaffari // Equation of State:
29477841947SLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
29577841947SLeila Ghaffari //
29677841947SLeila Ghaffari // Constants:
29777841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
29877841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
29977841947SLeila Ghaffari //   g               ,  Gravity
30077841947SLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
30177841947SLeila Ghaffari // *****************************************************************************
30277841947SLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q,
30377841947SLeila Ghaffari                       const CeedScalar *const *in, CeedScalar *const *out) {
30477841947SLeila Ghaffari   // *INDENT-OFF*
30577841947SLeila Ghaffari   // Inputs
30677841947SLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
307e6225c47SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
30877841947SLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
30977841947SLeila Ghaffari   // Outputs
31077841947SLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
31177841947SLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
31277841947SLeila Ghaffari 
313e6225c47SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
314932417b3SJed Brown   const CeedScalar c_tau = context->c_tau;
31577841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
31677841947SLeila Ghaffari 
31777841947SLeila Ghaffari   CeedPragmaSIMD
31877841947SLeila Ghaffari   // Quadrature Point Loop
31977841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
32077841947SLeila Ghaffari     // *INDENT-OFF*
32177841947SLeila Ghaffari     // Setup
32277841947SLeila Ghaffari     // -- Interp in
32377841947SLeila Ghaffari     const CeedScalar rho        =   q[0][i];
32477841947SLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
32577841947SLeila Ghaffari                                     q[2][i] / rho,
32677841947SLeila Ghaffari                                     q[3][i] / rho
32777841947SLeila Ghaffari                                    };
32877841947SLeila Ghaffari     const CeedScalar E          =   q[4][i];
329e6225c47SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
330e6225c47SLeila Ghaffari                                     dq[1][0][i],
331e6225c47SLeila Ghaffari                                     dq[2][0][i]
332e6225c47SLeila Ghaffari                                    };
333e6225c47SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
334e6225c47SLeila Ghaffari                                     dq[1][1][i],
335e6225c47SLeila Ghaffari                                     dq[2][1][i]},
336e6225c47SLeila Ghaffari                                    {dq[0][2][i],
337e6225c47SLeila Ghaffari                                     dq[1][2][i],
338e6225c47SLeila Ghaffari                                     dq[2][2][i]},
339e6225c47SLeila Ghaffari                                    {dq[0][3][i],
340e6225c47SLeila Ghaffari                                     dq[1][3][i],
341e6225c47SLeila Ghaffari                                     dq[2][3][i]}
342e6225c47SLeila Ghaffari                                   };
343e6225c47SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
344e6225c47SLeila Ghaffari                                     dq[1][4][i],
345e6225c47SLeila Ghaffari                                     dq[2][4][i]
346e6225c47SLeila Ghaffari                                    };
34777841947SLeila Ghaffari     // -- Interp-to-Interp q_data
34877841947SLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
34977841947SLeila Ghaffari     // -- Interp-to-Grad q_data
35077841947SLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
35177841947SLeila Ghaffari     // *INDENT-OFF*
35277841947SLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
35377841947SLeila Ghaffari                                     q_data[2][i],
35477841947SLeila Ghaffari                                     q_data[3][i]},
35577841947SLeila Ghaffari                                    {q_data[4][i],
35677841947SLeila Ghaffari                                     q_data[5][i],
35777841947SLeila Ghaffari                                     q_data[6][i]},
35877841947SLeila Ghaffari                                    {q_data[7][i],
35977841947SLeila Ghaffari                                     q_data[8][i],
36077841947SLeila Ghaffari                                     q_data[9][i]}
36177841947SLeila Ghaffari                                   };
36277841947SLeila Ghaffari     // *INDENT-ON*
363e6225c47SLeila Ghaffari     // dU/dx
364e6225c47SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
365e6225c47SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
366e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
367e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
368e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
369e6225c47SLeila Ghaffari       for (int k=0; k<3; k++) {
370e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
371e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
372e6225c47SLeila Ghaffari         for (int l=0; l<3; l++) {
373e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
374e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
375e6225c47SLeila Ghaffari         }
376e6225c47SLeila Ghaffari       }
377e6225c47SLeila Ghaffari     }
378e6225c47SLeila Ghaffari     // Pressure
37977841947SLeila Ghaffari     const CeedScalar
38077841947SLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
38177841947SLeila Ghaffari     E_internal = E - E_kinetic,
382e6225c47SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
38377841947SLeila Ghaffari 
38477841947SLeila Ghaffari     // The Physics
38577841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
38677841947SLeila Ghaffari     for (int j=0; j<5; j++) {
387e6225c47SLeila Ghaffari       v[j][i] = 0.;
38877841947SLeila Ghaffari       for (int k=0; k<3; k++)
389e6225c47SLeila Ghaffari         dv[k][j][i] = 0.;
39077841947SLeila Ghaffari     }
39177841947SLeila Ghaffari 
39277841947SLeila Ghaffari     // -- Density
39377841947SLeila Ghaffari     // ---- u rho
39477841947SLeila Ghaffari     for (int j=0; j<3; j++)
39577841947SLeila Ghaffari       dv[j][0][i]  += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
39677841947SLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
39777841947SLeila Ghaffari     // -- Momentum
39877841947SLeila Ghaffari     // ---- rho (u x u) + P I3
39977841947SLeila Ghaffari     for (int j=0; j<3; j++)
40077841947SLeila Ghaffari       for (int k=0; k<3; k++)
401e6225c47SLeila Ghaffari         dv[k][j+1][i]  += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
402e6225c47SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
403e6225c47SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
40477841947SLeila Ghaffari     // -- Total Energy Density
40577841947SLeila Ghaffari     // ---- (E + P) u
40677841947SLeila Ghaffari     for (int j=0; j<3; j++)
40777841947SLeila Ghaffari       dv[j][4][i]  += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
40877841947SLeila Ghaffari                                          u[2]*dXdx[j][2]);
409e6225c47SLeila Ghaffari 
410e6225c47SLeila Ghaffari     // --Stabilization terms
411e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
412e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
413932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
414e6225c47SLeila Ghaffari 
415e6225c47SLeila Ghaffari     // ---- Transpose of the Jacobian
416e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv_T[3][5][5];
417e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
418e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
419e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
420e6225c47SLeila Ghaffari           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
421e6225c47SLeila Ghaffari 
422e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
423e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
424e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
425e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
426e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
427e6225c47SLeila Ghaffari       for (int k=0; k<3; k++)
428e6225c47SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
429e6225c47SLeila Ghaffari     }
430e6225c47SLeila Ghaffari 
431e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
432e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
433e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
434e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
435e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
436e6225c47SLeila Ghaffari           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
437e6225c47SLeila Ghaffari 
438932417b3SJed Brown     // Stabilization
439932417b3SJed Brown     // -- Tau elements
440932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
441932417b3SJed Brown     CeedScalar Tau_x[3] = {0.};
442932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
443e6225c47SLeila Ghaffari 
444932417b3SJed Brown     // -- Stabilization method: none or SU
445e6225c47SLeila Ghaffari     CeedScalar stab[5][3];
446e6225c47SLeila Ghaffari     switch (context->stabilization) {
447e6225c47SLeila Ghaffari     case 0:        // Galerkin
448e6225c47SLeila Ghaffari       break;
449e6225c47SLeila Ghaffari     case 1:        // SU
450e6225c47SLeila Ghaffari       for (int j=0; j<3; j++)
451e6225c47SLeila Ghaffari         for (int k=0; k<5; k++)
452e6225c47SLeila Ghaffari           for (int l=0; l<5; l++)
453932417b3SJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
454e6225c47SLeila Ghaffari 
455e6225c47SLeila Ghaffari       for (int j=0; j<5; j++)
456e6225c47SLeila Ghaffari         for (int k=0; k<3; k++)
457e6225c47SLeila Ghaffari           dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] +
458e6225c47SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
459e6225c47SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
460e6225c47SLeila Ghaffari       break;
461e6225c47SLeila Ghaffari     case 2:        // SUPG is not implemented for explicit scheme
462e6225c47SLeila Ghaffari       break;
463e6225c47SLeila Ghaffari     }
464e6225c47SLeila Ghaffari 
46577841947SLeila Ghaffari   } // End Quadrature Point Loop
46677841947SLeila Ghaffari 
46777841947SLeila Ghaffari   // Return
46877841947SLeila Ghaffari   return 0;
46977841947SLeila Ghaffari }
47077841947SLeila Ghaffari 
47177841947SLeila Ghaffari // *****************************************************************************
47277841947SLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above)
47377841947SLeila Ghaffari //   with implicit time stepping method
47477841947SLeila Ghaffari //
47577841947SLeila Ghaffari // *****************************************************************************
47677841947SLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q,
47777841947SLeila Ghaffari                                 const CeedScalar *const *in, CeedScalar *const *out) {
47877841947SLeila Ghaffari   // *INDENT-OFF*
47977841947SLeila Ghaffari   // Inputs
48077841947SLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
481e6225c47SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
48277841947SLeila Ghaffari                    (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
48377841947SLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
48477841947SLeila Ghaffari   // Outputs
48577841947SLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
48677841947SLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
48777841947SLeila Ghaffari 
488e6225c47SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
489932417b3SJed Brown   const CeedScalar c_tau = context->c_tau;
49077841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
49177841947SLeila Ghaffari 
49277841947SLeila Ghaffari   CeedPragmaSIMD
49377841947SLeila Ghaffari   // Quadrature Point Loop
49477841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
49577841947SLeila Ghaffari     // *INDENT-OFF*
49677841947SLeila Ghaffari     // Setup
49777841947SLeila Ghaffari     // -- Interp in
49877841947SLeila Ghaffari     const CeedScalar rho        =   q[0][i];
49977841947SLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
50077841947SLeila Ghaffari                                     q[2][i] / rho,
50177841947SLeila Ghaffari                                     q[3][i] / rho
50277841947SLeila Ghaffari                                    };
50377841947SLeila Ghaffari     const CeedScalar E          =   q[4][i];
504e6225c47SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
505e6225c47SLeila Ghaffari                                     dq[1][0][i],
506e6225c47SLeila Ghaffari                                     dq[2][0][i]
507e6225c47SLeila Ghaffari                                    };
508e6225c47SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
509e6225c47SLeila Ghaffari                                     dq[1][1][i],
510e6225c47SLeila Ghaffari                                     dq[2][1][i]},
511e6225c47SLeila Ghaffari                                    {dq[0][2][i],
512e6225c47SLeila Ghaffari                                     dq[1][2][i],
513e6225c47SLeila Ghaffari                                     dq[2][2][i]},
514e6225c47SLeila Ghaffari                                    {dq[0][3][i],
515e6225c47SLeila Ghaffari                                     dq[1][3][i],
516e6225c47SLeila Ghaffari                                     dq[2][3][i]}
517e6225c47SLeila Ghaffari                                   };
518e6225c47SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
519e6225c47SLeila Ghaffari                                     dq[1][4][i],
520e6225c47SLeila Ghaffari                                     dq[2][4][i]
521e6225c47SLeila Ghaffari                                    };
52277841947SLeila Ghaffari     // -- Interp-to-Interp q_data
52377841947SLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
52477841947SLeila Ghaffari     // -- Interp-to-Grad q_data
52577841947SLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
52677841947SLeila Ghaffari     // *INDENT-OFF*
52777841947SLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
52877841947SLeila Ghaffari                                     q_data[2][i],
52977841947SLeila Ghaffari                                     q_data[3][i]},
53077841947SLeila Ghaffari                                    {q_data[4][i],
53177841947SLeila Ghaffari                                     q_data[5][i],
53277841947SLeila Ghaffari                                     q_data[6][i]},
53377841947SLeila Ghaffari                                    {q_data[7][i],
53477841947SLeila Ghaffari                                     q_data[8][i],
53577841947SLeila Ghaffari                                     q_data[9][i]}
53677841947SLeila Ghaffari                                   };
53777841947SLeila Ghaffari     // *INDENT-ON*
538e6225c47SLeila Ghaffari     // dU/dx
539e6225c47SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
540e6225c47SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
541e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
542e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
543e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
544e6225c47SLeila Ghaffari       for (int k=0; k<3; k++) {
545e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
546e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
547e6225c47SLeila Ghaffari         for (int l=0; l<3; l++) {
548e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
549e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
550e6225c47SLeila Ghaffari         }
551e6225c47SLeila Ghaffari       }
552e6225c47SLeila Ghaffari     }
55377841947SLeila Ghaffari     const CeedScalar
55477841947SLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
55577841947SLeila Ghaffari     E_internal = E - E_kinetic,
556e6225c47SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
55777841947SLeila Ghaffari 
55877841947SLeila Ghaffari     // The Physics
55977841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
56077841947SLeila Ghaffari     for (int j=0; j<5; j++) {
561e6225c47SLeila Ghaffari       v[j][i] = 0.;
56277841947SLeila Ghaffari       for (int k=0; k<3; k++)
563e6225c47SLeila Ghaffari         dv[k][j][i] = 0.;
56477841947SLeila Ghaffari     }
56577841947SLeila Ghaffari     //-----mass matrix
56677841947SLeila Ghaffari     for (int j=0; j<5; j++)
56777841947SLeila Ghaffari       v[j][i] += wdetJ*q_dot[j][i];
56877841947SLeila Ghaffari 
56977841947SLeila Ghaffari     // -- Density
57077841947SLeila Ghaffari     // ---- u rho
57177841947SLeila Ghaffari     for (int j=0; j<3; j++)
57277841947SLeila Ghaffari       dv[j][0][i]  -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
57377841947SLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
57477841947SLeila Ghaffari     // -- Momentum
57577841947SLeila Ghaffari     // ---- rho (u x u) + P I3
57677841947SLeila Ghaffari     for (int j=0; j<3; j++)
57777841947SLeila Ghaffari       for (int k=0; k<3; k++)
578e6225c47SLeila Ghaffari         dv[k][j+1][i]  -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
579e6225c47SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
580e6225c47SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
58177841947SLeila Ghaffari     // -- Total Energy Density
58277841947SLeila Ghaffari     // ---- (E + P) u
58377841947SLeila Ghaffari     for (int j=0; j<3; j++)
58477841947SLeila Ghaffari       dv[j][4][i]  -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
58577841947SLeila Ghaffari                                          u[2]*dXdx[j][2]);
586e6225c47SLeila Ghaffari 
587e6225c47SLeila Ghaffari     // -- Stabilization terms
588e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
589e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
590932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
591e6225c47SLeila Ghaffari 
592e6225c47SLeila Ghaffari     // ---- Transpose of the Jacobian
593e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv_T[3][5][5];
594e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
595e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
596e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
597e6225c47SLeila Ghaffari           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
598e6225c47SLeila Ghaffari 
599e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
600e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
601e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
602e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
603e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
604e6225c47SLeila Ghaffari       for (int k=0; k<3; k++)
605e6225c47SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
606e6225c47SLeila Ghaffari     }
607e6225c47SLeila Ghaffari 
608e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
609e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
610e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
611e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
612e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
613e6225c47SLeila Ghaffari           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
614e6225c47SLeila Ghaffari 
615e6225c47SLeila Ghaffari     // ---- Strong residual
616e6225c47SLeila Ghaffari     CeedScalar strong_res[5];
617e6225c47SLeila Ghaffari     for (int j=0; j<5; j++)
618e6225c47SLeila Ghaffari       strong_res[j] = q_dot[j][i] + strong_conv[j];
619e6225c47SLeila Ghaffari 
620932417b3SJed Brown     // Stabilization
621932417b3SJed Brown     // -- Tau elements
622932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
623932417b3SJed Brown     CeedScalar Tau_x[3] = {0.};
624932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
625e6225c47SLeila Ghaffari 
626932417b3SJed Brown     // -- Stabilization method: none, SU, or SUPG
627e6225c47SLeila Ghaffari     CeedScalar stab[5][3];
628e6225c47SLeila Ghaffari     switch (context->stabilization) {
629e6225c47SLeila Ghaffari     case 0:        // Galerkin
630e6225c47SLeila Ghaffari       break;
631e6225c47SLeila Ghaffari     case 1:        // SU
632e6225c47SLeila Ghaffari       for (int j=0; j<3; j++)
633e6225c47SLeila Ghaffari         for (int k=0; k<5; k++)
634e6225c47SLeila Ghaffari           for (int l=0; l<5; l++)
635932417b3SJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
636e6225c47SLeila Ghaffari 
637e6225c47SLeila Ghaffari       for (int j=0; j<5; j++)
638e6225c47SLeila Ghaffari         for (int k=0; k<3; k++)
639e6225c47SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
640e6225c47SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
641e6225c47SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
642e6225c47SLeila Ghaffari       break;
643e6225c47SLeila Ghaffari     case 2:        // SUPG
644e6225c47SLeila Ghaffari       for (int j=0; j<3; j++)
645e6225c47SLeila Ghaffari         for (int k=0; k<5; k++)
646e6225c47SLeila Ghaffari           for (int l=0; l<5; l++)
647932417b3SJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_res[l];
648e6225c47SLeila Ghaffari 
649e6225c47SLeila Ghaffari       for (int j=0; j<5; j++)
650e6225c47SLeila Ghaffari         for (int k=0; k<3; k++)
651e6225c47SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
652e6225c47SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
653e6225c47SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
654e6225c47SLeila Ghaffari       break;
655e6225c47SLeila Ghaffari     }
65677841947SLeila Ghaffari   } // End Quadrature Point Loop
65777841947SLeila Ghaffari 
65877841947SLeila Ghaffari   // Return
65977841947SLeila Ghaffari   return 0;
66077841947SLeila Ghaffari }
66177841947SLeila Ghaffari // *****************************************************************************
6622fe7aee7SLeila Ghaffari // This QFunction sets the inflow boundary conditions for
6632fe7aee7SLeila Ghaffari //   the traveling vortex problem.
66477841947SLeila Ghaffari //
66577841947SLeila Ghaffari //  Prescribed T_inlet and P_inlet are converted to conservative variables
66677841947SLeila Ghaffari //      and applied weakly.
66777841947SLeila Ghaffari //
66877841947SLeila Ghaffari // *****************************************************************************
6692fe7aee7SLeila Ghaffari CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q,
67077841947SLeila Ghaffari                                        const CeedScalar *const *in,
67177841947SLeila Ghaffari                                        CeedScalar *const *out) {
67277841947SLeila Ghaffari   // *INDENT-OFF*
67377841947SLeila Ghaffari   // Inputs
6742fe7aee7SLeila Ghaffari   const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
67577841947SLeila Ghaffari   // Outputs
67677841947SLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
67777841947SLeila Ghaffari   // *INDENT-ON*
67877841947SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
67977841947SLeila Ghaffari   const int euler_test      = context->euler_test;
68077841947SLeila Ghaffari   const bool implicit       = context->implicit;
68177841947SLeila Ghaffari   CeedScalar *mean_velocity = context->mean_velocity;
68277841947SLeila Ghaffari   const CeedScalar cv    = 2.5;
68377841947SLeila Ghaffari   const CeedScalar R     = 1.;
68477841947SLeila Ghaffari   CeedScalar T_inlet;
68577841947SLeila Ghaffari   CeedScalar P_inlet;
68677841947SLeila Ghaffari 
68777841947SLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
68877841947SLeila Ghaffari   if (euler_test == 1 || euler_test == 3)
68977841947SLeila Ghaffari     for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.;
69077841947SLeila Ghaffari 
69177841947SLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
69277841947SLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
69377841947SLeila Ghaffari   else T_inlet = P_inlet = 1.;
69477841947SLeila Ghaffari 
69577841947SLeila Ghaffari   CeedPragmaSIMD
69677841947SLeila Ghaffari   // Quadrature Point Loop
69777841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
69877841947SLeila Ghaffari     // Setup
69977841947SLeila Ghaffari     // -- Interp-to-Interp q_data
70077841947SLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
70177841947SLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
70277841947SLeila Ghaffari     // We can effect this by swapping the sign on this weight
70377841947SLeila Ghaffari     const CeedScalar wdetJb  = (implicit ? -1. : 1.) * q_data_sur[0][i];
7042fe7aee7SLeila Ghaffari     // ---- Normal vect
70577841947SLeila Ghaffari     const CeedScalar norm[3] = {q_data_sur[1][i],
70677841947SLeila Ghaffari                                 q_data_sur[2][i],
70777841947SLeila Ghaffari                                 q_data_sur[3][i]
70877841947SLeila Ghaffari                                };
70977841947SLeila Ghaffari 
71077841947SLeila Ghaffari     // face_normal = Normal vector of the face
71177841947SLeila Ghaffari     const CeedScalar face_normal = norm[0]*mean_velocity[0] +
71277841947SLeila Ghaffari                                    norm[1]*mean_velocity[1] +
71377841947SLeila Ghaffari                                    norm[2]*mean_velocity[2];
71477841947SLeila Ghaffari     // The Physics
71577841947SLeila Ghaffari     // Zero v so all future terms can safely sum into it
716e6225c47SLeila Ghaffari     for (int j=0; j<5; j++) v[j][i] = 0.;
71777841947SLeila Ghaffari 
71877841947SLeila Ghaffari     // Implementing in/outflow BCs
7192fe7aee7SLeila Ghaffari     if (face_normal > 0) {
72077841947SLeila Ghaffari     } else { // inflow
72177841947SLeila Ghaffari       const CeedScalar rho_inlet = P_inlet/(R*T_inlet);
72277841947SLeila Ghaffari       const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] +
72377841947SLeila Ghaffari                                           mean_velocity[1]*mean_velocity[1]) / 2.;
72477841947SLeila Ghaffari       // incoming total energy
72577841947SLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
72677841947SLeila Ghaffari 
72777841947SLeila Ghaffari       // The Physics
72877841947SLeila Ghaffari       // -- Density
72977841947SLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
73077841947SLeila Ghaffari 
73177841947SLeila Ghaffari       // -- Momentum
73277841947SLeila Ghaffari       for (int j=0; j<3; j++)
73377841947SLeila Ghaffari         v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] +
73477841947SLeila Ghaffari                               norm[j] * P_inlet);
73577841947SLeila Ghaffari 
73677841947SLeila Ghaffari       // -- Total Energy Density
73777841947SLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
73877841947SLeila Ghaffari     }
73977841947SLeila Ghaffari 
74077841947SLeila Ghaffari   } // End Quadrature Point Loop
74177841947SLeila Ghaffari   return 0;
74277841947SLeila Ghaffari }
74377841947SLeila Ghaffari 
74477841947SLeila Ghaffari // *****************************************************************************
74555e76554SLeila Ghaffari // This QFunction sets the outflow boundary conditions for
74655e76554SLeila Ghaffari //   the Euler solver.
74755e76554SLeila Ghaffari //
74855e76554SLeila Ghaffari //  Outflow BCs:
74955e76554SLeila Ghaffari //    The validity of the weak form of the governing equations is
75055e76554SLeila Ghaffari //      extended to the outflow.
75155e76554SLeila Ghaffari //
75255e76554SLeila Ghaffari // *****************************************************************************
75355e76554SLeila Ghaffari CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q,
75455e76554SLeila Ghaffari                               const CeedScalar *const *in,
75555e76554SLeila Ghaffari                               CeedScalar *const *out) {
75655e76554SLeila Ghaffari   // *INDENT-OFF*
75755e76554SLeila Ghaffari   // Inputs
75855e76554SLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
75955e76554SLeila Ghaffari                    (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
76055e76554SLeila Ghaffari   // Outputs
76155e76554SLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
76255e76554SLeila Ghaffari   // *INDENT-ON*
76355e76554SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
76455e76554SLeila Ghaffari   const bool implicit       = context->implicit;
76555e76554SLeila Ghaffari   CeedScalar *mean_velocity = context->mean_velocity;
76655e76554SLeila Ghaffari 
76755e76554SLeila Ghaffari   const CeedScalar gamma = 1.4;
76855e76554SLeila Ghaffari 
76955e76554SLeila Ghaffari   CeedPragmaSIMD
77055e76554SLeila Ghaffari   // Quadrature Point Loop
77155e76554SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
77255e76554SLeila Ghaffari     // Setup
77355e76554SLeila Ghaffari     // -- Interp in
77455e76554SLeila Ghaffari     const CeedScalar rho      =  q[0][i];
77555e76554SLeila Ghaffari     const CeedScalar u[3]     = {q[1][i] / rho,
77655e76554SLeila Ghaffari                                  q[2][i] / rho,
77755e76554SLeila Ghaffari                                  q[3][i] / rho
77855e76554SLeila Ghaffari                                 };
77955e76554SLeila Ghaffari     const CeedScalar E        =  q[4][i];
78055e76554SLeila Ghaffari 
78155e76554SLeila Ghaffari     // -- Interp-to-Interp q_data
78255e76554SLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
78355e76554SLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
78455e76554SLeila Ghaffari     // We can effect this by swapping the sign on this weight
78555e76554SLeila Ghaffari     const CeedScalar wdetJb     =   (implicit ? -1. : 1.) * q_data_sur[0][i];
78655e76554SLeila Ghaffari     // ---- Normal vectors
78755e76554SLeila Ghaffari     const CeedScalar norm[3]    =   {q_data_sur[1][i],
78855e76554SLeila Ghaffari                                      q_data_sur[2][i],
78955e76554SLeila Ghaffari                                      q_data_sur[3][i]
79055e76554SLeila Ghaffari                                     };
79155e76554SLeila Ghaffari 
79255e76554SLeila Ghaffari     // face_normal = Normal vector of the face
79355e76554SLeila Ghaffari     const CeedScalar face_normal = norm[0]*mean_velocity[0] +
79455e76554SLeila Ghaffari                                    norm[1]*mean_velocity[1] +
79555e76554SLeila Ghaffari                                    norm[2]*mean_velocity[2];
79655e76554SLeila Ghaffari     // The Physics
79755e76554SLeila Ghaffari     // Zero v so all future terms can safely sum into it
79855e76554SLeila Ghaffari     for (int j=0; j<5; j++) v[j][i] = 0;
79955e76554SLeila Ghaffari 
80055e76554SLeila Ghaffari     // Implementing in/outflow BCs
80155e76554SLeila Ghaffari     if (face_normal > 0) { // outflow
80255e76554SLeila Ghaffari       const CeedScalar E_kinetic = (u[0]*u[0] + u[1]*u[1]) / 2.;
80355e76554SLeila Ghaffari       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.); // pressure
80455e76554SLeila Ghaffari       const CeedScalar u_normal  = norm[0]*u[0] + norm[1]*u[1] +
80555e76554SLeila Ghaffari                                    norm[2]*u[2]; // Normal velocity
80655e76554SLeila Ghaffari       // The Physics
80755e76554SLeila Ghaffari       // -- Density
80855e76554SLeila Ghaffari       v[0][i] -= wdetJb * rho * u_normal;
80955e76554SLeila Ghaffari 
81055e76554SLeila Ghaffari       // -- Momentum
81155e76554SLeila Ghaffari       for (int j=0; j<3; j++)
81255e76554SLeila Ghaffari         v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + norm[j] * P);
81355e76554SLeila Ghaffari 
81455e76554SLeila Ghaffari       // -- Total Energy Density
81555e76554SLeila Ghaffari       v[4][i] -= wdetJb * u_normal * (E + P);
81655e76554SLeila Ghaffari     }
81755e76554SLeila Ghaffari   } // End Quadrature Point Loop
81855e76554SLeila Ghaffari   return 0;
81955e76554SLeila Ghaffari }
82055e76554SLeila Ghaffari 
82155e76554SLeila Ghaffari // *****************************************************************************
82277841947SLeila Ghaffari 
82377841947SLeila Ghaffari #endif // eulervortex_h
824