xref: /libCEED/examples/fluids/qfunctions/newtonian.h (revision 3d8e882215d238700cdceb37404f76ca7fa24eaa)
1*3d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*3d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
388b783a1SJames Wright //
4*3d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
588b783a1SJames Wright //
6*3d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
788b783a1SJames Wright 
888b783a1SJames Wright /// @file
988b783a1SJames Wright /// Operator for Navier-Stokes example using PETSc
1088b783a1SJames Wright 
1188b783a1SJames Wright 
1288b783a1SJames Wright #ifndef newtonian_h
1388b783a1SJames Wright #define newtonian_h
1488b783a1SJames Wright 
1588b783a1SJames Wright #include <math.h>
1688b783a1SJames Wright #include <ceed.h>
1788b783a1SJames Wright 
1888b783a1SJames Wright #ifndef M_PI
1988b783a1SJames Wright #define M_PI    3.14159265358979323846
2088b783a1SJames Wright #endif
2188b783a1SJames Wright 
2288b783a1SJames Wright #ifndef setup_context_struct
2388b783a1SJames Wright #define setup_context_struct
2488b783a1SJames Wright typedef struct SetupContext_ *SetupContext;
2588b783a1SJames Wright struct SetupContext_ {
2688b783a1SJames Wright   CeedScalar theta0;
2788b783a1SJames Wright   CeedScalar thetaC;
2888b783a1SJames Wright   CeedScalar P0;
2988b783a1SJames Wright   CeedScalar N;
3088b783a1SJames Wright   CeedScalar cv;
3188b783a1SJames Wright   CeedScalar cp;
3288b783a1SJames Wright   CeedScalar g;
3388b783a1SJames Wright   CeedScalar rc;
3488b783a1SJames Wright   CeedScalar lx;
3588b783a1SJames Wright   CeedScalar ly;
3688b783a1SJames Wright   CeedScalar lz;
3788b783a1SJames Wright   CeedScalar center[3];
3888b783a1SJames Wright   CeedScalar dc_axis[3];
3988b783a1SJames Wright   CeedScalar wind[3];
4088b783a1SJames Wright   CeedScalar time;
4188b783a1SJames Wright   int wind_type;              // See WindType: 0=ROTATION, 1=TRANSLATION
4288b783a1SJames Wright   int bubble_type;            // See BubbleType: 0=SPHERE, 1=CYLINDER
4388b783a1SJames Wright   int bubble_continuity_type; // See BubbleContinuityType: 0=SMOOTH, 1=BACK_SHARP 2=THICK
4488b783a1SJames Wright };
4588b783a1SJames Wright #endif
4688b783a1SJames Wright 
4788b783a1SJames Wright #ifndef newtonian_context_struct
4888b783a1SJames Wright #define newtonian_context_struct
4988b783a1SJames Wright typedef enum {
5088b783a1SJames Wright   STAB_NONE = 0,
5188b783a1SJames Wright   STAB_SU   = 1, // Streamline Upwind
5288b783a1SJames Wright   STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin
5388b783a1SJames Wright } StabilizationType;
5488b783a1SJames Wright 
5588b783a1SJames Wright typedef struct NewtonianIdealGasContext_ *NewtonianIdealGasContext;
5688b783a1SJames Wright struct NewtonianIdealGasContext_ {
5788b783a1SJames Wright   CeedScalar lambda;
5888b783a1SJames Wright   CeedScalar mu;
5988b783a1SJames Wright   CeedScalar k;
6088b783a1SJames Wright   CeedScalar cv;
6188b783a1SJames Wright   CeedScalar cp;
6288b783a1SJames Wright   CeedScalar g;
6388b783a1SJames Wright   CeedScalar c_tau;
6488b783a1SJames Wright   StabilizationType stabilization;
6588b783a1SJames Wright };
6688b783a1SJames Wright #endif
6788b783a1SJames Wright 
6888b783a1SJames Wright // *****************************************************************************
6988b783a1SJames Wright // Helper function for computing flux Jacobian
7088b783a1SJames Wright // *****************************************************************************
7188b783a1SJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NS(CeedScalar dF[3][5][5],
7288b783a1SJames Wright     const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
7388b783a1SJames Wright     const CeedScalar gamma, const CeedScalar g, CeedScalar z) {
7488b783a1SJames Wright   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square
7588b783a1SJames Wright   for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions
7688b783a1SJames Wright     for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix
7788b783a1SJames Wright       dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2. - g*z)) : 0.) - u[i]*u[j];
7888b783a1SJames Wright       for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix
7988b783a1SJames Wright         dF[i][0][k+1]   = ((i==k) ? 1. : 0.);
8088b783a1SJames Wright         dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) +
8188b783a1SJames Wright                           ((i==k) ? u[j] : 0.) -
8288b783a1SJames Wright                           ((i==j) ? u[k] : 0.) * (gamma-1.);
8388b783a1SJames Wright         dF[i][4][k+1]   = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) -
8488b783a1SJames Wright                           (gamma-1.)*u[i]*u[k];
8588b783a1SJames Wright       }
8688b783a1SJames Wright       dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.);
8788b783a1SJames Wright     }
8888b783a1SJames Wright     dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho);
8988b783a1SJames Wright     dF[i][4][4] = u[i] * gamma;
9088b783a1SJames Wright   }
9188b783a1SJames Wright }
9288b783a1SJames Wright 
9388b783a1SJames Wright // *****************************************************************************
9488b783a1SJames Wright // Helper function for computing Tau elements (stabilization constant)
9588b783a1SJames Wright //   Model from:
9688b783a1SJames Wright //     Stabilized Methods for Compressible Flows, Hughes et al 2010
9788b783a1SJames Wright //
9888b783a1SJames Wright //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
9988b783a1SJames Wright //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
10088b783a1SJames Wright //
10188b783a1SJames Wright // Where
10288b783a1SJames Wright //   c_tau     = stabilization constant (0.5 is reported as "optimal")
10388b783a1SJames Wright //   h[i]      = 2 length(dxdX[i])
10488b783a1SJames Wright //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
10588b783a1SJames Wright //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
10688b783a1SJames Wright //   rho(A[i]) = spectral radius of the convective flux Jacobian i,
10788b783a1SJames Wright //               wave speed in direction i
10888b783a1SJames Wright // *****************************************************************************
10988b783a1SJames Wright CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3],
11088b783a1SJames Wright                                        const CeedScalar dXdx[3][3], const CeedScalar u[3],
11188b783a1SJames Wright                                        const CeedScalar sound_speed, const CeedScalar c_tau) {
11288b783a1SJames Wright   for (int i=0; i<3; i++) {
11388b783a1SJames Wright     // length of element in direction i
11488b783a1SJames Wright     CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] +
11588b783a1SJames Wright                             dXdx[2][i]*dXdx[2][i]);
11688b783a1SJames Wright     // fastest wave in direction i
11788b783a1SJames Wright     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
11888b783a1SJames Wright     Tau_x[i] = c_tau * h / fastest_wave;
11988b783a1SJames Wright   }
12088b783a1SJames Wright }
12188b783a1SJames Wright 
12288b783a1SJames Wright // *****************************************************************************
12388b783a1SJames Wright // This QFunction sets a "still" initial condition for generic Newtonian IG problems
12488b783a1SJames Wright // *****************************************************************************
12588b783a1SJames Wright CEED_QFUNCTION(ICsNewtonianIG)(void *ctx, CeedInt Q,
12688b783a1SJames Wright                                const CeedScalar *const *in, CeedScalar *const *out) {
12788b783a1SJames Wright   // Inputs
12888b783a1SJames Wright   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
12988b783a1SJames Wright 
13088b783a1SJames Wright   // Outputs
13188b783a1SJames Wright   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
13288b783a1SJames Wright 
13388b783a1SJames Wright   // Quadrature Point Loop
13488b783a1SJames Wright   CeedPragmaSIMD
13588b783a1SJames Wright   for (CeedInt i=0; i<Q; i++) {
13688b783a1SJames Wright     CeedScalar q[5] = {0.};
13788b783a1SJames Wright 
13888b783a1SJames Wright     // Context
13988b783a1SJames Wright     const SetupContext context = (SetupContext)ctx;
14088b783a1SJames Wright     const CeedScalar theta0    = context->theta0;
14188b783a1SJames Wright     const CeedScalar P0        = context->P0;
14288b783a1SJames Wright     const CeedScalar N         = context->N;
14388b783a1SJames Wright     const CeedScalar cv        = context->cv;
14488b783a1SJames Wright     const CeedScalar cp        = context->cp;
14588b783a1SJames Wright     const CeedScalar g         = context->g;
14688b783a1SJames Wright     const CeedScalar Rd        = cp - cv;
14788b783a1SJames Wright 
14888b783a1SJames Wright     // Setup
14988b783a1SJames Wright     // -- Coordinates
15088b783a1SJames Wright     const CeedScalar z = X[2][i];
15188b783a1SJames Wright 
15288b783a1SJames Wright     // -- Exner pressure, hydrostatic balance
15388b783a1SJames Wright     const CeedScalar Pi = 1. + g*g*(exp(-N*N*z/g) - 1.) / (cp*theta0*N*N);
15488b783a1SJames Wright 
15588b783a1SJames Wright     // -- Density
15688b783a1SJames Wright     const CeedScalar rho = P0 * pow(Pi, cv/Rd) / (Rd*theta0);
15788b783a1SJames Wright 
15888b783a1SJames Wright     // Initial Conditions
15988b783a1SJames Wright     q[0] = rho;
16088b783a1SJames Wright     q[1] = 0.0;
16188b783a1SJames Wright     q[2] = 0.0;
16288b783a1SJames Wright     q[3] = 0.0;
16388b783a1SJames Wright     q[4] = rho * (cv*theta0*Pi + g*z);
16488b783a1SJames Wright 
16588b783a1SJames Wright     for (CeedInt j=0; j<5; j++)
16688b783a1SJames Wright       q0[j][i] = q[j];
16788b783a1SJames Wright   } // End of Quadrature Point Loop
16888b783a1SJames Wright   return 0;
16988b783a1SJames Wright }
17088b783a1SJames Wright 
17188b783a1SJames Wright // *****************************************************************************
17288b783a1SJames Wright // This QFunction implements the following formulation of Navier-Stokes with
17388b783a1SJames Wright //   explicit time stepping method
17488b783a1SJames Wright //
17588b783a1SJames Wright // This is 3D compressible Navier-Stokes in conservation form with state
17688b783a1SJames Wright //   variables of density, momentum density, and total energy density.
17788b783a1SJames Wright //
17888b783a1SJames Wright // State Variables: q = ( rho, U1, U2, U3, E )
17988b783a1SJames Wright //   rho - Mass Density
18088b783a1SJames Wright //   Ui  - Momentum Density,      Ui = rho ui
18188b783a1SJames Wright //   E   - Total Energy Density,  E  = rho (cv T + (u u)/2 + g z)
18288b783a1SJames Wright //
18388b783a1SJames Wright // Navier-Stokes Equations:
18488b783a1SJames Wright //   drho/dt + div( U )                               = 0
18588b783a1SJames Wright //   dU/dt   + div( rho (u x u) + P I3 ) + rho g khat = div( Fu )
18688b783a1SJames Wright //   dE/dt   + div( (E + P) u )                       = div( Fe )
18788b783a1SJames Wright //
18888b783a1SJames Wright // Viscous Stress:
18988b783a1SJames Wright //   Fu = mu (grad( u ) + grad( u )^T + lambda div ( u ) I3)
19088b783a1SJames Wright //
19188b783a1SJames Wright // Thermal Stress:
19288b783a1SJames Wright //   Fe = u Fu + k grad( T )
19388b783a1SJames Wright //
19488b783a1SJames Wright // Equation of State:
19588b783a1SJames Wright //   P = (gamma - 1) (E - rho (u u) / 2 - rho g z)
19688b783a1SJames Wright //
19788b783a1SJames Wright // Stabilization:
19888b783a1SJames Wright //   Tau = diag(TauC, TauM, TauM, TauM, TauE)
19988b783a1SJames Wright //     f1 = rho  sqrt(ui uj gij)
20088b783a1SJames Wright //     gij = dXi/dX * dXi/dX
20188b783a1SJames Wright //     TauC = Cc f1 / (8 gii)
20288b783a1SJames Wright //     TauM = min( 1 , 1 / f1 )
20388b783a1SJames Wright //     TauE = TauM / (Ce cv)
20488b783a1SJames Wright //
20588b783a1SJames Wright //  SU   = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) )
20688b783a1SJames Wright //
20788b783a1SJames Wright // Constants:
20888b783a1SJames Wright //   lambda = - 2 / 3,  From Stokes hypothesis
20988b783a1SJames Wright //   mu              ,  Dynamic viscosity
21088b783a1SJames Wright //   k               ,  Thermal conductivity
21188b783a1SJames Wright //   cv              ,  Specific heat, constant volume
21288b783a1SJames Wright //   cp              ,  Specific heat, constant pressure
21388b783a1SJames Wright //   g               ,  Gravity
21488b783a1SJames Wright //   gamma  = cp / cv,  Specific heat ratio
21588b783a1SJames Wright //
21688b783a1SJames Wright // We require the product of the inverse of the Jacobian (dXdx_j,k) and
21788b783a1SJames Wright // its transpose (dXdx_k,j) to properly compute integrals of the form:
21888b783a1SJames Wright // int( gradv gradu )
21988b783a1SJames Wright //
22088b783a1SJames Wright // *****************************************************************************
22188b783a1SJames Wright CEED_QFUNCTION(Newtonian)(void *ctx, CeedInt Q,
22288b783a1SJames Wright                           const CeedScalar *const *in, CeedScalar *const *out) {
22388b783a1SJames Wright   // *INDENT-OFF*
22488b783a1SJames Wright   // Inputs
22588b783a1SJames Wright   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
22688b783a1SJames Wright                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
22788b783a1SJames Wright                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
22888b783a1SJames Wright                    (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
22988b783a1SJames Wright   // Outputs
23088b783a1SJames Wright   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
23188b783a1SJames Wright              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
23288b783a1SJames Wright   // *INDENT-ON*
23388b783a1SJames Wright 
23488b783a1SJames Wright   // Context
23588b783a1SJames Wright   NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx;
23688b783a1SJames Wright   const CeedScalar lambda = context->lambda;
23788b783a1SJames Wright   const CeedScalar mu     = context->mu;
23888b783a1SJames Wright   const CeedScalar k      = context->k;
23988b783a1SJames Wright   const CeedScalar cv     = context->cv;
24088b783a1SJames Wright   const CeedScalar cp     = context->cp;
24188b783a1SJames Wright   const CeedScalar g      = context->g;
24288b783a1SJames Wright   const CeedScalar c_tau  = context->c_tau;
24388b783a1SJames Wright   const CeedScalar gamma  = cp / cv;
24488b783a1SJames Wright 
24588b783a1SJames Wright   CeedPragmaSIMD
24688b783a1SJames Wright   // Quadrature Point Loop
24788b783a1SJames Wright   for (CeedInt i=0; i<Q; i++) {
24888b783a1SJames Wright     // *INDENT-OFF*
24988b783a1SJames Wright     // Setup
25088b783a1SJames Wright     // -- Interp in
25188b783a1SJames Wright     const CeedScalar rho        =   q[0][i];
25288b783a1SJames Wright     const CeedScalar u[3]       =  {q[1][i] / rho,
25388b783a1SJames Wright                                     q[2][i] / rho,
25488b783a1SJames Wright                                     q[3][i] / rho
25588b783a1SJames Wright                                    };
25688b783a1SJames Wright     const CeedScalar E          =   q[4][i];
25788b783a1SJames Wright     // -- Grad in
25888b783a1SJames Wright     const CeedScalar drho[3]    =  {dq[0][0][i],
25988b783a1SJames Wright                                     dq[1][0][i],
26088b783a1SJames Wright                                     dq[2][0][i]
26188b783a1SJames Wright                                    };
26288b783a1SJames Wright     const CeedScalar dU[3][3]   = {{dq[0][1][i],
26388b783a1SJames Wright                                     dq[1][1][i],
26488b783a1SJames Wright                                     dq[2][1][i]},
26588b783a1SJames Wright                                    {dq[0][2][i],
26688b783a1SJames Wright                                     dq[1][2][i],
26788b783a1SJames Wright                                     dq[2][2][i]},
26888b783a1SJames Wright                                    {dq[0][3][i],
26988b783a1SJames Wright                                     dq[1][3][i],
27088b783a1SJames Wright                                     dq[2][3][i]}
27188b783a1SJames Wright                                   };
27288b783a1SJames Wright     const CeedScalar dE[3]      =  {dq[0][4][i],
27388b783a1SJames Wright                                     dq[1][4][i],
27488b783a1SJames Wright                                     dq[2][4][i]
27588b783a1SJames Wright                                    };
27688b783a1SJames Wright     // -- Interp-to-Interp q_data
27788b783a1SJames Wright     const CeedScalar wdetJ      =   q_data[0][i];
27888b783a1SJames Wright     // -- Interp-to-Grad q_data
27988b783a1SJames Wright     // ---- Inverse of change of coordinate matrix: X_i,j
28088b783a1SJames Wright     // *INDENT-OFF*
28188b783a1SJames Wright     const CeedScalar dXdx[3][3] = {{q_data[1][i],
28288b783a1SJames Wright                                     q_data[2][i],
28388b783a1SJames Wright                                     q_data[3][i]},
28488b783a1SJames Wright                                    {q_data[4][i],
28588b783a1SJames Wright                                     q_data[5][i],
28688b783a1SJames Wright                                     q_data[6][i]},
28788b783a1SJames Wright                                    {q_data[7][i],
28888b783a1SJames Wright                                     q_data[8][i],
28988b783a1SJames Wright                                     q_data[9][i]}
29088b783a1SJames Wright                                   };
29188b783a1SJames Wright     // *INDENT-ON*
29288b783a1SJames Wright     // -- Grad-to-Grad q_data
29388b783a1SJames Wright     // dU/dx
29488b783a1SJames Wright     CeedScalar du[3][3] = {{0}};
29588b783a1SJames Wright     CeedScalar drhodx[3] = {0};
29688b783a1SJames Wright     CeedScalar dEdx[3] = {0};
29788b783a1SJames Wright     CeedScalar dUdx[3][3] = {{0}};
29888b783a1SJames Wright     CeedScalar dXdxdXdxT[3][3] = {{0}};
29988b783a1SJames Wright     for (int j=0; j<3; j++) {
30088b783a1SJames Wright       for (int k=0; k<3; k++) {
30188b783a1SJames Wright         du[j][k] = (dU[j][k] - drho[k]*u[j]) / rho;
30288b783a1SJames Wright         drhodx[j] += drho[k] * dXdx[k][j];
30388b783a1SJames Wright         dEdx[j] += dE[k] * dXdx[k][j];
30488b783a1SJames Wright         for (int l=0; l<3; l++) {
30588b783a1SJames Wright           dUdx[j][k] += dU[j][l] * dXdx[l][k];
30688b783a1SJames Wright           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
30788b783a1SJames Wright         }
30888b783a1SJames Wright       }
30988b783a1SJames Wright     }
31088b783a1SJames Wright     CeedScalar dudx[3][3] = {{0}};
31188b783a1SJames Wright     for (int j=0; j<3; j++)
31288b783a1SJames Wright       for (int k=0; k<3; k++)
31388b783a1SJames Wright         for (int l=0; l<3; l++)
31488b783a1SJames Wright           dudx[j][k] += du[j][l] * dXdx[l][k];
31588b783a1SJames Wright     // -- grad_T
31688b783a1SJames Wright     const CeedScalar grad_T[3]  = {(dEdx[0]/rho - E*drhodx[0]/(rho*rho) - /* *NOPAD* */
31788b783a1SJames Wright                                     (u[0]*dudx[0][0] + u[1]*dudx[1][0] + u[2]*dudx[2][0]))/cv,
31888b783a1SJames Wright                                    (dEdx[1]/rho - E*drhodx[1]/(rho*rho) - /* *NOPAD* */
31988b783a1SJames Wright                                     (u[0]*dudx[0][1] + u[1]*dudx[1][1] + u[2]*dudx[2][1]))/cv,
32088b783a1SJames Wright                                    (dEdx[2]/rho - E*drhodx[2]/(rho*rho) - /* *NOPAD* */
32188b783a1SJames Wright                                     (u[0]*dudx[0][2] + u[1]*dudx[1][2] + u[2]*dudx[2][2]) - g)/cv
32288b783a1SJames Wright                                   };
32388b783a1SJames Wright 
32488b783a1SJames Wright     // -- Fuvisc
32588b783a1SJames Wright     // ---- Symmetric 3x3 matrix
32688b783a1SJames Wright     const CeedScalar Fu[6]     =  {mu*(dudx[0][0] * (2 + lambda) + /* *NOPAD* */
32788b783a1SJames Wright                                        lambda * (dudx[1][1] + dudx[2][2])),
32888b783a1SJames Wright                                    mu*(dudx[0][1] + dudx[1][0]), /* *NOPAD* */
32988b783a1SJames Wright                                    mu*(dudx[0][2] + dudx[2][0]), /* *NOPAD* */
33088b783a1SJames Wright                                    mu*(dudx[1][1] * (2 + lambda) + /* *NOPAD* */
33188b783a1SJames Wright                                        lambda * (dudx[0][0] + dudx[2][2])),
33288b783a1SJames Wright                                    mu*(dudx[1][2] + dudx[2][1]), /* *NOPAD* */
33388b783a1SJames Wright                                    mu*(dudx[2][2] * (2 + lambda) + /* *NOPAD* */
33488b783a1SJames Wright                                        lambda * (dudx[0][0] + dudx[1][1]))
33588b783a1SJames Wright                                   };
33688b783a1SJames Wright     // -- Fevisc
33788b783a1SJames Wright     const CeedScalar Fe[3]     =  {u[0]*Fu[0] + u[1]*Fu[1] + u[2]*Fu[2] + /* *NOPAD* */
33888b783a1SJames Wright                                    k*grad_T[0], /* *NOPAD* */
33988b783a1SJames Wright                                    u[0]*Fu[1] + u[1]*Fu[3] + u[2]*Fu[4] + /* *NOPAD* */
34088b783a1SJames Wright                                    k*grad_T[1], /* *NOPAD* */
34188b783a1SJames Wright                                    u[0]*Fu[2] + u[1]*Fu[4] + u[2]*Fu[5] + /* *NOPAD* */
34288b783a1SJames Wright                                    k*grad_T[2] /* *NOPAD* */
34388b783a1SJames Wright                                   };
34488b783a1SJames Wright     // Pressure
34588b783a1SJames Wright     const CeedScalar
34688b783a1SJames Wright     E_kinetic   = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
34788b783a1SJames Wright     E_potential = rho*g*x[2][i],
34888b783a1SJames Wright     E_internal  = E - E_kinetic - E_potential,
34988b783a1SJames Wright     P           = E_internal * (gamma - 1.); // P = pressure
35088b783a1SJames Wright 
35188b783a1SJames Wright     // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
35288b783a1SJames Wright     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
35388b783a1SJames Wright     computeFluxJacobian_NS(jacob_F_conv, rho, u, E, gamma, g, x[2][i]);
35488b783a1SJames Wright 
35588b783a1SJames Wright     // jacob_F_conv_T = jacob_F_conv^T
35688b783a1SJames Wright     CeedScalar jacob_F_conv_T[3][5][5];
35788b783a1SJames Wright     for (int j=0; j<3; j++)
35888b783a1SJames Wright       for (int k=0; k<5; k++)
35988b783a1SJames Wright         for (int l=0; l<5; l++)
36088b783a1SJames Wright           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
36188b783a1SJames Wright 
36288b783a1SJames Wright     // dqdx collects drhodx, dUdx and dEdx in one vector
36388b783a1SJames Wright     CeedScalar dqdx[5][3];
36488b783a1SJames Wright     for (int j=0; j<3; j++) {
36588b783a1SJames Wright       dqdx[0][j] = drhodx[j];
36688b783a1SJames Wright       dqdx[4][j] = dEdx[j];
36788b783a1SJames Wright       for (int k=0; k<3; k++)
36888b783a1SJames Wright         dqdx[k+1][j] = dUdx[k][j];
36988b783a1SJames Wright     }
37088b783a1SJames Wright 
37188b783a1SJames Wright     // strong_conv = dF/dq * dq/dx    (Strong convection)
37288b783a1SJames Wright     CeedScalar strong_conv[5] = {0};
37388b783a1SJames Wright     for (int j=0; j<3; j++)
37488b783a1SJames Wright       for (int k=0; k<5; k++)
37588b783a1SJames Wright         for (int l=0; l<5; l++)
37688b783a1SJames Wright           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
37788b783a1SJames Wright 
37888b783a1SJames Wright     // Body force
37988b783a1SJames Wright     const CeedScalar body_force[5] = {0, 0, 0, -rho*g, 0};
38088b783a1SJames Wright 
38188b783a1SJames Wright     // The Physics
38288b783a1SJames Wright     // Zero dv so all future terms can safely sum into it
38388b783a1SJames Wright     for (int j=0; j<5; j++)
38488b783a1SJames Wright       for (int k=0; k<3; k++)
38588b783a1SJames Wright         dv[k][j][i] = 0;
38688b783a1SJames Wright 
38788b783a1SJames Wright     // -- Density
38888b783a1SJames Wright     // ---- u rho
38988b783a1SJames Wright     for (int j=0; j<3; j++)
39088b783a1SJames Wright       dv[j][0][i]  += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
39188b783a1SJames Wright                              rho*u[2]*dXdx[j][2]);
39288b783a1SJames Wright     // -- Momentum
39388b783a1SJames Wright     // ---- rho (u x u) + P I3
39488b783a1SJames Wright     for (int j=0; j<3; j++)
39588b783a1SJames Wright       for (int k=0; k<3; k++)
39688b783a1SJames Wright         dv[k][j+1][i]  += wdetJ*((rho*u[j]*u[0] + (j==0?P:0))*dXdx[k][0] +
39788b783a1SJames Wright                                  (rho*u[j]*u[1] + (j==1?P:0))*dXdx[k][1] +
39888b783a1SJames Wright                                  (rho*u[j]*u[2] + (j==2?P:0))*dXdx[k][2]);
39988b783a1SJames Wright     // ---- Fuvisc
40088b783a1SJames Wright     const CeedInt Fuviscidx[3][3] = {{0, 1, 2}, {1, 3, 4}, {2, 4, 5}}; // symmetric matrix indices
40188b783a1SJames Wright     for (int j=0; j<3; j++)
40288b783a1SJames Wright       for (int k=0; k<3; k++)
40388b783a1SJames Wright         dv[k][j+1][i] -= wdetJ*(Fu[Fuviscidx[j][0]]*dXdx[k][0] +
40488b783a1SJames Wright                                 Fu[Fuviscidx[j][1]]*dXdx[k][1] +
40588b783a1SJames Wright                                 Fu[Fuviscidx[j][2]]*dXdx[k][2]);
40688b783a1SJames Wright     // -- Total Energy Density
40788b783a1SJames Wright     // ---- (E + P) u
40888b783a1SJames Wright     for (int j=0; j<3; j++)
40988b783a1SJames Wright       dv[j][4][i]  += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
41088b783a1SJames Wright                                          u[2]*dXdx[j][2]);
41188b783a1SJames Wright     // ---- Fevisc
41288b783a1SJames Wright     for (int j=0; j<3; j++)
41388b783a1SJames Wright       dv[j][4][i] -= wdetJ * (Fe[0]*dXdx[j][0] + Fe[1]*dXdx[j][1] +
41488b783a1SJames Wright                               Fe[2]*dXdx[j][2]);
41588b783a1SJames Wright     // Body Force
41688b783a1SJames Wright     for (int j=0; j<5; j++)
41788b783a1SJames Wright       v[j][i] = wdetJ * body_force[j];
41888b783a1SJames Wright 
41988b783a1SJames Wright     // Stabilization
42088b783a1SJames Wright     // -- Tau elements
42188b783a1SJames Wright     const CeedScalar sound_speed = sqrt(gamma * P / rho);
42288b783a1SJames Wright     CeedScalar Tau_x[3] = {0.};
42388b783a1SJames Wright     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
42488b783a1SJames Wright 
42588b783a1SJames Wright     // -- Stabilization method: none or SU
42688b783a1SJames Wright     CeedScalar stab[5][3];
42788b783a1SJames Wright     switch (context->stabilization) {
42888b783a1SJames Wright     case STAB_NONE:        // Galerkin
42988b783a1SJames Wright       break;
43088b783a1SJames Wright     case STAB_SU:        // SU
43188b783a1SJames Wright       for (int j=0; j<3; j++)
43288b783a1SJames Wright         for (int k=0; k<5; k++)
43388b783a1SJames Wright           for (int l=0; l<5; l++)
43488b783a1SJames Wright             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
43588b783a1SJames Wright 
43688b783a1SJames Wright       for (int j=0; j<5; j++)
43788b783a1SJames Wright         for (int k=0; k<3; k++)
43888b783a1SJames Wright           dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] +
43988b783a1SJames Wright                                 stab[j][1] * dXdx[k][1] +
44088b783a1SJames Wright                                 stab[j][2] * dXdx[k][2]);
44188b783a1SJames Wright       break;
44288b783a1SJames Wright     case STAB_SUPG:        // SUPG is not implemented for explicit scheme
44388b783a1SJames Wright       break;
44488b783a1SJames Wright     }
44588b783a1SJames Wright 
44688b783a1SJames Wright   } // End Quadrature Point Loop
44788b783a1SJames Wright 
44888b783a1SJames Wright   // Return
44988b783a1SJames Wright   return 0;
45088b783a1SJames Wright }
45188b783a1SJames Wright 
45288b783a1SJames Wright // *****************************************************************************
45388b783a1SJames Wright // This QFunction implements the Navier-Stokes equations (mentioned above) with
45488b783a1SJames Wright //   implicit time stepping method
45588b783a1SJames Wright //
45688b783a1SJames Wright //  SU   = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) )
45788b783a1SJames Wright //  SUPG = Galerkin + grad(v) . ( Ai^T * Tau * (q_dot + Aj q,j - body force) )
45888b783a1SJames Wright //                                       (diffussive terms will be added later)
45988b783a1SJames Wright //
46088b783a1SJames Wright // *****************************************************************************
46188b783a1SJames Wright CEED_QFUNCTION(IFunction_Newtonian)(void *ctx, CeedInt Q,
46288b783a1SJames Wright                                     const CeedScalar *const *in,
46388b783a1SJames Wright                                     CeedScalar *const *out) {
46488b783a1SJames Wright   // *INDENT-OFF*
46588b783a1SJames Wright   // Inputs
46688b783a1SJames Wright   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
46788b783a1SJames Wright                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
46888b783a1SJames Wright                    (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
46988b783a1SJames Wright                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3],
47088b783a1SJames Wright                    (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4];
47188b783a1SJames Wright   // Outputs
47288b783a1SJames Wright   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
47388b783a1SJames Wright              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
47488b783a1SJames Wright   // *INDENT-ON*
47588b783a1SJames Wright   // Context
47688b783a1SJames Wright   NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx;
47788b783a1SJames Wright   const CeedScalar lambda = context->lambda;
47888b783a1SJames Wright   const CeedScalar mu     = context->mu;
47988b783a1SJames Wright   const CeedScalar k      = context->k;
48088b783a1SJames Wright   const CeedScalar cv     = context->cv;
48188b783a1SJames Wright   const CeedScalar cp     = context->cp;
48288b783a1SJames Wright   const CeedScalar g      = context->g;
48388b783a1SJames Wright   const CeedScalar c_tau  = context->c_tau;
48488b783a1SJames Wright   const CeedScalar gamma  = cp / cv;
48588b783a1SJames Wright 
48688b783a1SJames Wright   CeedPragmaSIMD
48788b783a1SJames Wright   // Quadrature Point Loop
48888b783a1SJames Wright   for (CeedInt i=0; i<Q; i++) {
48988b783a1SJames Wright     // Setup
49088b783a1SJames Wright     // -- Interp in
49188b783a1SJames Wright     const CeedScalar rho        =   q[0][i];
49288b783a1SJames Wright     const CeedScalar u[3]       =  {q[1][i] / rho,
49388b783a1SJames Wright                                     q[2][i] / rho,
49488b783a1SJames Wright                                     q[3][i] / rho
49588b783a1SJames Wright                                    };
49688b783a1SJames Wright     const CeedScalar E          =   q[4][i];
49788b783a1SJames Wright     // -- Grad in
49888b783a1SJames Wright     const CeedScalar drho[3]    =  {dq[0][0][i],
49988b783a1SJames Wright                                     dq[1][0][i],
50088b783a1SJames Wright                                     dq[2][0][i]
50188b783a1SJames Wright                                    };
50288b783a1SJames Wright     // *INDENT-OFF*
50388b783a1SJames Wright     const CeedScalar dU[3][3]   = {{dq[0][1][i],
50488b783a1SJames Wright                                     dq[1][1][i],
50588b783a1SJames Wright                                     dq[2][1][i]},
50688b783a1SJames Wright                                    {dq[0][2][i],
50788b783a1SJames Wright                                     dq[1][2][i],
50888b783a1SJames Wright                                     dq[2][2][i]},
50988b783a1SJames Wright                                    {dq[0][3][i],
51088b783a1SJames Wright                                     dq[1][3][i],
51188b783a1SJames Wright                                     dq[2][3][i]}
51288b783a1SJames Wright                                   };
51388b783a1SJames Wright     // *INDENT-ON*
51488b783a1SJames Wright     const CeedScalar dE[3]      =  {dq[0][4][i],
51588b783a1SJames Wright                                     dq[1][4][i],
51688b783a1SJames Wright                                     dq[2][4][i]
51788b783a1SJames Wright                                    };
51888b783a1SJames Wright     // -- Interp-to-Interp q_data
51988b783a1SJames Wright     const CeedScalar wdetJ      =   q_data[0][i];
52088b783a1SJames Wright     // -- Interp-to-Grad q_data
52188b783a1SJames Wright     // ---- Inverse of change of coordinate matrix: X_i,j
52288b783a1SJames Wright     // *INDENT-OFF*
52388b783a1SJames Wright     const CeedScalar dXdx[3][3] = {{q_data[1][i],
52488b783a1SJames Wright                                     q_data[2][i],
52588b783a1SJames Wright                                     q_data[3][i]},
52688b783a1SJames Wright                                    {q_data[4][i],
52788b783a1SJames Wright                                     q_data[5][i],
52888b783a1SJames Wright                                     q_data[6][i]},
52988b783a1SJames Wright                                    {q_data[7][i],
53088b783a1SJames Wright                                     q_data[8][i],
53188b783a1SJames Wright                                     q_data[9][i]}
53288b783a1SJames Wright                                   };
53388b783a1SJames Wright     // *INDENT-ON*
53488b783a1SJames Wright     // -- Grad-to-Grad q_data
53588b783a1SJames Wright     // dU/dx
53688b783a1SJames Wright     CeedScalar du[3][3] = {{0}};
53788b783a1SJames Wright     CeedScalar drhodx[3] = {0};
53888b783a1SJames Wright     CeedScalar dEdx[3] = {0};
53988b783a1SJames Wright     CeedScalar dUdx[3][3] = {{0}};
54088b783a1SJames Wright     CeedScalar dXdxdXdxT[3][3] = {{0}};
54188b783a1SJames Wright     for (int j=0; j<3; j++) {
54288b783a1SJames Wright       for (int k=0; k<3; k++) {
54388b783a1SJames Wright         du[j][k] = (dU[j][k] - drho[k]*u[j]) / rho;
54488b783a1SJames Wright         drhodx[j] += drho[k] * dXdx[k][j];
54588b783a1SJames Wright         dEdx[j] += dE[k] * dXdx[k][j];
54688b783a1SJames Wright         for (int l=0; l<3; l++) {
54788b783a1SJames Wright           dUdx[j][k] += dU[j][l] * dXdx[l][k];
54888b783a1SJames Wright           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
54988b783a1SJames Wright         }
55088b783a1SJames Wright       }
55188b783a1SJames Wright     }
55288b783a1SJames Wright     CeedScalar dudx[3][3] = {{0}};
55388b783a1SJames Wright     for (int j=0; j<3; j++)
55488b783a1SJames Wright       for (int k=0; k<3; k++)
55588b783a1SJames Wright         for (int l=0; l<3; l++)
55688b783a1SJames Wright           dudx[j][k] += du[j][l] * dXdx[l][k];
55788b783a1SJames Wright     // -- grad_T
55888b783a1SJames Wright     const CeedScalar grad_T[3]  = {(dEdx[0]/rho - E*drhodx[0]/(rho*rho) - /* *NOPAD* */
55988b783a1SJames Wright                                     (u[0]*dudx[0][0] + u[1]*dudx[1][0] + u[2]*dudx[2][0]))/cv,
56088b783a1SJames Wright                                    (dEdx[1]/rho - E*drhodx[1]/(rho*rho) - /* *NOPAD* */
56188b783a1SJames Wright                                     (u[0]*dudx[0][1] + u[1]*dudx[1][1] + u[2]*dudx[2][1]))/cv,
56288b783a1SJames Wright                                    (dEdx[2]/rho - E*drhodx[2]/(rho*rho) - /* *NOPAD* */
56388b783a1SJames Wright                                     (u[0]*dudx[0][2] + u[1]*dudx[1][2] + u[2]*dudx[2][2]) - g)/cv
56488b783a1SJames Wright                                   };
56588b783a1SJames Wright     // -- Fuvisc
56688b783a1SJames Wright     // ---- Symmetric 3x3 matrix
56788b783a1SJames Wright     const CeedScalar Fu[6]     =  {mu*(dudx[0][0] * (2 + lambda) + /* *NOPAD* */
56888b783a1SJames Wright                                        lambda * (dudx[1][1] + dudx[2][2])),
56988b783a1SJames Wright                                    mu*(dudx[0][1] + dudx[1][0]), /* *NOPAD* */
57088b783a1SJames Wright                                    mu*(dudx[0][2] + dudx[2][0]), /* *NOPAD* */
57188b783a1SJames Wright                                    mu*(dudx[1][1] * (2 + lambda) + /* *NOPAD* */
57288b783a1SJames Wright                                        lambda * (dudx[0][0] + dudx[2][2])),
57388b783a1SJames Wright                                    mu*(dudx[1][2] + dudx[2][1]), /* *NOPAD* */
57488b783a1SJames Wright                                    mu*(dudx[2][2] * (2 + lambda) + /* *NOPAD* */
57588b783a1SJames Wright                                        lambda * (dudx[0][0] + dudx[1][1]))
57688b783a1SJames Wright                                   };
57788b783a1SJames Wright     // -- Fevisc
57888b783a1SJames Wright     const CeedScalar Fe[3]     =  {u[0]*Fu[0] + u[1]*Fu[1] + u[2]*Fu[2] + /* *NOPAD* */
57988b783a1SJames Wright                                    k*grad_T[0], /* *NOPAD* */
58088b783a1SJames Wright                                    u[0]*Fu[1] + u[1]*Fu[3] + u[2]*Fu[4] + /* *NOPAD* */
58188b783a1SJames Wright                                    k*grad_T[1], /* *NOPAD* */
58288b783a1SJames Wright                                    u[0]*Fu[2] + u[1]*Fu[4] + u[2]*Fu[5] + /* *NOPAD* */
58388b783a1SJames Wright                                    k*grad_T[2] /* *NOPAD* */
58488b783a1SJames Wright                                   };
58588b783a1SJames Wright     // Pressure
58688b783a1SJames Wright     const CeedScalar
58788b783a1SJames Wright     E_kinetic   = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
58888b783a1SJames Wright     E_potential = rho*g*x[2][i],
58988b783a1SJames Wright     E_internal  = E - E_kinetic - E_potential,
59088b783a1SJames Wright     P           = E_internal * (gamma - 1.); // P = pressure
59188b783a1SJames Wright 
59288b783a1SJames Wright     // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
59388b783a1SJames Wright     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
59488b783a1SJames Wright     computeFluxJacobian_NS(jacob_F_conv, rho, u, E, gamma, g, x[2][i]);
59588b783a1SJames Wright 
59688b783a1SJames Wright     // jacob_F_conv_T = jacob_F_conv^T
59788b783a1SJames Wright     CeedScalar jacob_F_conv_T[3][5][5];
59888b783a1SJames Wright     for (int j=0; j<3; j++)
59988b783a1SJames Wright       for (int k=0; k<5; k++)
60088b783a1SJames Wright         for (int l=0; l<5; l++)
60188b783a1SJames Wright           jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k];
60288b783a1SJames Wright     // dqdx collects drhodx, dUdx and dEdx in one vector
60388b783a1SJames Wright     CeedScalar dqdx[5][3];
60488b783a1SJames Wright     for (int j=0; j<3; j++) {
60588b783a1SJames Wright       dqdx[0][j] = drhodx[j];
60688b783a1SJames Wright       dqdx[4][j] = dEdx[j];
60788b783a1SJames Wright       for (int k=0; k<3; k++)
60888b783a1SJames Wright         dqdx[k+1][j] = dUdx[k][j];
60988b783a1SJames Wright     }
61088b783a1SJames Wright     // strong_conv = dF/dq * dq/dx    (Strong convection)
61188b783a1SJames Wright     CeedScalar strong_conv[5] = {0};
61288b783a1SJames Wright     for (int j=0; j<3; j++)
61388b783a1SJames Wright       for (int k=0; k<5; k++)
61488b783a1SJames Wright         for (int l=0; l<5; l++)
61588b783a1SJames Wright           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
61688b783a1SJames Wright 
61788b783a1SJames Wright     // Body force
61888b783a1SJames Wright     const CeedScalar body_force[5] = {0, 0, 0, -rho*g, 0};
61988b783a1SJames Wright 
62088b783a1SJames Wright     // Strong residual
62188b783a1SJames Wright     CeedScalar strong_res[5];
62288b783a1SJames Wright     for (int j=0; j<5; j++)
62388b783a1SJames Wright       strong_res[j] = q_dot[j][i] + strong_conv[j] - body_force[j];
62488b783a1SJames Wright 
62588b783a1SJames Wright     // The Physics
62688b783a1SJames Wright     //-----mass matrix
62788b783a1SJames Wright     for (int j=0; j<5; j++)
62888b783a1SJames Wright       v[j][i] = wdetJ*q_dot[j][i];
62988b783a1SJames Wright 
63088b783a1SJames Wright     // Zero dv so all future terms can safely sum into it
63188b783a1SJames Wright     for (int j=0; j<5; j++)
63288b783a1SJames Wright       for (int k=0; k<3; k++)
63388b783a1SJames Wright         dv[k][j][i] = 0;
63488b783a1SJames Wright 
63588b783a1SJames Wright     // -- Density
63688b783a1SJames Wright     // ---- u rho
63788b783a1SJames Wright     for (int j=0; j<3; j++)
63888b783a1SJames Wright       dv[j][0][i]  -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
63988b783a1SJames Wright                              rho*u[2]*dXdx[j][2]);
64088b783a1SJames Wright     // -- Momentum
64188b783a1SJames Wright     // ---- rho (u x u) + P I3
64288b783a1SJames Wright     for (int j=0; j<3; j++)
64388b783a1SJames Wright       for (int k=0; k<3; k++)
64488b783a1SJames Wright         dv[k][j+1][i]  -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0))*dXdx[k][0] +
64588b783a1SJames Wright                                  (rho*u[j]*u[1] + (j==1?P:0))*dXdx[k][1] +
64688b783a1SJames Wright                                  (rho*u[j]*u[2] + (j==2?P:0))*dXdx[k][2]);
64788b783a1SJames Wright     // ---- Fuvisc
64888b783a1SJames Wright     const CeedInt Fuviscidx[3][3] = {{0, 1, 2}, {1, 3, 4}, {2, 4, 5}}; // symmetric matrix indices
64988b783a1SJames Wright     for (int j=0; j<3; j++)
65088b783a1SJames Wright       for (int k=0; k<3; k++)
65188b783a1SJames Wright         dv[k][j+1][i] += wdetJ*(Fu[Fuviscidx[j][0]]*dXdx[k][0] +
65288b783a1SJames Wright                                 Fu[Fuviscidx[j][1]]*dXdx[k][1] +
65388b783a1SJames Wright                                 Fu[Fuviscidx[j][2]]*dXdx[k][2]);
65488b783a1SJames Wright     // -- Total Energy Density
65588b783a1SJames Wright     // ---- (E + P) u
65688b783a1SJames Wright     for (int j=0; j<3; j++)
65788b783a1SJames Wright       dv[j][4][i]  -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
65888b783a1SJames Wright                                          u[2]*dXdx[j][2]);
65988b783a1SJames Wright     // ---- Fevisc
66088b783a1SJames Wright     for (int j=0; j<3; j++)
66188b783a1SJames Wright       dv[j][4][i] += wdetJ * (Fe[0]*dXdx[j][0] + Fe[1]*dXdx[j][1] +
66288b783a1SJames Wright                               Fe[2]*dXdx[j][2]);
66388b783a1SJames Wright     // Body Force
66488b783a1SJames Wright     for (int j=0; j<5; j++)
66588b783a1SJames Wright       v[j][i] -= wdetJ*body_force[j];
66688b783a1SJames Wright 
66788b783a1SJames Wright     // Stabilization
66888b783a1SJames Wright     // -- Tau elements
66988b783a1SJames Wright     const CeedScalar sound_speed = sqrt(gamma * P / rho);
67088b783a1SJames Wright     CeedScalar Tau_x[3] = {0.};
67188b783a1SJames Wright     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
67288b783a1SJames Wright 
67388b783a1SJames Wright     // -- Stabilization method: none, SU, or SUPG
67488b783a1SJames Wright     CeedScalar stab[5][3];
67588b783a1SJames Wright     switch (context->stabilization) {
67688b783a1SJames Wright     case STAB_NONE:        // Galerkin
67788b783a1SJames Wright       break;
67888b783a1SJames Wright     case STAB_SU:        // SU
67988b783a1SJames Wright       for (int j=0; j<3; j++)
68088b783a1SJames Wright         for (int k=0; k<5; k++)
68188b783a1SJames Wright           for (int l=0; l<5; l++)
68288b783a1SJames Wright             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l];
68388b783a1SJames Wright 
68488b783a1SJames Wright       for (int j=0; j<5; j++)
68588b783a1SJames Wright         for (int k=0; k<3; k++)
68688b783a1SJames Wright           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
68788b783a1SJames Wright                                 stab[j][1] * dXdx[k][1] +
68888b783a1SJames Wright                                 stab[j][2] * dXdx[k][2]);
68988b783a1SJames Wright       break;
69088b783a1SJames Wright     case STAB_SUPG:        // SUPG
69188b783a1SJames Wright       for (int j=0; j<3; j++)
69288b783a1SJames Wright         for (int k=0; k<5; k++)
69388b783a1SJames Wright           for (int l=0; l<5; l++)
69488b783a1SJames Wright             stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_res[l];
69588b783a1SJames Wright 
69688b783a1SJames Wright       for (int j=0; j<5; j++)
69788b783a1SJames Wright         for (int k=0; k<3; k++)
69888b783a1SJames Wright           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
69988b783a1SJames Wright                                 stab[j][1] * dXdx[k][1] +
70088b783a1SJames Wright                                 stab[j][2] * dXdx[k][2]);
70188b783a1SJames Wright       break;
70288b783a1SJames Wright     }
70388b783a1SJames Wright 
70488b783a1SJames Wright   } // End Quadrature Point Loop
70588b783a1SJames Wright 
70688b783a1SJames Wright   // Return
70788b783a1SJames Wright   return 0;
70888b783a1SJames Wright }
70988b783a1SJames Wright // *****************************************************************************
71088b783a1SJames Wright #endif // newtonian_h
711