xref: /libCEED/examples/fluids/qfunctions/eulervortex.h (revision 932417b3091c5f226feee63990ab54384e85e3d6)
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;
41*932417b3SJed 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;
18277841947SLeila Ghaffari   }
18377841947SLeila Ghaffari   // Return
18477841947SLeila Ghaffari   return 0;
18577841947SLeila Ghaffari }
18677841947SLeila Ghaffari 
18777841947SLeila Ghaffari // *****************************************************************************
188e6225c47SLeila Ghaffari // Helper function for computing flux Jacobian
189e6225c47SLeila Ghaffari // *****************************************************************************
190*932417b3SJed Brown CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5],
191e6225c47SLeila Ghaffari     const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
192e6225c47SLeila Ghaffari     const CeedScalar gamma) {
193e6225c47SLeila Ghaffari   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square
194e6225c47SLeila Ghaffari   for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions
195e6225c47SLeila Ghaffari     for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix
196e6225c47SLeila Ghaffari       dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j];
197e6225c47SLeila Ghaffari       for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix
198e6225c47SLeila Ghaffari         dF[i][0][k+1]   = ((i==k) ? 1. : 0.);
199e6225c47SLeila Ghaffari         dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) +
200e6225c47SLeila Ghaffari                           ((i==k) ? u[j] : 0.) -
201e6225c47SLeila Ghaffari                           ((i==j) ? u[k] : 0.) * (gamma-1.);
202e6225c47SLeila Ghaffari         dF[i][4][k+1]   = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) -
203e6225c47SLeila Ghaffari                           (gamma-1.)*u[i]*u[k];
204e6225c47SLeila Ghaffari       }
205e6225c47SLeila Ghaffari       dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.);
206e6225c47SLeila Ghaffari     }
207e6225c47SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho);
208e6225c47SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
209e6225c47SLeila Ghaffari   }
210e6225c47SLeila Ghaffari }
211e6225c47SLeila Ghaffari 
212e6225c47SLeila Ghaffari // *****************************************************************************
213*932417b3SJed Brown // Helper function for computing Tau elements (stabilization constant)
214*932417b3SJed Brown //   Model from:
215*932417b3SJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
216*932417b3SJed Brown //
217*932417b3SJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
218*932417b3SJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
219*932417b3SJed Brown //
220*932417b3SJed Brown // Where
221*932417b3SJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
222*932417b3SJed Brown //   h[i]      = 2 length(dxdX[i])
223*932417b3SJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
224*932417b3SJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
225*932417b3SJed Brown //   rho(A[i]) = spectral radius of the convective flux Jacobian i,
226*932417b3SJed Brown //               wave speed in direction i
227*932417b3SJed Brown // *****************************************************************************
228*932417b3SJed Brown CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3],
229*932417b3SJed Brown                                        const CeedScalar dXdx[3][3], const CeedScalar u[3],
230*932417b3SJed Brown                                        const CeedScalar sound_speed, const CeedScalar c_tau) {
231*932417b3SJed Brown   for (int i=0; i<3; i++) {
232*932417b3SJed Brown     // length of element in direction i
233*932417b3SJed Brown     CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] +
234*932417b3SJed Brown                             dXdx[2][i]*dXdx[2][i]);
235*932417b3SJed Brown     // fastest wave in direction i
236*932417b3SJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
237*932417b3SJed Brown     Tau_x[i] = c_tau * h / fastest_wave;
238*932417b3SJed Brown   }
239*932417b3SJed Brown }
240*932417b3SJed Brown 
241*932417b3SJed Brown // *****************************************************************************
24277841947SLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
24377841947SLeila Ghaffari // *****************************************************************************
24477841947SLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q,
24577841947SLeila Ghaffari                          const CeedScalar *const *in, CeedScalar *const *out) {
24677841947SLeila Ghaffari   // Inputs
24777841947SLeila Ghaffari   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
24877841947SLeila Ghaffari 
24977841947SLeila Ghaffari   // Outputs
25077841947SLeila Ghaffari   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
25177841947SLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
25277841947SLeila Ghaffari 
25377841947SLeila Ghaffari   CeedPragmaSIMD
25477841947SLeila Ghaffari   // Quadrature Point Loop
25577841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
25677841947SLeila Ghaffari     const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
257e6225c47SLeila Ghaffari     CeedScalar q[5] = {0.};
25877841947SLeila Ghaffari 
25977841947SLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
26077841947SLeila Ghaffari 
26177841947SLeila Ghaffari     for (CeedInt j=0; j<5; j++)
26277841947SLeila Ghaffari       q0[j][i] = q[j];
26377841947SLeila Ghaffari   } // End of Quadrature Point Loop
26477841947SLeila Ghaffari 
26577841947SLeila Ghaffari   // Return
26677841947SLeila Ghaffari   return 0;
26777841947SLeila Ghaffari }
26877841947SLeila Ghaffari 
26977841947SLeila Ghaffari // *****************************************************************************
27077841947SLeila Ghaffari // This QFunction implements the following formulation of Euler equations
27177841947SLeila Ghaffari //   with explicit time stepping method
27277841947SLeila Ghaffari //
27377841947SLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation
27477841947SLeila Ghaffari //   form with state variables of density, momentum density, and total
27577841947SLeila Ghaffari //   energy density.
27677841947SLeila Ghaffari //
27777841947SLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
27877841947SLeila Ghaffari //   rho - Mass Density
27977841947SLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
28077841947SLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
28177841947SLeila Ghaffari //
28277841947SLeila Ghaffari // Euler Equations:
28377841947SLeila Ghaffari //   drho/dt + div( U )                   = 0
28477841947SLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
28577841947SLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
28677841947SLeila Ghaffari //
28777841947SLeila Ghaffari // Equation of State:
28877841947SLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
28977841947SLeila Ghaffari //
29077841947SLeila Ghaffari // Constants:
29177841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
29277841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
29377841947SLeila Ghaffari //   g               ,  Gravity
29477841947SLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
29577841947SLeila Ghaffari // *****************************************************************************
29677841947SLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q,
29777841947SLeila Ghaffari                       const CeedScalar *const *in, CeedScalar *const *out) {
29877841947SLeila Ghaffari   // *INDENT-OFF*
29977841947SLeila Ghaffari   // Inputs
30077841947SLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
301e6225c47SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
30277841947SLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
30377841947SLeila Ghaffari   // Outputs
30477841947SLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
30577841947SLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
30677841947SLeila Ghaffari 
307e6225c47SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
308*932417b3SJed Brown   const CeedScalar c_tau = context->c_tau;
30977841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
31077841947SLeila Ghaffari 
31177841947SLeila Ghaffari   CeedPragmaSIMD
31277841947SLeila Ghaffari   // Quadrature Point Loop
31377841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
31477841947SLeila Ghaffari     // *INDENT-OFF*
31577841947SLeila Ghaffari     // Setup
31677841947SLeila Ghaffari     // -- Interp in
31777841947SLeila Ghaffari     const CeedScalar rho        =   q[0][i];
31877841947SLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
31977841947SLeila Ghaffari                                     q[2][i] / rho,
32077841947SLeila Ghaffari                                     q[3][i] / rho
32177841947SLeila Ghaffari                                    };
32277841947SLeila Ghaffari     const CeedScalar E          =   q[4][i];
323e6225c47SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
324e6225c47SLeila Ghaffari                                     dq[1][0][i],
325e6225c47SLeila Ghaffari                                     dq[2][0][i]
326e6225c47SLeila Ghaffari                                    };
327e6225c47SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
328e6225c47SLeila Ghaffari                                     dq[1][1][i],
329e6225c47SLeila Ghaffari                                     dq[2][1][i]},
330e6225c47SLeila Ghaffari                                    {dq[0][2][i],
331e6225c47SLeila Ghaffari                                     dq[1][2][i],
332e6225c47SLeila Ghaffari                                     dq[2][2][i]},
333e6225c47SLeila Ghaffari                                    {dq[0][3][i],
334e6225c47SLeila Ghaffari                                     dq[1][3][i],
335e6225c47SLeila Ghaffari                                     dq[2][3][i]}
336e6225c47SLeila Ghaffari                                   };
337e6225c47SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
338e6225c47SLeila Ghaffari                                     dq[1][4][i],
339e6225c47SLeila Ghaffari                                     dq[2][4][i]
340e6225c47SLeila Ghaffari                                    };
34177841947SLeila Ghaffari     // -- Interp-to-Interp q_data
34277841947SLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
34377841947SLeila Ghaffari     // -- Interp-to-Grad q_data
34477841947SLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
34577841947SLeila Ghaffari     // *INDENT-OFF*
34677841947SLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
34777841947SLeila Ghaffari                                     q_data[2][i],
34877841947SLeila Ghaffari                                     q_data[3][i]},
34977841947SLeila Ghaffari                                    {q_data[4][i],
35077841947SLeila Ghaffari                                     q_data[5][i],
35177841947SLeila Ghaffari                                     q_data[6][i]},
35277841947SLeila Ghaffari                                    {q_data[7][i],
35377841947SLeila Ghaffari                                     q_data[8][i],
35477841947SLeila Ghaffari                                     q_data[9][i]}
35577841947SLeila Ghaffari                                   };
35677841947SLeila Ghaffari     // *INDENT-ON*
357e6225c47SLeila Ghaffari     // dU/dx
358e6225c47SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
359e6225c47SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
360e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
361e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
362e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
363e6225c47SLeila Ghaffari       for (int k=0; k<3; k++) {
364e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
365e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
366e6225c47SLeila Ghaffari         for (int l=0; l<3; l++) {
367e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
368e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
369e6225c47SLeila Ghaffari         }
370e6225c47SLeila Ghaffari       }
371e6225c47SLeila Ghaffari     }
372e6225c47SLeila Ghaffari     // Pressure
37377841947SLeila Ghaffari     const CeedScalar
37477841947SLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
37577841947SLeila Ghaffari     E_internal = E - E_kinetic,
376e6225c47SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
37777841947SLeila Ghaffari 
37877841947SLeila Ghaffari     // The Physics
37977841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
38077841947SLeila Ghaffari     for (int j=0; j<5; j++) {
381e6225c47SLeila Ghaffari       v[j][i] = 0.;
38277841947SLeila Ghaffari       for (int k=0; k<3; k++)
383e6225c47SLeila Ghaffari         dv[k][j][i] = 0.;
38477841947SLeila Ghaffari     }
38577841947SLeila Ghaffari 
38677841947SLeila Ghaffari     // -- Density
38777841947SLeila Ghaffari     // ---- u rho
38877841947SLeila Ghaffari     for (int j=0; j<3; j++)
38977841947SLeila Ghaffari       dv[j][0][i]  += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
39077841947SLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
39177841947SLeila Ghaffari     // -- Momentum
39277841947SLeila Ghaffari     // ---- rho (u x u) + P I3
39377841947SLeila Ghaffari     for (int j=0; j<3; j++)
39477841947SLeila Ghaffari       for (int k=0; k<3; k++)
395e6225c47SLeila Ghaffari         dv[k][j+1][i]  += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
396e6225c47SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
397e6225c47SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
39877841947SLeila Ghaffari     // -- Total Energy Density
39977841947SLeila Ghaffari     // ---- (E + P) u
40077841947SLeila Ghaffari     for (int j=0; j<3; j++)
40177841947SLeila Ghaffari       dv[j][4][i]  += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
40277841947SLeila Ghaffari                                          u[2]*dXdx[j][2]);
403e6225c47SLeila Ghaffari 
404e6225c47SLeila Ghaffari     // --Stabilization terms
405e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
406e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
407*932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
408e6225c47SLeila Ghaffari 
409e6225c47SLeila Ghaffari     // ---- Transpose of the Jacobian
410e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv_T[3][5][5];
411e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
412e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
413e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
414e6225c47SLeila Ghaffari           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
415e6225c47SLeila Ghaffari 
416e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
417e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
418e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
419e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
420e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
421e6225c47SLeila Ghaffari       for (int k=0; k<3; k++)
422e6225c47SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
423e6225c47SLeila Ghaffari     }
424e6225c47SLeila Ghaffari 
425e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
426e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
427e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
428e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
429e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
430e6225c47SLeila Ghaffari           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
431e6225c47SLeila Ghaffari 
432*932417b3SJed Brown     // Stabilization
433*932417b3SJed Brown     // -- Tau elements
434*932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
435*932417b3SJed Brown     CeedScalar Tau_x[3] = {0.};
436*932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
437e6225c47SLeila Ghaffari 
438*932417b3SJed Brown     // -- Stabilization method: none or SU
439e6225c47SLeila Ghaffari     CeedScalar stab[5][3];
440e6225c47SLeila Ghaffari     switch (context->stabilization) {
441e6225c47SLeila Ghaffari     case 0:        // Galerkin
442e6225c47SLeila Ghaffari       break;
443e6225c47SLeila Ghaffari     case 1:        // SU
444e6225c47SLeila Ghaffari       for (int j=0; j<3; j++)
445e6225c47SLeila Ghaffari         for (int k=0; k<5; k++)
446e6225c47SLeila Ghaffari           for (int l=0; l<5; l++)
447*932417b3SJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
448e6225c47SLeila Ghaffari 
449e6225c47SLeila Ghaffari       for (int j=0; j<5; j++)
450e6225c47SLeila Ghaffari         for (int k=0; k<3; k++)
451e6225c47SLeila Ghaffari           dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] +
452e6225c47SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
453e6225c47SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
454e6225c47SLeila Ghaffari       break;
455e6225c47SLeila Ghaffari     case 2:        // SUPG is not implemented for explicit scheme
456e6225c47SLeila Ghaffari       break;
457e6225c47SLeila Ghaffari     }
458e6225c47SLeila Ghaffari 
45977841947SLeila Ghaffari   } // End Quadrature Point Loop
46077841947SLeila Ghaffari 
46177841947SLeila Ghaffari   // Return
46277841947SLeila Ghaffari   return 0;
46377841947SLeila Ghaffari }
46477841947SLeila Ghaffari 
46577841947SLeila Ghaffari // *****************************************************************************
46677841947SLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above)
46777841947SLeila Ghaffari //   with implicit time stepping method
46877841947SLeila Ghaffari //
46977841947SLeila Ghaffari // *****************************************************************************
47077841947SLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q,
47177841947SLeila Ghaffari                                 const CeedScalar *const *in, CeedScalar *const *out) {
47277841947SLeila Ghaffari   // *INDENT-OFF*
47377841947SLeila Ghaffari   // Inputs
47477841947SLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
475e6225c47SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
47677841947SLeila Ghaffari                    (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
47777841947SLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
47877841947SLeila Ghaffari   // Outputs
47977841947SLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
48077841947SLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
48177841947SLeila Ghaffari 
482e6225c47SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
483*932417b3SJed Brown   const CeedScalar c_tau = context->c_tau;
48477841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
48577841947SLeila Ghaffari 
48677841947SLeila Ghaffari   CeedPragmaSIMD
48777841947SLeila Ghaffari   // Quadrature Point Loop
48877841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
48977841947SLeila Ghaffari     // *INDENT-OFF*
49077841947SLeila Ghaffari     // Setup
49177841947SLeila Ghaffari     // -- Interp in
49277841947SLeila Ghaffari     const CeedScalar rho        =   q[0][i];
49377841947SLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
49477841947SLeila Ghaffari                                     q[2][i] / rho,
49577841947SLeila Ghaffari                                     q[3][i] / rho
49677841947SLeila Ghaffari                                    };
49777841947SLeila Ghaffari     const CeedScalar E          =   q[4][i];
498e6225c47SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
499e6225c47SLeila Ghaffari                                     dq[1][0][i],
500e6225c47SLeila Ghaffari                                     dq[2][0][i]
501e6225c47SLeila Ghaffari                                    };
502e6225c47SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
503e6225c47SLeila Ghaffari                                     dq[1][1][i],
504e6225c47SLeila Ghaffari                                     dq[2][1][i]},
505e6225c47SLeila Ghaffari                                    {dq[0][2][i],
506e6225c47SLeila Ghaffari                                     dq[1][2][i],
507e6225c47SLeila Ghaffari                                     dq[2][2][i]},
508e6225c47SLeila Ghaffari                                    {dq[0][3][i],
509e6225c47SLeila Ghaffari                                     dq[1][3][i],
510e6225c47SLeila Ghaffari                                     dq[2][3][i]}
511e6225c47SLeila Ghaffari                                   };
512e6225c47SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
513e6225c47SLeila Ghaffari                                     dq[1][4][i],
514e6225c47SLeila Ghaffari                                     dq[2][4][i]
515e6225c47SLeila Ghaffari                                    };
51677841947SLeila Ghaffari     // -- Interp-to-Interp q_data
51777841947SLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
51877841947SLeila Ghaffari     // -- Interp-to-Grad q_data
51977841947SLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
52077841947SLeila Ghaffari     // *INDENT-OFF*
52177841947SLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
52277841947SLeila Ghaffari                                     q_data[2][i],
52377841947SLeila Ghaffari                                     q_data[3][i]},
52477841947SLeila Ghaffari                                    {q_data[4][i],
52577841947SLeila Ghaffari                                     q_data[5][i],
52677841947SLeila Ghaffari                                     q_data[6][i]},
52777841947SLeila Ghaffari                                    {q_data[7][i],
52877841947SLeila Ghaffari                                     q_data[8][i],
52977841947SLeila Ghaffari                                     q_data[9][i]}
53077841947SLeila Ghaffari                                   };
53177841947SLeila Ghaffari     // *INDENT-ON*
532e6225c47SLeila Ghaffari     // dU/dx
533e6225c47SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
534e6225c47SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
535e6225c47SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
536e6225c47SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
537e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
538e6225c47SLeila Ghaffari       for (int k=0; k<3; k++) {
539e6225c47SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
540e6225c47SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
541e6225c47SLeila Ghaffari         for (int l=0; l<3; l++) {
542e6225c47SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
543e6225c47SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
544e6225c47SLeila Ghaffari         }
545e6225c47SLeila Ghaffari       }
546e6225c47SLeila Ghaffari     }
54777841947SLeila Ghaffari     const CeedScalar
54877841947SLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
54977841947SLeila Ghaffari     E_internal = E - E_kinetic,
550e6225c47SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
55177841947SLeila Ghaffari 
55277841947SLeila Ghaffari     // The Physics
55377841947SLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
55477841947SLeila Ghaffari     for (int j=0; j<5; j++) {
555e6225c47SLeila Ghaffari       v[j][i] = 0.;
55677841947SLeila Ghaffari       for (int k=0; k<3; k++)
557e6225c47SLeila Ghaffari         dv[k][j][i] = 0.;
55877841947SLeila Ghaffari     }
55977841947SLeila Ghaffari     //-----mass matrix
56077841947SLeila Ghaffari     for (int j=0; j<5; j++)
56177841947SLeila Ghaffari       v[j][i] += wdetJ*q_dot[j][i];
56277841947SLeila Ghaffari 
56377841947SLeila Ghaffari     // -- Density
56477841947SLeila Ghaffari     // ---- u rho
56577841947SLeila Ghaffari     for (int j=0; j<3; j++)
56677841947SLeila Ghaffari       dv[j][0][i]  -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
56777841947SLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
56877841947SLeila Ghaffari     // -- Momentum
56977841947SLeila Ghaffari     // ---- rho (u x u) + P I3
57077841947SLeila Ghaffari     for (int j=0; j<3; j++)
57177841947SLeila Ghaffari       for (int k=0; k<3; k++)
572e6225c47SLeila Ghaffari         dv[k][j+1][i]  -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
573e6225c47SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
574e6225c47SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
57577841947SLeila Ghaffari     // -- Total Energy Density
57677841947SLeila Ghaffari     // ---- (E + P) u
57777841947SLeila Ghaffari     for (int j=0; j<3; j++)
57877841947SLeila Ghaffari       dv[j][4][i]  -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
57977841947SLeila Ghaffari                                          u[2]*dXdx[j][2]);
580e6225c47SLeila Ghaffari 
581e6225c47SLeila Ghaffari     // -- Stabilization terms
582e6225c47SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
583e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
584*932417b3SJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
585e6225c47SLeila Ghaffari 
586e6225c47SLeila Ghaffari     // ---- Transpose of the Jacobian
587e6225c47SLeila Ghaffari     CeedScalar jacob_F_conv_T[3][5][5];
588e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
589e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
590e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
591e6225c47SLeila Ghaffari           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
592e6225c47SLeila Ghaffari 
593e6225c47SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
594e6225c47SLeila Ghaffari     CeedScalar dqdx[5][3];
595e6225c47SLeila Ghaffari     for (int j=0; j<3; j++) {
596e6225c47SLeila Ghaffari       dqdx[0][j] = drhodx[j];
597e6225c47SLeila Ghaffari       dqdx[4][j] = dEdx[j];
598e6225c47SLeila Ghaffari       for (int k=0; k<3; k++)
599e6225c47SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
600e6225c47SLeila Ghaffari     }
601e6225c47SLeila Ghaffari 
602e6225c47SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
603e6225c47SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
604e6225c47SLeila Ghaffari     for (int j=0; j<3; j++)
605e6225c47SLeila Ghaffari       for (int k=0; k<5; k++)
606e6225c47SLeila Ghaffari         for (int l=0; l<5; l++)
607e6225c47SLeila Ghaffari           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
608e6225c47SLeila Ghaffari 
609e6225c47SLeila Ghaffari     // ---- Strong residual
610e6225c47SLeila Ghaffari     CeedScalar strong_res[5];
611e6225c47SLeila Ghaffari     for (int j=0; j<5; j++)
612e6225c47SLeila Ghaffari       strong_res[j] = q_dot[j][i] + strong_conv[j];
613e6225c47SLeila Ghaffari 
614*932417b3SJed Brown     // Stabilization
615*932417b3SJed Brown     // -- Tau elements
616*932417b3SJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
617*932417b3SJed Brown     CeedScalar Tau_x[3] = {0.};
618*932417b3SJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
619e6225c47SLeila Ghaffari 
620*932417b3SJed Brown     // -- Stabilization method: none, SU, or SUPG
621e6225c47SLeila Ghaffari     CeedScalar stab[5][3];
622e6225c47SLeila Ghaffari     switch (context->stabilization) {
623e6225c47SLeila Ghaffari     case 0:        // Galerkin
624e6225c47SLeila Ghaffari       break;
625e6225c47SLeila Ghaffari     case 1:        // SU
626e6225c47SLeila Ghaffari       for (int j=0; j<3; j++)
627e6225c47SLeila Ghaffari         for (int k=0; k<5; k++)
628e6225c47SLeila Ghaffari           for (int l=0; l<5; l++)
629*932417b3SJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
630e6225c47SLeila Ghaffari 
631e6225c47SLeila Ghaffari       for (int j=0; j<5; j++)
632e6225c47SLeila Ghaffari         for (int k=0; k<3; k++)
633e6225c47SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
634e6225c47SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
635e6225c47SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
636e6225c47SLeila Ghaffari       break;
637e6225c47SLeila Ghaffari     case 2:        // SUPG
638e6225c47SLeila Ghaffari       for (int j=0; j<3; j++)
639e6225c47SLeila Ghaffari         for (int k=0; k<5; k++)
640e6225c47SLeila Ghaffari           for (int l=0; l<5; l++)
641*932417b3SJed Brown             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_res[l];
642e6225c47SLeila Ghaffari 
643e6225c47SLeila Ghaffari       for (int j=0; j<5; j++)
644e6225c47SLeila Ghaffari         for (int k=0; k<3; k++)
645e6225c47SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
646e6225c47SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
647e6225c47SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
648e6225c47SLeila Ghaffari       break;
649e6225c47SLeila Ghaffari     }
65077841947SLeila Ghaffari   } // End Quadrature Point Loop
65177841947SLeila Ghaffari 
65277841947SLeila Ghaffari   // Return
65377841947SLeila Ghaffari   return 0;
65477841947SLeila Ghaffari }
65577841947SLeila Ghaffari // *****************************************************************************
65677841947SLeila Ghaffari // This QFunction sets the boundary conditions
65777841947SLeila Ghaffari //   In this problem, only in/outflow BCs are implemented
65877841947SLeila Ghaffari //
65977841947SLeila Ghaffari //  Inflow and outflow faces are determined based on
66077841947SLeila Ghaffari //    sign(dot(mean_velocity, normal)):
66177841947SLeila Ghaffari //      sign(dot(mean_velocity, normal)) > 0 : outflow BCs
66277841947SLeila Ghaffari //      sign(dot(mean_velocity, normal)) < 0 : inflow BCs
66377841947SLeila Ghaffari //
66477841947SLeila Ghaffari //  Outflow BCs:
66577841947SLeila Ghaffari //    The validity of the weak form of the governing equations is
66677841947SLeila Ghaffari //      extended to the outflow.
66777841947SLeila Ghaffari //
66877841947SLeila Ghaffari //  Inflow BCs:
66977841947SLeila Ghaffari //    Prescribed T_inlet and P_inlet are converted to conservative variables
67077841947SLeila Ghaffari //      and applied weakly.
67177841947SLeila Ghaffari //
67277841947SLeila Ghaffari // *****************************************************************************
67377841947SLeila Ghaffari CEED_QFUNCTION(Euler_Sur)(void *ctx, CeedInt Q,
67477841947SLeila Ghaffari                           const CeedScalar *const *in,
67577841947SLeila Ghaffari                           CeedScalar *const *out) {
67677841947SLeila Ghaffari   // *INDENT-OFF*
67777841947SLeila Ghaffari   // Inputs
67877841947SLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
67977841947SLeila Ghaffari                    (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
68077841947SLeila Ghaffari   // Outputs
68177841947SLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
68277841947SLeila Ghaffari   // *INDENT-ON*
68377841947SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
68477841947SLeila Ghaffari   const int euler_test      = context->euler_test;
68577841947SLeila Ghaffari   const bool implicit       = context->implicit;
68677841947SLeila Ghaffari   CeedScalar *mean_velocity = context->mean_velocity;
68777841947SLeila Ghaffari 
68877841947SLeila Ghaffari   const CeedScalar gamma = 1.4;
68977841947SLeila Ghaffari   const CeedScalar cv    = 2.5;
69077841947SLeila Ghaffari   const CeedScalar R     = 1.;
69177841947SLeila Ghaffari   CeedScalar T_inlet;
69277841947SLeila Ghaffari   CeedScalar P_inlet;
69377841947SLeila Ghaffari 
69477841947SLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
69577841947SLeila Ghaffari   if (euler_test == 1 || euler_test == 3)
69677841947SLeila Ghaffari     for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.;
69777841947SLeila Ghaffari 
69877841947SLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
69977841947SLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
70077841947SLeila Ghaffari   else T_inlet = P_inlet = 1.;
70177841947SLeila Ghaffari 
70277841947SLeila Ghaffari   CeedPragmaSIMD
70377841947SLeila Ghaffari   // Quadrature Point Loop
70477841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
70577841947SLeila Ghaffari     // Setup
70677841947SLeila Ghaffari     // -- Interp in
70777841947SLeila Ghaffari     const CeedScalar rho      =  q[0][i];
70877841947SLeila Ghaffari     const CeedScalar u[3]     = {q[1][i] / rho,
70977841947SLeila Ghaffari                                  q[2][i] / rho,
71077841947SLeila Ghaffari                                  q[3][i] / rho
71177841947SLeila Ghaffari                                 };
71277841947SLeila Ghaffari     const CeedScalar E        =  q[4][i];
71377841947SLeila Ghaffari 
71477841947SLeila Ghaffari     // -- Interp-to-Interp q_data
71577841947SLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
71677841947SLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
71777841947SLeila Ghaffari     // We can effect this by swapping the sign on this weight
71877841947SLeila Ghaffari     const CeedScalar wdetJb     =   (implicit ? -1. : 1.) * q_data_sur[0][i];
71977841947SLeila Ghaffari     // ---- Normal vectors
72077841947SLeila Ghaffari     const CeedScalar norm[3]    =   {q_data_sur[1][i],
72177841947SLeila Ghaffari                                      q_data_sur[2][i],
72277841947SLeila Ghaffari                                      q_data_sur[3][i]
72377841947SLeila Ghaffari                                     };
72477841947SLeila Ghaffari 
72577841947SLeila Ghaffari     // face_normal = Normal vector of the face
72677841947SLeila Ghaffari     const CeedScalar face_normal = norm[0]*mean_velocity[0] +
72777841947SLeila Ghaffari                                    norm[1]*mean_velocity[1] +
72877841947SLeila Ghaffari                                    norm[2]*mean_velocity[2];
72977841947SLeila Ghaffari     // The Physics
73077841947SLeila Ghaffari     // Zero v so all future terms can safely sum into it
731e6225c47SLeila Ghaffari     for (int j=0; j<5; j++) v[j][i] = 0.;
73277841947SLeila Ghaffari 
73377841947SLeila Ghaffari     // Implementing in/outflow BCs
73477841947SLeila Ghaffari     if (face_normal > 0) { // outflow
73577841947SLeila Ghaffari       const CeedScalar E_kinetic = (u[0]*u[0] + u[1]*u[1]) / 2.;
73677841947SLeila Ghaffari       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.); // pressure
73777841947SLeila Ghaffari       const CeedScalar u_normal  = norm[0]*u[0] + norm[1]*u[1] +
73877841947SLeila Ghaffari                                    norm[2]*u[2]; // Normal velocity
73977841947SLeila Ghaffari       // The Physics
74077841947SLeila Ghaffari       // -- Density
74177841947SLeila Ghaffari       v[0][i] -= wdetJb * rho * u_normal;
74277841947SLeila Ghaffari 
74377841947SLeila Ghaffari       // -- Momentum
74477841947SLeila Ghaffari       for (int j=0; j<3; j++)
74577841947SLeila Ghaffari         v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + norm[j] * P);
74677841947SLeila Ghaffari 
74777841947SLeila Ghaffari       // -- Total Energy Density
74877841947SLeila Ghaffari       v[4][i] -= wdetJb * u_normal * (E + P);
74977841947SLeila Ghaffari 
75077841947SLeila Ghaffari     } else { // inflow
75177841947SLeila Ghaffari       const CeedScalar rho_inlet = P_inlet/(R*T_inlet);
75277841947SLeila Ghaffari       const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] +
75377841947SLeila Ghaffari                                           mean_velocity[1]*mean_velocity[1]) / 2.;
75477841947SLeila Ghaffari       // incoming total energy
75577841947SLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
75677841947SLeila Ghaffari 
75777841947SLeila Ghaffari       // The Physics
75877841947SLeila Ghaffari       // -- Density
75977841947SLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
76077841947SLeila Ghaffari 
76177841947SLeila Ghaffari       // -- Momentum
76277841947SLeila Ghaffari       for (int j=0; j<3; j++)
76377841947SLeila Ghaffari         v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] +
76477841947SLeila Ghaffari                               norm[j] * P_inlet);
76577841947SLeila Ghaffari 
76677841947SLeila Ghaffari       // -- Total Energy Density
76777841947SLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
76877841947SLeila Ghaffari     }
76977841947SLeila Ghaffari 
77077841947SLeila Ghaffari   } // End Quadrature Point Loop
77177841947SLeila Ghaffari   return 0;
77277841947SLeila Ghaffari }
77377841947SLeila Ghaffari 
77477841947SLeila Ghaffari // *****************************************************************************
77577841947SLeila Ghaffari 
77677841947SLeila Ghaffari #endif // eulervortex_h
777