xref: /libCEED/examples/fluids/qfunctions/eulervortex.h (revision 2fe7aee7a837446ac12fe008ea12313e81769bbb)
177841947SLeila Ghaffari // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
277841947SLeila Ghaffari // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
377841947SLeila Ghaffari // reserved. See files LICENSE and NOTICE for details.
477841947SLeila Ghaffari //
577841947SLeila Ghaffari // This file is part of CEED, a collection of benchmarks, miniapps, software
677841947SLeila Ghaffari // libraries and APIs for efficient high-order finite element and spectral
777841947SLeila Ghaffari // element discretizations for exascale applications. For more information and
877841947SLeila Ghaffari // source code availability see http://github.com/ceed.
977841947SLeila Ghaffari //
1077841947SLeila Ghaffari // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
1177841947SLeila Ghaffari // a collaborative effort of two U.S. Department of Energy organizations (Office
1277841947SLeila Ghaffari // of Science and the National Nuclear Security Administration) responsible for
1377841947SLeila Ghaffari // the planning and preparation of a capable exascale ecosystem, including
1477841947SLeila Ghaffari // software, applications, hardware, advanced system engineering and early
1577841947SLeila Ghaffari // testbed platforms, in support of the nation's exascale computing imperative.
1677841947SLeila Ghaffari 
1777841947SLeila Ghaffari /// @file
1877841947SLeila Ghaffari /// Euler traveling vortex initial condition and operator for Navier-Stokes
1977841947SLeila Ghaffari /// example using PETSc
2077841947SLeila Ghaffari 
2177841947SLeila Ghaffari // Model from:
2277841947SLeila Ghaffari //   On the Order of Accuracy and Numerical Performance of Two Classes of
2377841947SLeila Ghaffari //   Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011).
2477841947SLeila Ghaffari 
2577841947SLeila Ghaffari #ifndef eulervortex_h
2677841947SLeila Ghaffari #define eulervortex_h
2777841947SLeila Ghaffari 
2877841947SLeila Ghaffari #include <math.h>
2977841947SLeila Ghaffari 
3077841947SLeila Ghaffari #ifndef M_PI
3177841947SLeila Ghaffari #define M_PI    3.14159265358979323846
3277841947SLeila Ghaffari #endif
3377841947SLeila Ghaffari 
3477841947SLeila Ghaffari #ifndef euler_context_struct
3577841947SLeila Ghaffari #define euler_context_struct
3677841947SLeila Ghaffari typedef struct EulerContext_ *EulerContext;
3777841947SLeila Ghaffari struct EulerContext_ {
3877841947SLeila Ghaffari   CeedScalar center[3];
3977841947SLeila Ghaffari   CeedScalar curr_time;
4077841947SLeila Ghaffari   CeedScalar vortex_strength;
41932417b3SJed Brown   CeedScalar c_tau;
4277841947SLeila Ghaffari   CeedScalar mean_velocity[3];
4377841947SLeila Ghaffari   bool implicit;
44e6225c47SLeila Ghaffari   int euler_test;
45e6225c47SLeila Ghaffari   int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG
4677841947SLeila Ghaffari };
4777841947SLeila Ghaffari #endif
4877841947SLeila Ghaffari 
4977841947SLeila Ghaffari // *****************************************************************************
5077841947SLeila Ghaffari // This function sets the initial conditions
5177841947SLeila Ghaffari //
5277841947SLeila Ghaffari //   Temperature:
5377841947SLeila Ghaffari //     T   = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
5477841947SLeila Ghaffari //   Density:
5577841947SLeila Ghaffari //     rho = (T/S_vortex)^(1 / (gamma - 1))
5677841947SLeila Ghaffari //   Pressure:
5777841947SLeila Ghaffari //     P   = rho * T
5877841947SLeila Ghaffari //   Velocity:
5977841947SLeila Ghaffari //     ui  = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
6077841947SLeila Ghaffari //     r   = sqrt( (x - xc)**2 + (y - yc)**2 )
6177841947SLeila Ghaffari //   Velocity/Momentum Density:
6277841947SLeila Ghaffari //     Ui  = rho ui
6377841947SLeila Ghaffari //   Total Energy:
6477841947SLeila Ghaffari //     E   = P / (gamma - 1) + rho (u u)/2
6577841947SLeila Ghaffari //
6677841947SLeila Ghaffari // Constants:
6777841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
6877841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
6977841947SLeila Ghaffari //   vortex_strength ,  Strength of vortex
7077841947SLeila Ghaffari //   center          ,  Location of bubble center
7177841947SLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
7277841947SLeila Ghaffari //
7377841947SLeila Ghaffari // *****************************************************************************
7477841947SLeila Ghaffari 
7577841947SLeila Ghaffari // *****************************************************************************
7677841947SLeila Ghaffari // This helper function provides support for the exact, time-dependent solution
7777841947SLeila Ghaffari //   (currently not implemented) and IC formulation for Euler traveling vortex
7877841947SLeila Ghaffari // *****************************************************************************
7977841947SLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time,
8077841947SLeila Ghaffari                                       const CeedScalar X[], CeedInt Nf, CeedScalar q[],
8177841947SLeila Ghaffari                                       void *ctx) {
8277841947SLeila Ghaffari   // Context
8377841947SLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
8477841947SLeila Ghaffari   const CeedScalar vortex_strength    = context->vortex_strength;
8577841947SLeila Ghaffari   const CeedScalar *center            = context->center; // Center of the domain
8677841947SLeila Ghaffari   const CeedScalar *mean_velocity = context->mean_velocity;
8777841947SLeila Ghaffari 
8877841947SLeila Ghaffari   // Setup
8977841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
9077841947SLeila Ghaffari   const CeedScalar cv    = 2.5;
9177841947SLeila Ghaffari   const CeedScalar R     = 1.;
9277841947SLeila Ghaffari   const CeedScalar x     = X[0], y = X[1]; // Coordinates
9377841947SLeila Ghaffari   // Vortex center
9477841947SLeila Ghaffari   const CeedScalar xc = center[0] + mean_velocity[0] * time;
9577841947SLeila Ghaffari   const CeedScalar yc = center[1] + mean_velocity[1] * time;
9677841947SLeila Ghaffari 
9777841947SLeila Ghaffari   const CeedScalar x0       = x - xc;
9877841947SLeila Ghaffari   const CeedScalar y0       = y - yc;
9977841947SLeila Ghaffari   const CeedScalar r        = sqrt( x0*x0 + y0*y0 );
10077841947SLeila Ghaffari   const CeedScalar C        = vortex_strength * exp((1. - r*r)/2.) / (2. * M_PI);
101e6225c47SLeila Ghaffari   const CeedScalar delta_T  = - (gamma - 1.) * vortex_strength * vortex_strength *
102e6225c47SLeila Ghaffari                               exp(1 - r*r) / (8. * gamma * M_PI * M_PI);
10377841947SLeila Ghaffari   const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma
10477841947SLeila Ghaffari   const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength /
10577841947SLeila Ghaffari                               (8.*gamma*M_PI*M_PI);
10677841947SLeila Ghaffari   CeedScalar rho, P, T, E, u[3] = {0.};
10777841947SLeila Ghaffari 
10877841947SLeila Ghaffari   // Initial Conditions
10977841947SLeila Ghaffari   switch (context->euler_test) {
11077841947SLeila Ghaffari   case 0: // Traveling vortex
11177841947SLeila Ghaffari     T = 1 + delta_T;
11277841947SLeila Ghaffari     // P = rho * T
11377841947SLeila Ghaffari     // P = S * rho^gamma
11477841947SLeila Ghaffari     // Solve for rho, then substitute for P
115e6225c47SLeila Ghaffari     rho  = pow(T/S_vortex, 1 / (gamma - 1.));
11677841947SLeila Ghaffari     P    = rho * T;
11777841947SLeila Ghaffari     u[0] = mean_velocity[0] - C*y0;
11877841947SLeila Ghaffari     u[1] = mean_velocity[1] + C*x0;
11977841947SLeila Ghaffari 
12077841947SLeila Ghaffari     // Assign exact solution
12177841947SLeila Ghaffari     q[0] = rho;
12277841947SLeila Ghaffari     q[1] = rho * u[0];
12377841947SLeila Ghaffari     q[2] = rho * u[1];
12477841947SLeila Ghaffari     q[3] = rho * u[2];
12577841947SLeila Ghaffari     q[4] = P / (gamma - 1.) + rho * (u[0]*u[0] + u[1]*u[1]) / 2.;
12677841947SLeila Ghaffari     break;
12777841947SLeila Ghaffari   case 1: // Constant zero velocity, density constant, total energy constant
12877841947SLeila Ghaffari     rho  = 1.;
12977841947SLeila Ghaffari     E    = 2.;
13077841947SLeila Ghaffari 
13177841947SLeila Ghaffari     // Assign exact solution
13277841947SLeila Ghaffari     q[0] = rho;
13377841947SLeila Ghaffari     q[1] = rho * u[0];
13477841947SLeila Ghaffari     q[2] = rho * u[1];
13577841947SLeila Ghaffari     q[3] = rho * u[2];
13677841947SLeila Ghaffari     q[4] = E;
13777841947SLeila Ghaffari     break;
13877841947SLeila Ghaffari   case 2: // Constant nonzero velocity, density constant, total energy constant
13977841947SLeila Ghaffari     rho  = 1.;
14077841947SLeila Ghaffari     E    = 2.;
14177841947SLeila Ghaffari     u[0] = mean_velocity[0];
14277841947SLeila Ghaffari     u[1] = mean_velocity[1];
14377841947SLeila Ghaffari 
14477841947SLeila Ghaffari     // Assign exact solution
14577841947SLeila Ghaffari     q[0] = rho;
14677841947SLeila Ghaffari     q[1] = rho * u[0];
14777841947SLeila Ghaffari     q[2] = rho * u[1];
14877841947SLeila Ghaffari     q[3] = rho * u[2];
14977841947SLeila Ghaffari     q[4] = E;
15077841947SLeila Ghaffari     break;
15177841947SLeila Ghaffari   case 3: // Velocity zero, pressure constant
15277841947SLeila Ghaffari     // (so density and internal energy will be non-constant),
15377841947SLeila Ghaffari     // but the velocity should stay zero and the bubble won't diffuse
15477841947SLeila Ghaffari     // (for Euler, where there is no thermal conductivity)
15577841947SLeila Ghaffari     P    = 1.;
15677841947SLeila Ghaffari     T    = 1. - S_bubble * exp(1. - r*r);
15777841947SLeila Ghaffari     rho  = P / (R*T);
15877841947SLeila Ghaffari 
15977841947SLeila Ghaffari     // Assign exact solution
16077841947SLeila Ghaffari     q[0] = rho;
16177841947SLeila Ghaffari     q[1] = rho * u[0];
16277841947SLeila Ghaffari     q[2] = rho * u[1];
16377841947SLeila Ghaffari     q[3] = rho * u[2];
16477841947SLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
16577841947SLeila Ghaffari     break;
16677841947SLeila Ghaffari   case 4: // Constant nonzero velocity, pressure constant
16777841947SLeila Ghaffari     // (so density and internal energy will be non-constant),
16877841947SLeila Ghaffari     // it should be transported across the domain, but velocity stays constant
16977841947SLeila Ghaffari     P    = 1.;
17077841947SLeila Ghaffari     T    = 1. - S_bubble * exp(1. - r*r);
17177841947SLeila Ghaffari     rho  = P / (R*T);
17277841947SLeila Ghaffari     u[0] = mean_velocity[0];
17377841947SLeila Ghaffari     u[1] = mean_velocity[1];
17477841947SLeila Ghaffari 
17577841947SLeila Ghaffari     // Assign exact solution
17677841947SLeila Ghaffari     q[0] = rho;
17777841947SLeila Ghaffari     q[1] = rho * u[0];
17877841947SLeila Ghaffari     q[2] = rho * u[1];
17977841947SLeila Ghaffari     q[3] = rho * u[2];
18077841947SLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
18177841947SLeila Ghaffari     break;
18232f166c6SLeila Ghaffari   case 5: // non-smooth thermal bubble - cylinder
18332f166c6SLeila Ghaffari     P    = 1.;
18432f166c6SLeila Ghaffari     T = 1. - (r < 1. ? S_bubble : 0.);
18532f166c6SLeila Ghaffari     rho  = P / (R*T);
18632f166c6SLeila Ghaffari     u[0] = mean_velocity[0];
18732f166c6SLeila Ghaffari     u[1] = mean_velocity[1];
18832f166c6SLeila Ghaffari 
18932f166c6SLeila Ghaffari     // Assign exact solution
19032f166c6SLeila Ghaffari     q[0] = rho;
19132f166c6SLeila Ghaffari     q[1] = rho * u[0];
19232f166c6SLeila Ghaffari     q[2] = rho * u[1];
19332f166c6SLeila Ghaffari     q[3] = rho * u[2];
19432f166c6SLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
19532f166c6SLeila Ghaffari     break;
19677841947SLeila Ghaffari   }
19777841947SLeila Ghaffari   // Return
19877841947SLeila Ghaffari   return 0;
19977841947SLeila Ghaffari }
20077841947SLeila Ghaffari 
20177841947SLeila Ghaffari // *****************************************************************************
202e6225c47SLeila Ghaffari // Helper function for computing flux Jacobian
203e6225c47SLeila Ghaffari // *****************************************************************************
204932417b3SJed Brown CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5],
205e6225c47SLeila Ghaffari     const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
206e6225c47SLeila Ghaffari     const CeedScalar gamma) {
207e6225c47SLeila Ghaffari   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square
208e6225c47SLeila Ghaffari   for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions
209e6225c47SLeila Ghaffari     for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix
210e6225c47SLeila Ghaffari       dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j];
211e6225c47SLeila Ghaffari       for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix
212e6225c47SLeila Ghaffari         dF[i][0][k+1]   = ((i==k) ? 1. : 0.);
213e6225c47SLeila Ghaffari         dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) +
214e6225c47SLeila Ghaffari                           ((i==k) ? u[j] : 0.) -
215e6225c47SLeila Ghaffari                           ((i==j) ? u[k] : 0.) * (gamma-1.);
216e6225c47SLeila Ghaffari         dF[i][4][k+1]   = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) -
217e6225c47SLeila Ghaffari                           (gamma-1.)*u[i]*u[k];
218e6225c47SLeila Ghaffari       }
219e6225c47SLeila Ghaffari       dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.);
220e6225c47SLeila Ghaffari     }
221e6225c47SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho);
222e6225c47SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
223e6225c47SLeila Ghaffari   }
224e6225c47SLeila Ghaffari }
225e6225c47SLeila Ghaffari 
226e6225c47SLeila Ghaffari // *****************************************************************************
227932417b3SJed Brown // Helper function for computing Tau elements (stabilization constant)
228932417b3SJed Brown //   Model from:
229932417b3SJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
230932417b3SJed Brown //
231932417b3SJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
232932417b3SJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
233932417b3SJed Brown //
234932417b3SJed Brown // Where
235932417b3SJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
236932417b3SJed Brown //   h[i]      = 2 length(dxdX[i])
237932417b3SJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
238932417b3SJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
239932417b3SJed Brown //   rho(A[i]) = spectral radius of the convective flux Jacobian i,
240932417b3SJed Brown //               wave speed in direction i
241932417b3SJed Brown // *****************************************************************************
242932417b3SJed Brown CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3],
243932417b3SJed Brown                                        const CeedScalar dXdx[3][3], const CeedScalar u[3],
244932417b3SJed Brown                                        const CeedScalar sound_speed, const CeedScalar c_tau) {
245932417b3SJed Brown   for (int i=0; i<3; i++) {
246932417b3SJed Brown     // length of element in direction i
247932417b3SJed Brown     CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] +
248932417b3SJed Brown                             dXdx[2][i]*dXdx[2][i]);
249932417b3SJed Brown     // fastest wave in direction i
250932417b3SJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
251932417b3SJed Brown     Tau_x[i] = c_tau * h / fastest_wave;
252932417b3SJed Brown   }
253932417b3SJed Brown }
254932417b3SJed Brown 
255932417b3SJed Brown // *****************************************************************************
25677841947SLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
25777841947SLeila Ghaffari // *****************************************************************************
25877841947SLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q,
25977841947SLeila Ghaffari                          const CeedScalar *const *in, CeedScalar *const *out) {
26077841947SLeila Ghaffari   // Inputs
26177841947SLeila Ghaffari   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
26277841947SLeila Ghaffari 
26377841947SLeila Ghaffari   // Outputs
26477841947SLeila Ghaffari   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
26577841947SLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
26677841947SLeila Ghaffari 
26777841947SLeila Ghaffari   CeedPragmaSIMD
26877841947SLeila Ghaffari   // Quadrature Point Loop
26977841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
27077841947SLeila Ghaffari     const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
271e6225c47SLeila Ghaffari     CeedScalar q[5] = {0.};
27277841947SLeila Ghaffari 
27377841947SLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
27477841947SLeila Ghaffari 
27577841947SLeila Ghaffari     for (CeedInt j=0; j<5; j++)
27677841947SLeila Ghaffari       q0[j][i] = q[j];
27777841947SLeila Ghaffari   } // End of Quadrature Point Loop
27877841947SLeila Ghaffari 
27977841947SLeila Ghaffari   // Return
28077841947SLeila Ghaffari   return 0;
28177841947SLeila Ghaffari }
28277841947SLeila Ghaffari 
28377841947SLeila Ghaffari // *****************************************************************************
28477841947SLeila Ghaffari // This QFunction implements the following formulation of Euler equations
28577841947SLeila Ghaffari //   with explicit time stepping method
28677841947SLeila Ghaffari //
28777841947SLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation
28877841947SLeila Ghaffari //   form with state variables of density, momentum density, and total
28977841947SLeila Ghaffari //   energy density.
29077841947SLeila Ghaffari //
29177841947SLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
29277841947SLeila Ghaffari //   rho - Mass Density
29377841947SLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
29477841947SLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
29577841947SLeila Ghaffari //
29677841947SLeila Ghaffari // Euler Equations:
29777841947SLeila Ghaffari //   drho/dt + div( U )                   = 0
29877841947SLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
29977841947SLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
30077841947SLeila Ghaffari //
30177841947SLeila Ghaffari // Equation of State:
30277841947SLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
30377841947SLeila Ghaffari //
30477841947SLeila Ghaffari // Constants:
30577841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
30677841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
30777841947SLeila Ghaffari //   g               ,  Gravity
30877841947SLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
30977841947SLeila Ghaffari // *****************************************************************************
31077841947SLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q,
31177841947SLeila Ghaffari                       const CeedScalar *const *in, CeedScalar *const *out) {
31277841947SLeila Ghaffari   // *INDENT-OFF*
31377841947SLeila Ghaffari   // Inputs
31477841947SLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
315e6225c47SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
31677841947SLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
31777841947SLeila Ghaffari   // Outputs
31877841947SLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
31977841947SLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
32077841947SLeila Ghaffari 
321e6225c47SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
322932417b3SJed Brown   const CeedScalar c_tau = context->c_tau;
32377841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
32477841947SLeila Ghaffari 
32577841947SLeila Ghaffari   CeedPragmaSIMD
32677841947SLeila Ghaffari   // Quadrature Point Loop
32777841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
32877841947SLeila Ghaffari     // *INDENT-OFF*
32977841947SLeila Ghaffari     // Setup
33077841947SLeila Ghaffari     // -- Interp in
33177841947SLeila Ghaffari     const CeedScalar rho        =   q[0][i];
33277841947SLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
33377841947SLeila Ghaffari                                     q[2][i] / rho,
33477841947SLeila Ghaffari                                     q[3][i] / rho
33577841947SLeila Ghaffari                                    };
33677841947SLeila Ghaffari     const CeedScalar E          =   q[4][i];
337e6225c47SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
338e6225c47SLeila Ghaffari                                     dq[1][0][i],
339e6225c47SLeila Ghaffari                                     dq[2][0][i]
340e6225c47SLeila Ghaffari                                    };
341e6225c47SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
342e6225c47SLeila Ghaffari                                     dq[1][1][i],
343e6225c47SLeila Ghaffari                                     dq[2][1][i]},
344e6225c47SLeila Ghaffari                                    {dq[0][2][i],
345e6225c47SLeila Ghaffari                                     dq[1][2][i],
346e6225c47SLeila Ghaffari                                     dq[2][2][i]},
347e6225c47SLeila Ghaffari                                    {dq[0][3][i],
348e6225c47SLeila Ghaffari                                     dq[1][3][i],
349e6225c47SLeila Ghaffari                                     dq[2][3][i]}
350e6225c47SLeila Ghaffari                                   };
351e6225c47SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
352e6225c47SLeila Ghaffari                                     dq[1][4][i],
353e6225c47SLeila Ghaffari                                     dq[2][4][i]
354e6225c47SLeila Ghaffari                                    };
35577841947SLeila Ghaffari     // -- Interp-to-Interp q_data
35677841947SLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
35777841947SLeila Ghaffari     // -- Interp-to-Grad q_data
35877841947SLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
35977841947SLeila Ghaffari     // *INDENT-OFF*
36077841947SLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
36177841947SLeila Ghaffari                                     q_data[2][i],
36277841947SLeila Ghaffari                                     q_data[3][i]},
36377841947SLeila Ghaffari                                    {q_data[4][i],
36477841947SLeila Ghaffari                                     q_data[5][i],
36577841947SLeila Ghaffari                                     q_data[6][i]},
36677841947SLeila Ghaffari                                    {q_data[7][i],
36777841947SLeila Ghaffari                                     q_data[8][i],
36877841947SLeila Ghaffari                                     q_data[9][i]}
36977841947SLeila Ghaffari                                   };
37077841947SLeila Ghaffari     // *INDENT-ON*
371e6225c47SLeila Ghaffari     // dU/dx
372e6225c47SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
373e6225c47SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
374e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
375e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
376e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
377e6225c47SLeila Ghaffari       for (int k=0; k<3; k++) {
378e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
379e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
380e6225c47SLeila Ghaffari         for (int l=0; l<3; l++) {
381e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
382e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
383e6225c47SLeila Ghaffari         }
384e6225c47SLeila Ghaffari       }
385e6225c47SLeila Ghaffari     }
386e6225c47SLeila Ghaffari     // Pressure
38777841947SLeila Ghaffari     const CeedScalar
38877841947SLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
38977841947SLeila Ghaffari     E_internal = E - E_kinetic,
390e6225c47SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
39177841947SLeila Ghaffari 
39277841947SLeila Ghaffari     // The Physics
39377841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
39477841947SLeila Ghaffari     for (int j=0; j<5; j++) {
395e6225c47SLeila Ghaffari       v[j][i] = 0.;
39677841947SLeila Ghaffari       for (int k=0; k<3; k++)
397e6225c47SLeila Ghaffari         dv[k][j][i] = 0.;
39877841947SLeila Ghaffari     }
39977841947SLeila Ghaffari 
40077841947SLeila Ghaffari     // -- Density
40177841947SLeila Ghaffari     // ---- u rho
40277841947SLeila Ghaffari     for (int j=0; j<3; j++)
40377841947SLeila Ghaffari       dv[j][0][i]  += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
40477841947SLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
40577841947SLeila Ghaffari     // -- Momentum
40677841947SLeila Ghaffari     // ---- rho (u x u) + P I3
40777841947SLeila Ghaffari     for (int j=0; j<3; j++)
40877841947SLeila Ghaffari       for (int k=0; k<3; k++)
409e6225c47SLeila Ghaffari         dv[k][j+1][i]  += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
410e6225c47SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
411e6225c47SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
41277841947SLeila Ghaffari     // -- Total Energy Density
41377841947SLeila Ghaffari     // ---- (E + P) u
41477841947SLeila Ghaffari     for (int j=0; j<3; j++)
41577841947SLeila Ghaffari       dv[j][4][i]  += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
41677841947SLeila Ghaffari                                          u[2]*dXdx[j][2]);
417e6225c47SLeila Ghaffari 
418e6225c47SLeila Ghaffari     // --Stabilization terms
419e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
420e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
421932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
422e6225c47SLeila Ghaffari 
423e6225c47SLeila Ghaffari     // ---- Transpose of the Jacobian
424e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv_T[3][5][5];
425e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
426e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
427e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
428e6225c47SLeila Ghaffari           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
429e6225c47SLeila Ghaffari 
430e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
431e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
432e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
433e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
434e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
435e6225c47SLeila Ghaffari       for (int k=0; k<3; k++)
436e6225c47SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
437e6225c47SLeila Ghaffari     }
438e6225c47SLeila Ghaffari 
439e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
440e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
441e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
442e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
443e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
444e6225c47SLeila Ghaffari           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
445e6225c47SLeila Ghaffari 
446932417b3SJed Brown     // Stabilization
447932417b3SJed Brown     // -- Tau elements
448932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
449932417b3SJed Brown     CeedScalar Tau_x[3] = {0.};
450932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
451e6225c47SLeila Ghaffari 
452932417b3SJed Brown     // -- Stabilization method: none or SU
453e6225c47SLeila Ghaffari     CeedScalar stab[5][3];
454e6225c47SLeila Ghaffari     switch (context->stabilization) {
455e6225c47SLeila Ghaffari     case 0:        // Galerkin
456e6225c47SLeila Ghaffari       break;
457e6225c47SLeila Ghaffari     case 1:        // SU
458e6225c47SLeila Ghaffari       for (int j=0; j<3; j++)
459e6225c47SLeila Ghaffari         for (int k=0; k<5; k++)
460e6225c47SLeila Ghaffari           for (int l=0; l<5; l++)
461932417b3SJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
462e6225c47SLeila Ghaffari 
463e6225c47SLeila Ghaffari       for (int j=0; j<5; j++)
464e6225c47SLeila Ghaffari         for (int k=0; k<3; k++)
465e6225c47SLeila Ghaffari           dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] +
466e6225c47SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
467e6225c47SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
468e6225c47SLeila Ghaffari       break;
469e6225c47SLeila Ghaffari     case 2:        // SUPG is not implemented for explicit scheme
470e6225c47SLeila Ghaffari       break;
471e6225c47SLeila Ghaffari     }
472e6225c47SLeila Ghaffari 
47377841947SLeila Ghaffari   } // End Quadrature Point Loop
47477841947SLeila Ghaffari 
47577841947SLeila Ghaffari   // Return
47677841947SLeila Ghaffari   return 0;
47777841947SLeila Ghaffari }
47877841947SLeila Ghaffari 
47977841947SLeila Ghaffari // *****************************************************************************
48077841947SLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above)
48177841947SLeila Ghaffari //   with implicit time stepping method
48277841947SLeila Ghaffari //
48377841947SLeila Ghaffari // *****************************************************************************
48477841947SLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q,
48577841947SLeila Ghaffari                                 const CeedScalar *const *in, CeedScalar *const *out) {
48677841947SLeila Ghaffari   // *INDENT-OFF*
48777841947SLeila Ghaffari   // Inputs
48877841947SLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
489e6225c47SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
49077841947SLeila Ghaffari                    (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
49177841947SLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
49277841947SLeila Ghaffari   // Outputs
49377841947SLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
49477841947SLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
49577841947SLeila Ghaffari 
496e6225c47SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
497932417b3SJed Brown   const CeedScalar c_tau = context->c_tau;
49877841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
49977841947SLeila Ghaffari 
50077841947SLeila Ghaffari   CeedPragmaSIMD
50177841947SLeila Ghaffari   // Quadrature Point Loop
50277841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
50377841947SLeila Ghaffari     // *INDENT-OFF*
50477841947SLeila Ghaffari     // Setup
50577841947SLeila Ghaffari     // -- Interp in
50677841947SLeila Ghaffari     const CeedScalar rho        =   q[0][i];
50777841947SLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
50877841947SLeila Ghaffari                                     q[2][i] / rho,
50977841947SLeila Ghaffari                                     q[3][i] / rho
51077841947SLeila Ghaffari                                    };
51177841947SLeila Ghaffari     const CeedScalar E          =   q[4][i];
512e6225c47SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
513e6225c47SLeila Ghaffari                                     dq[1][0][i],
514e6225c47SLeila Ghaffari                                     dq[2][0][i]
515e6225c47SLeila Ghaffari                                    };
516e6225c47SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
517e6225c47SLeila Ghaffari                                     dq[1][1][i],
518e6225c47SLeila Ghaffari                                     dq[2][1][i]},
519e6225c47SLeila Ghaffari                                    {dq[0][2][i],
520e6225c47SLeila Ghaffari                                     dq[1][2][i],
521e6225c47SLeila Ghaffari                                     dq[2][2][i]},
522e6225c47SLeila Ghaffari                                    {dq[0][3][i],
523e6225c47SLeila Ghaffari                                     dq[1][3][i],
524e6225c47SLeila Ghaffari                                     dq[2][3][i]}
525e6225c47SLeila Ghaffari                                   };
526e6225c47SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
527e6225c47SLeila Ghaffari                                     dq[1][4][i],
528e6225c47SLeila Ghaffari                                     dq[2][4][i]
529e6225c47SLeila Ghaffari                                    };
53077841947SLeila Ghaffari     // -- Interp-to-Interp q_data
53177841947SLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
53277841947SLeila Ghaffari     // -- Interp-to-Grad q_data
53377841947SLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
53477841947SLeila Ghaffari     // *INDENT-OFF*
53577841947SLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
53677841947SLeila Ghaffari                                     q_data[2][i],
53777841947SLeila Ghaffari                                     q_data[3][i]},
53877841947SLeila Ghaffari                                    {q_data[4][i],
53977841947SLeila Ghaffari                                     q_data[5][i],
54077841947SLeila Ghaffari                                     q_data[6][i]},
54177841947SLeila Ghaffari                                    {q_data[7][i],
54277841947SLeila Ghaffari                                     q_data[8][i],
54377841947SLeila Ghaffari                                     q_data[9][i]}
54477841947SLeila Ghaffari                                   };
54577841947SLeila Ghaffari     // *INDENT-ON*
546e6225c47SLeila Ghaffari     // dU/dx
547e6225c47SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
548e6225c47SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
549e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
550e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
551e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
552e6225c47SLeila Ghaffari       for (int k=0; k<3; k++) {
553e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
554e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
555e6225c47SLeila Ghaffari         for (int l=0; l<3; l++) {
556e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
557e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
558e6225c47SLeila Ghaffari         }
559e6225c47SLeila Ghaffari       }
560e6225c47SLeila Ghaffari     }
56177841947SLeila Ghaffari     const CeedScalar
56277841947SLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
56377841947SLeila Ghaffari     E_internal = E - E_kinetic,
564e6225c47SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
56577841947SLeila Ghaffari 
56677841947SLeila Ghaffari     // The Physics
56777841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
56877841947SLeila Ghaffari     for (int j=0; j<5; j++) {
569e6225c47SLeila Ghaffari       v[j][i] = 0.;
57077841947SLeila Ghaffari       for (int k=0; k<3; k++)
571e6225c47SLeila Ghaffari         dv[k][j][i] = 0.;
57277841947SLeila Ghaffari     }
57377841947SLeila Ghaffari     //-----mass matrix
57477841947SLeila Ghaffari     for (int j=0; j<5; j++)
57577841947SLeila Ghaffari       v[j][i] += wdetJ*q_dot[j][i];
57677841947SLeila Ghaffari 
57777841947SLeila Ghaffari     // -- Density
57877841947SLeila Ghaffari     // ---- u rho
57977841947SLeila Ghaffari     for (int j=0; j<3; j++)
58077841947SLeila Ghaffari       dv[j][0][i]  -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
58177841947SLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
58277841947SLeila Ghaffari     // -- Momentum
58377841947SLeila Ghaffari     // ---- rho (u x u) + P I3
58477841947SLeila Ghaffari     for (int j=0; j<3; j++)
58577841947SLeila Ghaffari       for (int k=0; k<3; k++)
586e6225c47SLeila Ghaffari         dv[k][j+1][i]  -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
587e6225c47SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
588e6225c47SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
58977841947SLeila Ghaffari     // -- Total Energy Density
59077841947SLeila Ghaffari     // ---- (E + P) u
59177841947SLeila Ghaffari     for (int j=0; j<3; j++)
59277841947SLeila Ghaffari       dv[j][4][i]  -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
59377841947SLeila Ghaffari                                          u[2]*dXdx[j][2]);
594e6225c47SLeila Ghaffari 
595e6225c47SLeila Ghaffari     // -- Stabilization terms
596e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
597e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
598932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
599e6225c47SLeila Ghaffari 
600e6225c47SLeila Ghaffari     // ---- Transpose of the Jacobian
601e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv_T[3][5][5];
602e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
603e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
604e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
605e6225c47SLeila Ghaffari           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
606e6225c47SLeila Ghaffari 
607e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
608e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
609e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
610e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
611e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
612e6225c47SLeila Ghaffari       for (int k=0; k<3; k++)
613e6225c47SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
614e6225c47SLeila Ghaffari     }
615e6225c47SLeila Ghaffari 
616e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
617e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
618e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
619e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
620e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
621e6225c47SLeila Ghaffari           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
622e6225c47SLeila Ghaffari 
623e6225c47SLeila Ghaffari     // ---- Strong residual
624e6225c47SLeila Ghaffari     CeedScalar strong_res[5];
625e6225c47SLeila Ghaffari     for (int j=0; j<5; j++)
626e6225c47SLeila Ghaffari       strong_res[j] = q_dot[j][i] + strong_conv[j];
627e6225c47SLeila Ghaffari 
628932417b3SJed Brown     // Stabilization
629932417b3SJed Brown     // -- Tau elements
630932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
631932417b3SJed Brown     CeedScalar Tau_x[3] = {0.};
632932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
633e6225c47SLeila Ghaffari 
634932417b3SJed Brown     // -- Stabilization method: none, SU, or SUPG
635e6225c47SLeila Ghaffari     CeedScalar stab[5][3];
636e6225c47SLeila Ghaffari     switch (context->stabilization) {
637e6225c47SLeila Ghaffari     case 0:        // Galerkin
638e6225c47SLeila Ghaffari       break;
639e6225c47SLeila Ghaffari     case 1:        // SU
640e6225c47SLeila Ghaffari       for (int j=0; j<3; j++)
641e6225c47SLeila Ghaffari         for (int k=0; k<5; k++)
642e6225c47SLeila Ghaffari           for (int l=0; l<5; l++)
643932417b3SJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
644e6225c47SLeila Ghaffari 
645e6225c47SLeila Ghaffari       for (int j=0; j<5; j++)
646e6225c47SLeila Ghaffari         for (int k=0; k<3; k++)
647e6225c47SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
648e6225c47SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
649e6225c47SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
650e6225c47SLeila Ghaffari       break;
651e6225c47SLeila Ghaffari     case 2:        // SUPG
652e6225c47SLeila Ghaffari       for (int j=0; j<3; j++)
653e6225c47SLeila Ghaffari         for (int k=0; k<5; k++)
654e6225c47SLeila Ghaffari           for (int l=0; l<5; l++)
655932417b3SJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_res[l];
656e6225c47SLeila Ghaffari 
657e6225c47SLeila Ghaffari       for (int j=0; j<5; j++)
658e6225c47SLeila Ghaffari         for (int k=0; k<3; k++)
659e6225c47SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
660e6225c47SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
661e6225c47SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
662e6225c47SLeila Ghaffari       break;
663e6225c47SLeila Ghaffari     }
66477841947SLeila Ghaffari   } // End Quadrature Point Loop
66577841947SLeila Ghaffari 
66677841947SLeila Ghaffari   // Return
66777841947SLeila Ghaffari   return 0;
66877841947SLeila Ghaffari }
66977841947SLeila Ghaffari // *****************************************************************************
670*2fe7aee7SLeila Ghaffari // This QFunction sets the inflow boundary conditions for
671*2fe7aee7SLeila Ghaffari //   the traveling vortex problem.
67277841947SLeila Ghaffari //
67377841947SLeila Ghaffari //  Prescribed T_inlet and P_inlet are converted to conservative variables
67477841947SLeila Ghaffari //      and applied weakly.
67577841947SLeila Ghaffari //
67677841947SLeila Ghaffari // *****************************************************************************
677*2fe7aee7SLeila Ghaffari CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q,
67877841947SLeila Ghaffari                                        const CeedScalar *const *in,
67977841947SLeila Ghaffari                                        CeedScalar *const *out) {
68077841947SLeila Ghaffari   // *INDENT-OFF*
68177841947SLeila Ghaffari   // Inputs
682*2fe7aee7SLeila Ghaffari   const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
68377841947SLeila Ghaffari   // Outputs
68477841947SLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
68577841947SLeila Ghaffari   // *INDENT-ON*
68677841947SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
68777841947SLeila Ghaffari   const int euler_test      = context->euler_test;
68877841947SLeila Ghaffari   const bool implicit       = context->implicit;
68977841947SLeila Ghaffari   CeedScalar *mean_velocity = context->mean_velocity;
69077841947SLeila Ghaffari   const CeedScalar cv    = 2.5;
69177841947SLeila Ghaffari   const CeedScalar R     = 1.;
69277841947SLeila Ghaffari   CeedScalar T_inlet;
69377841947SLeila Ghaffari   CeedScalar P_inlet;
69477841947SLeila Ghaffari 
69577841947SLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
69677841947SLeila Ghaffari   if (euler_test == 1 || euler_test == 3)
69777841947SLeila Ghaffari     for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.;
69877841947SLeila Ghaffari 
69977841947SLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
70077841947SLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
70177841947SLeila Ghaffari   else T_inlet = P_inlet = 1.;
70277841947SLeila Ghaffari 
70377841947SLeila Ghaffari   CeedPragmaSIMD
70477841947SLeila Ghaffari   // Quadrature Point Loop
70577841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
70677841947SLeila Ghaffari     // Setup
70777841947SLeila Ghaffari     // -- Interp-to-Interp q_data
70877841947SLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
70977841947SLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
71077841947SLeila Ghaffari     // We can effect this by swapping the sign on this weight
71177841947SLeila Ghaffari     const CeedScalar wdetJb  = (implicit ? -1. : 1.) * q_data_sur[0][i];
712*2fe7aee7SLeila Ghaffari     // ---- Normal vect
71377841947SLeila Ghaffari     const CeedScalar norm[3] = {q_data_sur[1][i],
71477841947SLeila Ghaffari                                 q_data_sur[2][i],
71577841947SLeila Ghaffari                                 q_data_sur[3][i]
71677841947SLeila Ghaffari                                };
71777841947SLeila Ghaffari 
71877841947SLeila Ghaffari     // face_normal = Normal vector of the face
71977841947SLeila Ghaffari     const CeedScalar face_normal = norm[0]*mean_velocity[0] +
72077841947SLeila Ghaffari                                    norm[1]*mean_velocity[1] +
72177841947SLeila Ghaffari                                    norm[2]*mean_velocity[2];
72277841947SLeila Ghaffari     // The Physics
72377841947SLeila Ghaffari     // Zero v so all future terms can safely sum into it
724e6225c47SLeila Ghaffari     for (int j=0; j<5; j++) v[j][i] = 0.;
72577841947SLeila Ghaffari 
72677841947SLeila Ghaffari     // Implementing in/outflow BCs
727*2fe7aee7SLeila Ghaffari     if (face_normal > 0) {
72877841947SLeila Ghaffari     } else { // inflow
72977841947SLeila Ghaffari       const CeedScalar rho_inlet = P_inlet/(R*T_inlet);
73077841947SLeila Ghaffari       const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] +
73177841947SLeila Ghaffari                                           mean_velocity[1]*mean_velocity[1]) / 2.;
73277841947SLeila Ghaffari       // incoming total energy
73377841947SLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
73477841947SLeila Ghaffari 
73577841947SLeila Ghaffari       // The Physics
73677841947SLeila Ghaffari       // -- Density
73777841947SLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
73877841947SLeila Ghaffari 
73977841947SLeila Ghaffari       // -- Momentum
74077841947SLeila Ghaffari       for (int j=0; j<3; j++)
74177841947SLeila Ghaffari         v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] +
74277841947SLeila Ghaffari                               norm[j] * P_inlet);
74377841947SLeila Ghaffari 
74477841947SLeila Ghaffari       // -- Total Energy Density
74577841947SLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
74677841947SLeila Ghaffari     }
74777841947SLeila Ghaffari 
74877841947SLeila Ghaffari   } // End Quadrature Point Loop
74977841947SLeila Ghaffari   return 0;
75077841947SLeila Ghaffari }
75177841947SLeila Ghaffari 
75277841947SLeila Ghaffari // *****************************************************************************
75377841947SLeila Ghaffari 
75477841947SLeila Ghaffari #endif // eulervortex_h
755