xref: /libCEED/examples/fluids/qfunctions/newtonian.h (revision 6f00d0e672188ac01d0f151e04000885fb0358d9)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
388b783a1SJames Wright //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
588b783a1SJames Wright //
63d8e8822SJeremy 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>
17841e4c73SJed Brown #include "newtonian_types.h"
1888b783a1SJames Wright 
1988b783a1SJames Wright #ifndef M_PI
2088b783a1SJames Wright #define M_PI    3.14159265358979323846
2188b783a1SJames Wright #endif
2288b783a1SJames Wright 
235c677226SJed Brown typedef struct {
245c677226SJed Brown   CeedScalar pressure;
255c677226SJed Brown   CeedScalar velocity[3];
265c677226SJed Brown   CeedScalar temperature;
275c677226SJed Brown } StatePrimitive;
285c677226SJed Brown 
295c677226SJed Brown typedef struct {
305c677226SJed Brown   CeedScalar density;
315c677226SJed Brown   CeedScalar momentum[3];
325c677226SJed Brown   CeedScalar E_total;
335c677226SJed Brown } StateConservative;
345c677226SJed Brown 
355c677226SJed Brown typedef struct {
365c677226SJed Brown   StateConservative U;
375c677226SJed Brown   StatePrimitive Y;
385c677226SJed Brown } State;
395c677226SJed Brown 
405c677226SJed Brown CEED_QFUNCTION_HELPER CeedScalar Dot3(const CeedScalar u[3],
415c677226SJed Brown                                       const CeedScalar v[3]) {
425c677226SJed Brown   return u[0]*v[0] + u[1]*v[1] + u[2]*v[2];
435c677226SJed Brown }
445c677226SJed Brown 
455c677226SJed Brown CEED_QFUNCTION_HELPER StatePrimitive StatePrimitiveFromConservative(
465c677226SJed Brown   NewtonianIdealGasContext gas, StateConservative U, const CeedScalar x[3]) {
475c677226SJed Brown   StatePrimitive Y;
485c677226SJed Brown   for (int i=0; i<3; i++) Y.velocity[i] = U.momentum[i] / U.density;
495c677226SJed Brown   CeedScalar e_kinetic = .5 * Dot3(Y.velocity, Y.velocity);
505c677226SJed Brown   CeedScalar e_potential = -Dot3(gas->g, x);
515c677226SJed Brown   CeedScalar e_total = U.E_total / U.density;
525c677226SJed Brown   CeedScalar e_internal = e_total - e_kinetic - e_potential;
535c677226SJed Brown   Y.temperature = e_internal / gas->cv;
545c677226SJed Brown   Y.pressure = (gas->cp / gas->cv - 1) * U.density * e_internal;
555c677226SJed Brown   return Y;
565c677226SJed Brown }
575c677226SJed Brown 
585c677226SJed Brown CEED_QFUNCTION_HELPER StatePrimitive StatePrimitiveFromConservative_fwd(
595c677226SJed Brown   NewtonianIdealGasContext gas, State s, StateConservative dU,
605c677226SJed Brown   const CeedScalar x[3], const CeedScalar dx[3]) {
615c677226SJed Brown   StatePrimitive dY;
625c677226SJed Brown   for (int i=0; i<3; i++) {
635c677226SJed Brown     dY.velocity[i] = (dU.momentum[i] - s.Y.velocity[i] * dU.density) / s.U.density;
645c677226SJed Brown   }
655c677226SJed Brown   CeedScalar e_kinetic = .5 * Dot3(s.Y.velocity, s.Y.velocity);
665c677226SJed Brown   CeedScalar de_kinetic = Dot3(dY.velocity, s.Y.velocity);
675c677226SJed Brown   CeedScalar e_potential = -Dot3(gas->g, x);
685c677226SJed Brown   CeedScalar de_potential = -Dot3(gas->g, dx);
695c677226SJed Brown   CeedScalar e_total = s.U.E_total / s.U.density;
705c677226SJed Brown   CeedScalar de_total = (dU.E_total - e_total * dU.density) / s.U.density;
715c677226SJed Brown   CeedScalar e_internal = e_total - e_kinetic - e_potential;
725c677226SJed Brown   CeedScalar de_internal = de_total - de_kinetic - de_potential;
735c677226SJed Brown   dY.temperature = de_internal / gas->cv;
745c677226SJed Brown   dY.pressure = (gas->cp / gas->cv - 1)
755c677226SJed Brown                 * (dU.density * e_internal + s.U.density * de_internal);
765c677226SJed Brown   return dY;
775c677226SJed Brown }
785c677226SJed Brown 
795c677226SJed Brown CEED_QFUNCTION_HELPER State StateFromU(NewtonianIdealGasContext gas,
805c677226SJed Brown                                        const CeedScalar U[5], const CeedScalar x[3]) {
815c677226SJed Brown   State s;
825c677226SJed Brown   s.U.density = U[0];
835c677226SJed Brown   s.U.momentum[0] = U[1];
845c677226SJed Brown   s.U.momentum[1] = U[2];
855c677226SJed Brown   s.U.momentum[2] = U[3];
865c677226SJed Brown   s.U.E_total = U[4];
875c677226SJed Brown   s.Y = StatePrimitiveFromConservative(gas, s.U, x);
885c677226SJed Brown   return s;
895c677226SJed Brown }
905c677226SJed Brown 
91*6f00d0e6SJed Brown CEED_QFUNCTION_HELPER State StateFromU_fwd(NewtonianIdealGasContext gas,
92*6f00d0e6SJed Brown     State s, const CeedScalar dU[5],
93*6f00d0e6SJed Brown     const CeedScalar x[3], const CeedScalar dx[3]) {
94*6f00d0e6SJed Brown   State ds;
95*6f00d0e6SJed Brown   ds.U.density = dU[0];
96*6f00d0e6SJed Brown   ds.U.momentum[0] = dU[1];
97*6f00d0e6SJed Brown   ds.U.momentum[1] = dU[2];
98*6f00d0e6SJed Brown   ds.U.momentum[2] = dU[3];
99*6f00d0e6SJed Brown   ds.U.E_total = dU[4];
100*6f00d0e6SJed Brown   ds.Y = StatePrimitiveFromConservative_fwd(gas, s, ds.U, x, dx);
101*6f00d0e6SJed Brown   return ds;
102*6f00d0e6SJed Brown }
103*6f00d0e6SJed Brown 
1045c677226SJed Brown CEED_QFUNCTION_HELPER void FluxInviscid(NewtonianIdealGasContext gas, State s,
1055c677226SJed Brown                                         StateConservative Flux[3]) {
1065c677226SJed Brown   for (int i=0; i<3; i++) {
1075c677226SJed Brown     Flux[i].density = s.U.momentum[i];
1085c677226SJed Brown     for (int j=0; j<3; j++)
1095c677226SJed Brown       Flux[i].momentum[j] = s.U.momentum[i] * s.Y.velocity[j]
1105c677226SJed Brown                             + s.Y.pressure * (i == j);
1115c677226SJed Brown     Flux[i].E_total = (s.U.E_total + s.Y.pressure) * s.Y.velocity[i];
1125c677226SJed Brown   }
1135c677226SJed Brown }
1145c677226SJed Brown 
1155c677226SJed Brown CEED_QFUNCTION_HELPER void FluxInviscid_fwd(NewtonianIdealGasContext gas,
1165c677226SJed Brown     State s, State ds, StateConservative dFlux[3]) {
1175c677226SJed Brown   for (int i=0; i<3; i++) {
1185c677226SJed Brown     dFlux[i].density = ds.U.momentum[i];
1195c677226SJed Brown     for (int j=0; j<3; j++)
1205c677226SJed Brown       dFlux[i].momentum[j] = ds.U.momentum[i] * s.Y.velocity[j] +
1215c677226SJed Brown                              s.U.momentum[i] * ds.Y.velocity[j] + ds.Y.pressure * (i == j);
1225c677226SJed Brown     dFlux[i].E_total = (ds.U.E_total + ds.Y.pressure) * s.Y.velocity[i] +
1235c677226SJed Brown                        (s.U.E_total + s.Y.pressure) * ds.Y.velocity[i];
1245c677226SJed Brown   }
1255c677226SJed Brown }
1265c677226SJed Brown 
1275c677226SJed Brown // Kelvin-Mandel notation
1285c677226SJed Brown CEED_QFUNCTION_HELPER void KMStrainRate(const State grad_s[3],
1295c677226SJed Brown                                         CeedScalar strain_rate[6]) {
1305c677226SJed Brown   const CeedScalar weight = 1 / sqrt(2.);
1315c677226SJed Brown   strain_rate[0] = grad_s[0].Y.velocity[0];
1325c677226SJed Brown   strain_rate[1] = grad_s[1].Y.velocity[1];
1335c677226SJed Brown   strain_rate[2] = grad_s[2].Y.velocity[2];
1345c677226SJed Brown   strain_rate[3] = weight * (grad_s[2].Y.velocity[1] + grad_s[1].Y.velocity[2]);
1355c677226SJed Brown   strain_rate[4] = weight * (grad_s[2].Y.velocity[0] + grad_s[0].Y.velocity[2]);
1365c677226SJed Brown   strain_rate[5] = weight * (grad_s[1].Y.velocity[0] + grad_s[0].Y.velocity[1]);
1375c677226SJed Brown }
1385c677226SJed Brown 
1395c677226SJed Brown CEED_QFUNCTION_HELPER void KMUnpack(const CeedScalar v[6], CeedScalar A[3][3]) {
1405c677226SJed Brown   const CeedScalar weight = 1 / sqrt(2.);
1415c677226SJed Brown   A[0][0] = v[0];
1425c677226SJed Brown   A[1][1] = v[1];
1435c677226SJed Brown   A[2][2] = v[2];
1445c677226SJed Brown   A[2][1] = A[1][2] = weight * v[3];
1455c677226SJed Brown   A[2][0] = A[0][2] = weight * v[4];
1465c677226SJed Brown   A[1][0] = A[0][1] = weight * v[5];
1475c677226SJed Brown }
1485c677226SJed Brown 
1495c677226SJed Brown CEED_QFUNCTION_HELPER void NewtonianStress(NewtonianIdealGasContext gas,
1505c677226SJed Brown     const CeedScalar strain_rate[6], CeedScalar stress[6]) {
1515c677226SJed Brown   CeedScalar div_u = strain_rate[0] + strain_rate[1] + strain_rate[2];
1525c677226SJed Brown   for (int i=0; i<6; i++) {
1535c677226SJed Brown     stress[i] = gas->mu * (2 * strain_rate[i] + gas->lambda * div_u * (i < 3));
1545c677226SJed Brown   }
1555c677226SJed Brown }
1565c677226SJed Brown 
1575c677226SJed Brown CEED_QFUNCTION_HELPER void ViscousEnergyFlux(NewtonianIdealGasContext gas,
1585c677226SJed Brown     StatePrimitive Y, const State grad_s[3], const CeedScalar stress[3][3],
1595c677226SJed Brown     CeedScalar Fe[3]) {
1605c677226SJed Brown   for (int i=0; i<3; i++) {
1615c677226SJed Brown     Fe[i] = - Y.velocity[0] * stress[0][i]
1625c677226SJed Brown             - Y.velocity[1] * stress[1][i]
1635c677226SJed Brown             - Y.velocity[2] * stress[2][i]
1645c677226SJed Brown             - gas->k * grad_s[i].Y.temperature;
1655c677226SJed Brown   }
1665c677226SJed Brown }
1675c677226SJed Brown 
16888b783a1SJames Wright // *****************************************************************************
16988b783a1SJames Wright // Helper function for computing flux Jacobian
17088b783a1SJames Wright // *****************************************************************************
17188b783a1SJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NS(CeedScalar dF[3][5][5],
17288b783a1SJames Wright     const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
17388626eedSJames Wright     const CeedScalar gamma, const CeedScalar g[3], const CeedScalar x[3]) {
17488b783a1SJames Wright   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square
17588626eedSJames Wright   CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]);
17688b783a1SJames Wright   for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions
17788b783a1SJames Wright     for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix
17888626eedSJames Wright       dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2. - e_potential)) : 0.) -
17988626eedSJames Wright                       u[i]*u[j];
18088b783a1SJames Wright       for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix
18188b783a1SJames Wright         dF[i][0][k+1]   = ((i==k) ? 1. : 0.);
18288b783a1SJames Wright         dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) +
18388b783a1SJames Wright                           ((i==k) ? u[j] : 0.) -
18488b783a1SJames Wright                           ((i==j) ? u[k] : 0.) * (gamma-1.);
18588b783a1SJames Wright         dF[i][4][k+1]   = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) -
18688b783a1SJames Wright                           (gamma-1.)*u[i]*u[k];
18788b783a1SJames Wright       }
18888b783a1SJames Wright       dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.);
18988b783a1SJames Wright     }
19088b783a1SJames Wright     dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho);
19188b783a1SJames Wright     dF[i][4][4] = u[i] * gamma;
19288b783a1SJames Wright   }
19388b783a1SJames Wright }
19488b783a1SJames Wright 
19588b783a1SJames Wright // *****************************************************************************
19688626eedSJames Wright // Helper function for computing flux Jacobian of Primitive variables
19788626eedSJames Wright // *****************************************************************************
19888626eedSJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NSp(CeedScalar dF[3][5][5],
19988626eedSJames Wright     const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
20088626eedSJames Wright     const CeedScalar Rd, const CeedScalar cv) {
20188626eedSJames Wright   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square
20288626eedSJames Wright   // TODO Add in gravity's contribution
20388626eedSJames Wright 
20488626eedSJames Wright   CeedScalar T    = ( E / rho - u_sq / 2. ) / cv;
20588626eedSJames Wright   CeedScalar drdT = -rho / T;
20688626eedSJames Wright   CeedScalar drdP = 1. / ( Rd * T);
20788626eedSJames Wright   CeedScalar etot =  E / rho ;
20888626eedSJames Wright   CeedScalar e2p  = drdP * etot + 1. ;
20988626eedSJames Wright   CeedScalar e3p  = ( E  + rho * Rd * T );
21088626eedSJames Wright   CeedScalar e4p  = drdT * etot + rho * cv ;
21188626eedSJames Wright 
21288626eedSJames Wright   for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions
21388626eedSJames Wright     for (CeedInt j=0; j<3; j++) { // j counts F^{m_j}
21488626eedSJames Wright //        [row][col] of A_i
21588626eedSJames Wright       dF[i][j+1][0] = drdP * u[i] * u[j] + ((i==j) ? 1. : 0.); // F^{{m_j} wrt p
21688626eedSJames Wright       for (CeedInt k=0; k<3; k++) { // k counts the wrt vel_k
217871db79fSKenneth E. Jansen         dF[i][0][k+1]   =  ((i==k) ? rho  : 0.);   // F^c wrt u_k
21888626eedSJames Wright         dF[i][j+1][k+1] = (((j==k) ? u[i] : 0.) +  // F^m_j wrt u_k
21988626eedSJames Wright                            ((i==k) ? u[j] : 0.) ) * rho;
22088626eedSJames Wright         dF[i][4][k+1]   = rho * u[i] * u[k]
22188626eedSJames Wright                           + ((i==k) ? e3p  : 0.) ; // F^e wrt u_k
22288626eedSJames Wright       }
22388626eedSJames Wright       dF[i][j+1][4] = drdT * u[i] * u[j]; // F^{m_j} wrt T
22488626eedSJames Wright     }
22588626eedSJames Wright     dF[i][4][0] = u[i] * e2p; // F^e wrt p
22688626eedSJames Wright     dF[i][4][4] = u[i] * e4p; // F^e wrt T
22788626eedSJames Wright     dF[i][0][0] = u[i] * drdP; // F^c wrt p
22888626eedSJames Wright     dF[i][0][4] = u[i] * drdT; // F^c wrt T
22988626eedSJames Wright   }
23088626eedSJames Wright }
23188626eedSJames Wright 
23288626eedSJames Wright CEED_QFUNCTION_HELPER void PrimitiveToConservative_fwd(const CeedScalar rho,
23388626eedSJames Wright     const CeedScalar u[3], const CeedScalar E, const CeedScalar Rd,
23488626eedSJames Wright     const CeedScalar cv, const CeedScalar dY[5], CeedScalar dU[5]) {
23588626eedSJames Wright   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2];
23688626eedSJames Wright   CeedScalar T    = ( E / rho - u_sq / 2. ) / cv;
23788626eedSJames Wright   CeedScalar drdT = -rho / T;
23888626eedSJames Wright   CeedScalar drdP = 1. / ( Rd * T);
23988626eedSJames Wright   dU[0] = drdP * dY[0] + drdT * dY[4];
24088626eedSJames Wright   CeedScalar de_kinetic = 0;
24188626eedSJames Wright   for (int i=0; i<3; i++) {
24288626eedSJames Wright     dU[1+i] = dU[0] * u[i] + rho * dY[1+i];
24388626eedSJames Wright     de_kinetic += u[i] * dY[1+i];
24488626eedSJames Wright   }
24588626eedSJames Wright   dU[4] = rho * cv * dY[4] + dU[0] * cv * T // internal energy: rho * e
24688626eedSJames Wright           + rho * de_kinetic + .5 * dU[0] * u_sq; // kinetic energy: .5 * rho * |u|^2
24788626eedSJames Wright }
24888626eedSJames Wright 
24988626eedSJames Wright // *****************************************************************************
25088626eedSJames Wright // Helper function for computing Tau elements (stabilization constant)
25188626eedSJames Wright //   Model from:
25288626eedSJames Wright //     PHASTA
25388626eedSJames Wright //
25488626eedSJames Wright //   Tau[i] = itau=0 which is diagonal-Shakib (3 values still but not spatial)
25588626eedSJames Wright //
25688626eedSJames Wright // Where NOT UPDATED YET
25788626eedSJames Wright // *****************************************************************************
25888626eedSJames Wright CEED_QFUNCTION_HELPER void Tau_diagPrim(CeedScalar Tau_d[3],
25988626eedSJames Wright                                         const CeedScalar dXdx[3][3], const CeedScalar u[3],
26088626eedSJames Wright                                         const CeedScalar cv, const NewtonianIdealGasContext newt_ctx,
26188626eedSJames Wright                                         const CeedScalar mu, const CeedScalar dt,
26288626eedSJames Wright                                         const CeedScalar rho) {
26388626eedSJames Wright   // Context
26488626eedSJames Wright   const CeedScalar Ctau_t = newt_ctx->Ctau_t;
26588626eedSJames Wright   const CeedScalar Ctau_v = newt_ctx->Ctau_v;
26688626eedSJames Wright   const CeedScalar Ctau_C = newt_ctx->Ctau_C;
26788626eedSJames Wright   const CeedScalar Ctau_M = newt_ctx->Ctau_M;
26888626eedSJames Wright   const CeedScalar Ctau_E = newt_ctx->Ctau_E;
26988626eedSJames Wright   CeedScalar gijd[6];
27088626eedSJames Wright   CeedScalar tau;
27188626eedSJames Wright   CeedScalar dts;
27288626eedSJames Wright   CeedScalar fact;
27388626eedSJames Wright 
27488626eedSJames Wright   //*INDENT-OFF*
27588626eedSJames Wright   gijd[0] =   dXdx[0][0] * dXdx[0][0]
27688626eedSJames Wright             + dXdx[1][0] * dXdx[1][0]
27788626eedSJames Wright             + dXdx[2][0] * dXdx[2][0];
27888626eedSJames Wright 
27988626eedSJames Wright   gijd[1] =   dXdx[0][0] * dXdx[0][1]
28088626eedSJames Wright             + dXdx[1][0] * dXdx[1][1]
28188626eedSJames Wright             + dXdx[2][0] * dXdx[2][1];
28288626eedSJames Wright 
28388626eedSJames Wright   gijd[2] =   dXdx[0][1] * dXdx[0][1]
28488626eedSJames Wright             + dXdx[1][1] * dXdx[1][1]
28588626eedSJames Wright             + dXdx[2][1] * dXdx[2][1];
28688626eedSJames Wright 
28788626eedSJames Wright   gijd[3] =   dXdx[0][0] * dXdx[0][2]
28888626eedSJames Wright             + dXdx[1][0] * dXdx[1][2]
28988626eedSJames Wright             + dXdx[2][0] * dXdx[2][2];
29088626eedSJames Wright 
29188626eedSJames Wright   gijd[4] =   dXdx[0][1] * dXdx[0][2]
29288626eedSJames Wright             + dXdx[1][1] * dXdx[1][2]
29388626eedSJames Wright             + dXdx[2][1] * dXdx[2][2];
29488626eedSJames Wright 
29588626eedSJames Wright   gijd[5] =   dXdx[0][2] * dXdx[0][2]
29688626eedSJames Wright             + dXdx[1][2] * dXdx[1][2]
29788626eedSJames Wright             + dXdx[2][2] * dXdx[2][2];
29888626eedSJames Wright   //*INDENT-ON*
29988626eedSJames Wright 
30088626eedSJames Wright   dts = Ctau_t / dt ;
30188626eedSJames Wright 
30288626eedSJames Wright   tau = rho*rho*((4. * dts * dts)
30388626eedSJames Wright                  + u[0] * ( u[0] * gijd[0] + 2. * ( u[1] * gijd[1] + u[2] * gijd[3]))
30488626eedSJames Wright                  + u[1] * ( u[1] * gijd[2] + 2. *   u[2] * gijd[4])
30588626eedSJames Wright                  + u[2] *   u[2] * gijd[5])
30688626eedSJames Wright         + Ctau_v* mu * mu *
30788626eedSJames Wright         (gijd[0]*gijd[0] + gijd[2]*gijd[2] + gijd[5]*gijd[5] +
30888626eedSJames Wright          + 2. * (gijd[1]*gijd[1] + gijd[3]*gijd[3] + gijd[4]*gijd[4]));
30988626eedSJames Wright 
31088626eedSJames Wright   fact=sqrt(tau);
31188626eedSJames Wright 
31288626eedSJames Wright   Tau_d[0] = Ctau_C * fact / (rho*(gijd[0] + gijd[2] + gijd[5]))*0.125;
31388626eedSJames Wright 
31488626eedSJames Wright   Tau_d[1] = Ctau_M / fact;
31588626eedSJames Wright   Tau_d[2] = Ctau_E / ( fact * cv );
31688626eedSJames Wright 
31788626eedSJames Wright // consider putting back the way I initially had it  Ctau_E * Tau_d[1] /cv
31888626eedSJames Wright //  to avoid a division if the compiler is smart enough to see that cv IS
31988626eedSJames Wright // a constant that it could invert once for all elements
32088626eedSJames Wright // but in that case energy tau is scaled by the product of Ctau_E * Ctau_M
32188626eedSJames Wright // OR we could absorb cv into Ctau_E but this puts more burden on user to
32288626eedSJames Wright // know how to change constants with a change of fluid or units.  Same for
32388626eedSJames Wright // Ctau_v * mu * mu IF AND ONLY IF we don't add viscosity law =f(T)
32488626eedSJames Wright }
32588626eedSJames Wright 
32688626eedSJames Wright // *****************************************************************************
32788b783a1SJames Wright // This QFunction sets a "still" initial condition for generic Newtonian IG problems
32888b783a1SJames Wright // *****************************************************************************
32988b783a1SJames Wright CEED_QFUNCTION(ICsNewtonianIG)(void *ctx, CeedInt Q,
33088b783a1SJames Wright                                const CeedScalar *const *in, CeedScalar *const *out) {
33188b783a1SJames Wright   // Inputs
33288b783a1SJames Wright   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
33388b783a1SJames Wright 
33488b783a1SJames Wright   // Outputs
33588b783a1SJames Wright   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
33688b783a1SJames Wright 
33788626eedSJames Wright   // Context
33888626eedSJames Wright   const SetupContext context = (SetupContext)ctx;
33988626eedSJames Wright   const CeedScalar theta0    = context->theta0;
34088626eedSJames Wright   const CeedScalar P0        = context->P0;
34188626eedSJames Wright   const CeedScalar cv        = context->cv;
34288626eedSJames Wright   const CeedScalar cp        = context->cp;
34388626eedSJames Wright   const CeedScalar *g        = context->g;
34488626eedSJames Wright   const CeedScalar Rd        = cp - cv;
34588626eedSJames Wright 
34688b783a1SJames Wright   // Quadrature Point Loop
34788b783a1SJames Wright   CeedPragmaSIMD
34888b783a1SJames Wright   for (CeedInt i=0; i<Q; i++) {
34988b783a1SJames Wright     CeedScalar q[5] = {0.};
35088b783a1SJames Wright 
35188b783a1SJames Wright     // Setup
35288b783a1SJames Wright     // -- Coordinates
35388626eedSJames Wright     const CeedScalar x[3] = {X[0][i], X[1][i], X[2][i]};
35488626eedSJames Wright     const CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]);
35588b783a1SJames Wright 
35688b783a1SJames Wright     // -- Density
35788626eedSJames Wright     const CeedScalar rho = P0 / (Rd*theta0);
35888b783a1SJames Wright 
35988b783a1SJames Wright     // Initial Conditions
36088b783a1SJames Wright     q[0] = rho;
36188b783a1SJames Wright     q[1] = 0.0;
36288b783a1SJames Wright     q[2] = 0.0;
36388b783a1SJames Wright     q[3] = 0.0;
36488626eedSJames Wright     q[4] = rho * (cv*theta0 + e_potential);
36588b783a1SJames Wright 
36688b783a1SJames Wright     for (CeedInt j=0; j<5; j++)
36788b783a1SJames Wright       q0[j][i] = q[j];
36888b783a1SJames Wright   } // End of Quadrature Point Loop
36988b783a1SJames Wright   return 0;
37088b783a1SJames Wright }
37188b783a1SJames Wright 
37288b783a1SJames Wright // *****************************************************************************
37388b783a1SJames Wright // This QFunction implements the following formulation of Navier-Stokes with
37488b783a1SJames Wright //   explicit time stepping method
37588b783a1SJames Wright //
37688b783a1SJames Wright // This is 3D compressible Navier-Stokes in conservation form with state
37788b783a1SJames Wright //   variables of density, momentum density, and total energy density.
37888b783a1SJames Wright //
37988b783a1SJames Wright // State Variables: q = ( rho, U1, U2, U3, E )
38088b783a1SJames Wright //   rho - Mass Density
38188b783a1SJames Wright //   Ui  - Momentum Density,      Ui = rho ui
38288b783a1SJames Wright //   E   - Total Energy Density,  E  = rho (cv T + (u u)/2 + g z)
38388b783a1SJames Wright //
38488b783a1SJames Wright // Navier-Stokes Equations:
38588b783a1SJames Wright //   drho/dt + div( U )                               = 0
38688b783a1SJames Wright //   dU/dt   + div( rho (u x u) + P I3 ) + rho g khat = div( Fu )
38788b783a1SJames Wright //   dE/dt   + div( (E + P) u )                       = div( Fe )
38888b783a1SJames Wright //
38988b783a1SJames Wright // Viscous Stress:
39088b783a1SJames Wright //   Fu = mu (grad( u ) + grad( u )^T + lambda div ( u ) I3)
39188b783a1SJames Wright //
39288b783a1SJames Wright // Thermal Stress:
39388b783a1SJames Wright //   Fe = u Fu + k grad( T )
39488626eedSJames Wright // Equation of State
39588b783a1SJames Wright //   P = (gamma - 1) (E - rho (u u) / 2 - rho g z)
39688b783a1SJames Wright //
39788b783a1SJames Wright // Stabilization:
39888b783a1SJames Wright //   Tau = diag(TauC, TauM, TauM, TauM, TauE)
39988b783a1SJames Wright //     f1 = rho  sqrt(ui uj gij)
40088b783a1SJames Wright //     gij = dXi/dX * dXi/dX
40188b783a1SJames Wright //     TauC = Cc f1 / (8 gii)
40288b783a1SJames Wright //     TauM = min( 1 , 1 / f1 )
40388b783a1SJames Wright //     TauE = TauM / (Ce cv)
40488b783a1SJames Wright //
40588b783a1SJames Wright //  SU   = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) )
40688b783a1SJames Wright //
40788b783a1SJames Wright // Constants:
40888b783a1SJames Wright //   lambda = - 2 / 3,  From Stokes hypothesis
40988b783a1SJames Wright //   mu              ,  Dynamic viscosity
41088b783a1SJames Wright //   k               ,  Thermal conductivity
41188b783a1SJames Wright //   cv              ,  Specific heat, constant volume
41288b783a1SJames Wright //   cp              ,  Specific heat, constant pressure
41388b783a1SJames Wright //   g               ,  Gravity
41488b783a1SJames Wright //   gamma  = cp / cv,  Specific heat ratio
41588b783a1SJames Wright //
41688b783a1SJames Wright // We require the product of the inverse of the Jacobian (dXdx_j,k) and
41788b783a1SJames Wright // its transpose (dXdx_k,j) to properly compute integrals of the form:
41888b783a1SJames Wright // int( gradv gradu )
41988b783a1SJames Wright //
42088b783a1SJames Wright // *****************************************************************************
4215c677226SJed Brown CEED_QFUNCTION(RHSFunction_Newtonian)(void *ctx, CeedInt Q,
42288b783a1SJames Wright                                       const CeedScalar *const *in, CeedScalar *const *out) {
42388b783a1SJames Wright   // *INDENT-OFF*
42488b783a1SJames Wright   // Inputs
42588b783a1SJames Wright   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
42688b783a1SJames Wright                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
42788b783a1SJames Wright                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
42888b783a1SJames Wright                    (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
42988b783a1SJames Wright   // Outputs
43088b783a1SJames Wright   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
43188b783a1SJames Wright              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
43288b783a1SJames Wright   // *INDENT-ON*
43388b783a1SJames Wright 
43488b783a1SJames Wright   // Context
43588b783a1SJames Wright   NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx;
43688b783a1SJames Wright   const CeedScalar mu     = context->mu;
43788b783a1SJames Wright   const CeedScalar cv     = context->cv;
43888b783a1SJames Wright   const CeedScalar cp     = context->cp;
43988626eedSJames Wright   const CeedScalar *g     = context->g;
44088626eedSJames Wright   const CeedScalar dt     = context->dt;
44188b783a1SJames Wright   const CeedScalar gamma  = cp / cv;
44288626eedSJames Wright   const CeedScalar Rd     = cp - cv;
44388b783a1SJames Wright 
44488b783a1SJames Wright   CeedPragmaSIMD
44588b783a1SJames Wright   // Quadrature Point Loop
44688b783a1SJames Wright   for (CeedInt i=0; i<Q; i++) {
4475c677226SJed Brown     CeedScalar U[5];
4485c677226SJed Brown     for (int j=0; j<5; j++) U[j] = q[j][i];
4495c677226SJed Brown     const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]};
4505c677226SJed Brown     State s = StateFromU(context, U, x_i);
4515c677226SJed Brown 
45288b783a1SJames Wright     // -- Interp-to-Interp q_data
45388b783a1SJames Wright     const CeedScalar wdetJ      =   q_data[0][i];
45488b783a1SJames Wright     // -- Interp-to-Grad q_data
45588b783a1SJames Wright     // ---- Inverse of change of coordinate matrix: X_i,j
45688b783a1SJames Wright     // *INDENT-OFF*
45788b783a1SJames Wright     const CeedScalar dXdx[3][3] = {{q_data[1][i],
45888b783a1SJames Wright                                     q_data[2][i],
45988b783a1SJames Wright                                     q_data[3][i]},
46088b783a1SJames Wright                                    {q_data[4][i],
46188b783a1SJames Wright                                     q_data[5][i],
46288b783a1SJames Wright                                     q_data[6][i]},
46388b783a1SJames Wright                                    {q_data[7][i],
46488b783a1SJames Wright                                     q_data[8][i],
46588b783a1SJames Wright                                     q_data[9][i]}
46688b783a1SJames Wright                                   };
46788b783a1SJames Wright     // *INDENT-ON*
46888b783a1SJames Wright 
4695c677226SJed Brown     State grad_s[3];
4705c677226SJed Brown     for (int j=0; j<3; j++) {
471*6f00d0e6SJed Brown       CeedScalar dx_i[3] = {0}, dU[5];
472*6f00d0e6SJed Brown       for (int k=0; k<5; k++) dU[k] = dq[0][k][i] * dXdx[0][j]
473*6f00d0e6SJed Brown                                         + dq[1][k][i] * dXdx[1][j]
474*6f00d0e6SJed Brown                                         + dq[2][k][i] * dXdx[2][j];
4755c677226SJed Brown       dx_i[j] = 1.;
476*6f00d0e6SJed Brown       grad_s[j] = StateFromU_fwd(context, s, dU, x_i, dx_i);
4775c677226SJed Brown     }
4785c677226SJed Brown 
4795c677226SJed Brown     CeedScalar strain_rate[6], kmstress[6], stress[3][3], Fe[3];
4805c677226SJed Brown     KMStrainRate(grad_s, strain_rate);
4815c677226SJed Brown     NewtonianStress(context, strain_rate, kmstress);
4825c677226SJed Brown     KMUnpack(kmstress, stress);
4835c677226SJed Brown     ViscousEnergyFlux(context, s.Y, grad_s, stress, Fe);
4845c677226SJed Brown 
4855c677226SJed Brown     StateConservative F_inviscid[3];
4865c677226SJed Brown     FluxInviscid(context, s, F_inviscid);
4875c677226SJed Brown 
4885c677226SJed Brown     // Total flux
4895c677226SJed Brown     CeedScalar Flux[5][3];
4905c677226SJed Brown     for (int j=0; j<3; j++) {
4915c677226SJed Brown       Flux[0][j] = F_inviscid[j].density;
4925c677226SJed Brown       for (int k=0; k<3; k++)
4935c677226SJed Brown         Flux[k+1][j] = F_inviscid[j].momentum[k] - stress[k][j];
4945c677226SJed Brown       Flux[4][j] = F_inviscid[j].E_total + Fe[j];
4955c677226SJed Brown     }
4965c677226SJed Brown 
4975c677226SJed Brown     for (int j=0; j<3; j++) {
4985c677226SJed Brown       for (int k=0; k<5; k++) {
4995c677226SJed Brown         dv[j][k][i] = wdetJ * (dXdx[j][0] * Flux[k][0] +
5005c677226SJed Brown                                dXdx[j][1] * Flux[k][1] +
5015c677226SJed Brown                                dXdx[j][2] * Flux[k][2]);
5025c677226SJed Brown       }
5035c677226SJed Brown     }
5045c677226SJed Brown 
5055c677226SJed Brown     const CeedScalar body_force[5] = {0, s.U.density *g[0], s.U.density *g[1], s.U.density *g[2], 0};
5065c677226SJed Brown     for (int j=0; j<5; j++)
5075c677226SJed Brown       v[j][i] = wdetJ * body_force[j];
50888b783a1SJames Wright 
50988b783a1SJames Wright     // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
5105c677226SJed Brown     CeedScalar jacob_F_conv[3][5][5] = {0};
5115c677226SJed Brown     computeFluxJacobian_NS(jacob_F_conv, s.U.density, s.Y.velocity, s.U.E_total,
5125c677226SJed Brown                            gamma, g, x_i);
5135c677226SJed Brown     CeedScalar grad_U[5][3];
51488b783a1SJames Wright     for (int j=0; j<3; j++) {
5155c677226SJed Brown       grad_U[0][j] = grad_s[j].U.density;
5165c677226SJed Brown       for (int k=0; k<3; k++) grad_U[k+1][j] = grad_s[j].U.momentum[k];
5175c677226SJed Brown       grad_U[4][j] = grad_s[j].U.E_total;
51888b783a1SJames Wright     }
51988b783a1SJames Wright 
52088b783a1SJames Wright     // strong_conv = dF/dq * dq/dx    (Strong convection)
52188b783a1SJames Wright     CeedScalar strong_conv[5] = {0};
52288b783a1SJames Wright     for (int j=0; j<3; j++)
52388b783a1SJames Wright       for (int k=0; k<5; k++)
52488b783a1SJames Wright         for (int l=0; l<5; l++)
5255c677226SJed Brown           strong_conv[k] += jacob_F_conv[j][k][l] * grad_U[l][j];
52688b783a1SJames Wright 
52788626eedSJames Wright     // -- Stabilization method: none, SU, or SUPG
52888626eedSJames Wright     CeedScalar stab[5][3] = {{0.}};
52988626eedSJames Wright     CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0};
53088626eedSJames Wright     CeedScalar Tau_d[3] = {0.};
53188b783a1SJames Wright     switch (context->stabilization) {
53288b783a1SJames Wright     case STAB_NONE:        // Galerkin
53388b783a1SJames Wright       break;
53488b783a1SJames Wright     case STAB_SU:        // SU
5355c677226SJed Brown       Tau_diagPrim(Tau_d, dXdx, s.Y.velocity, cv, context, mu, dt, s.U.density);
53688626eedSJames Wright       tau_strong_conv[0] = Tau_d[0] * strong_conv[0];
53788626eedSJames Wright       tau_strong_conv[1] = Tau_d[1] * strong_conv[1];
53888626eedSJames Wright       tau_strong_conv[2] = Tau_d[1] * strong_conv[2];
53988626eedSJames Wright       tau_strong_conv[3] = Tau_d[1] * strong_conv[3];
54088626eedSJames Wright       tau_strong_conv[4] = Tau_d[2] * strong_conv[4];
5415c677226SJed Brown       PrimitiveToConservative_fwd(s.U.density, s.Y.velocity, s.U.E_total, Rd, cv,
5425c677226SJed Brown                                   tau_strong_conv,
54388626eedSJames Wright                                   tau_strong_conv_conservative);
54488b783a1SJames Wright       for (int j=0; j<3; j++)
54588b783a1SJames Wright         for (int k=0; k<5; k++)
54688b783a1SJames Wright           for (int l=0; l<5; l++)
54788626eedSJames Wright             stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l];
54888b783a1SJames Wright 
54988b783a1SJames Wright       for (int j=0; j<5; j++)
55088b783a1SJames Wright         for (int k=0; k<3; k++)
55188b783a1SJames Wright           dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] +
55288b783a1SJames Wright                                 stab[j][1] * dXdx[k][1] +
55388b783a1SJames Wright                                 stab[j][2] * dXdx[k][2]);
55488b783a1SJames Wright       break;
55588b783a1SJames Wright     case STAB_SUPG:        // SUPG is not implemented for explicit scheme
55688b783a1SJames Wright       break;
55788b783a1SJames Wright     }
55888b783a1SJames Wright 
55988b783a1SJames Wright   } // End Quadrature Point Loop
56088b783a1SJames Wright 
56188b783a1SJames Wright   // Return
56288b783a1SJames Wright   return 0;
56388b783a1SJames Wright }
56488b783a1SJames Wright 
56588b783a1SJames Wright // *****************************************************************************
56688b783a1SJames Wright // This QFunction implements the Navier-Stokes equations (mentioned above) with
56788b783a1SJames Wright //   implicit time stepping method
56888b783a1SJames Wright //
56988b783a1SJames Wright //  SU   = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) )
57088b783a1SJames Wright //  SUPG = Galerkin + grad(v) . ( Ai^T * Tau * (q_dot + Aj q,j - body force) )
57188b783a1SJames Wright //                                       (diffussive terms will be added later)
57288b783a1SJames Wright //
57388b783a1SJames Wright // *****************************************************************************
57488b783a1SJames Wright CEED_QFUNCTION(IFunction_Newtonian)(void *ctx, CeedInt Q,
57588b783a1SJames Wright                                     const CeedScalar *const *in,
57688b783a1SJames Wright                                     CeedScalar *const *out) {
57788b783a1SJames Wright   // *INDENT-OFF*
57888b783a1SJames Wright   // Inputs
57988b783a1SJames Wright   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
58088b783a1SJames Wright                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
58188b783a1SJames Wright                    (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
58288b783a1SJames Wright                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3],
58388b783a1SJames Wright                    (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4];
58488b783a1SJames Wright   // Outputs
58588b783a1SJames Wright   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
58688b783a1SJames Wright              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
58788b783a1SJames Wright   // *INDENT-ON*
58888b783a1SJames Wright   // Context
58988b783a1SJames Wright   NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx;
59088b783a1SJames Wright   const CeedScalar mu     = context->mu;
59188b783a1SJames Wright   const CeedScalar cv     = context->cv;
59288b783a1SJames Wright   const CeedScalar cp     = context->cp;
59388626eedSJames Wright   const CeedScalar *g     = context->g;
59488626eedSJames Wright   const CeedScalar dt     = context->dt;
59588b783a1SJames Wright   const CeedScalar gamma  = cp / cv;
59688626eedSJames Wright   const CeedScalar Rd     = cp-cv;
59788b783a1SJames Wright 
59888b783a1SJames Wright   CeedPragmaSIMD
59988b783a1SJames Wright   // Quadrature Point Loop
60088b783a1SJames Wright   for (CeedInt i=0; i<Q; i++) {
6015c677226SJed Brown     CeedScalar U[5];
6025c677226SJed Brown     for (int j=0; j<5; j++) U[j] = q[j][i];
6035c677226SJed Brown     const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]};
6045c677226SJed Brown     State s = StateFromU(context, U, x_i);
6055c677226SJed Brown 
60688b783a1SJames Wright     // -- Interp-to-Interp q_data
60788b783a1SJames Wright     const CeedScalar wdetJ      =   q_data[0][i];
60888b783a1SJames Wright     // -- Interp-to-Grad q_data
60988b783a1SJames Wright     // ---- Inverse of change of coordinate matrix: X_i,j
61088b783a1SJames Wright     // *INDENT-OFF*
61188b783a1SJames Wright     const CeedScalar dXdx[3][3] = {{q_data[1][i],
61288b783a1SJames Wright                                     q_data[2][i],
61388b783a1SJames Wright                                     q_data[3][i]},
61488b783a1SJames Wright                                    {q_data[4][i],
61588b783a1SJames Wright                                     q_data[5][i],
61688b783a1SJames Wright                                     q_data[6][i]},
61788b783a1SJames Wright                                    {q_data[7][i],
61888b783a1SJames Wright                                     q_data[8][i],
61988b783a1SJames Wright                                     q_data[9][i]}
62088b783a1SJames Wright                                   };
62188b783a1SJames Wright     // *INDENT-ON*
6225c677226SJed Brown     State grad_s[3];
62388b783a1SJames Wright     for (int j=0; j<3; j++) {
624*6f00d0e6SJed Brown       CeedScalar dx_i[3] = {0}, dU[5];
625*6f00d0e6SJed Brown       for (int k=0; k<5; k++) dU[k] = dq[0][k][i] * dXdx[0][j]
626*6f00d0e6SJed Brown                                         + dq[1][k][i] * dXdx[1][j]
627*6f00d0e6SJed Brown                                         + dq[2][k][i] * dXdx[2][j];
6285c677226SJed Brown       dx_i[j] = 1.;
629*6f00d0e6SJed Brown       grad_s[j] = StateFromU_fwd(context, s, dU, x_i, dx_i);
63088b783a1SJames Wright     }
6315c677226SJed Brown 
6325c677226SJed Brown     CeedScalar strain_rate[6], kmstress[6], stress[3][3], Fe[3];
6335c677226SJed Brown     KMStrainRate(grad_s, strain_rate);
6345c677226SJed Brown     NewtonianStress(context, strain_rate, kmstress);
6355c677226SJed Brown     KMUnpack(kmstress, stress);
6365c677226SJed Brown     ViscousEnergyFlux(context, s.Y, grad_s, stress, Fe);
6375c677226SJed Brown 
6385c677226SJed Brown     StateConservative F_inviscid[3];
6395c677226SJed Brown     FluxInviscid(context, s, F_inviscid);
6405c677226SJed Brown 
6415c677226SJed Brown 
6425c677226SJed Brown     // Total flux
6435c677226SJed Brown     CeedScalar Flux[5][3];
6445c677226SJed Brown     for (int j=0; j<3; j++) {
6455c677226SJed Brown       Flux[0][j] = F_inviscid[j].density;
64688b783a1SJames Wright       for (int k=0; k<3; k++)
6475c677226SJed Brown         Flux[k+1][j] = F_inviscid[j].momentum[k] - stress[k][j];
6485c677226SJed Brown       Flux[4][j] = F_inviscid[j].E_total + Fe[j];
6495c677226SJed Brown     }
6505c677226SJed Brown 
6515c677226SJed Brown     for (int j=0; j<3; j++) {
6525c677226SJed Brown       for (int k=0; k<5; k++) {
6535c677226SJed Brown         dv[j][k][i] = -wdetJ * (dXdx[j][0] * Flux[k][0] +
6545c677226SJed Brown                                 dXdx[j][1] * Flux[k][1] +
6555c677226SJed Brown                                 dXdx[j][2] * Flux[k][2]);
6565c677226SJed Brown       }
6575c677226SJed Brown     }
6585c677226SJed Brown 
6595c677226SJed Brown     const CeedScalar body_force[5] = {0, s.U.density *g[0], s.U.density *g[1], s.U.density *g[2], 0};
6605c677226SJed Brown     for (int j=0; j<5; j++)
6615c677226SJed Brown       v[j][i] = wdetJ * (q_dot[j][i] - body_force[j]);
66288b783a1SJames Wright 
66388b783a1SJames Wright     // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
6645c677226SJed Brown     CeedScalar jacob_F_conv[3][5][5] = {0};
6655c677226SJed Brown     computeFluxJacobian_NS(jacob_F_conv, s.U.density, s.Y.velocity, s.U.E_total,
6665c677226SJed Brown                            gamma, g, x_i);
6675c677226SJed Brown     CeedScalar grad_U[5][3];
66888b783a1SJames Wright     for (int j=0; j<3; j++) {
6695c677226SJed Brown       grad_U[0][j] = grad_s[j].U.density;
6705c677226SJed Brown       for (int k=0; k<3; k++) grad_U[k+1][j] = grad_s[j].U.momentum[k];
6715c677226SJed Brown       grad_U[4][j] = grad_s[j].U.E_total;
67288b783a1SJames Wright     }
6735c677226SJed Brown 
67488b783a1SJames Wright     // strong_conv = dF/dq * dq/dx    (Strong convection)
67588b783a1SJames Wright     CeedScalar strong_conv[5] = {0};
67688b783a1SJames Wright     for (int j=0; j<3; j++)
67788b783a1SJames Wright       for (int k=0; k<5; k++)
67888b783a1SJames Wright         for (int l=0; l<5; l++)
6795c677226SJed Brown           strong_conv[k] += jacob_F_conv[j][k][l] * grad_U[l][j];
68088b783a1SJames Wright 
68188b783a1SJames Wright     // Strong residual
68288b783a1SJames Wright     CeedScalar strong_res[5];
68388b783a1SJames Wright     for (int j=0; j<5; j++)
68488b783a1SJames Wright       strong_res[j] = q_dot[j][i] + strong_conv[j] - body_force[j];
68588b783a1SJames Wright 
68688b783a1SJames Wright     // -- Stabilization method: none, SU, or SUPG
68788626eedSJames Wright     CeedScalar stab[5][3] = {{0.}};
68888626eedSJames Wright     CeedScalar tau_strong_res[5] = {0.}, tau_strong_res_conservative[5] = {0};
68988626eedSJames Wright     CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0};
69088626eedSJames Wright     CeedScalar Tau_d[3] = {0.};
69188b783a1SJames Wright     switch (context->stabilization) {
69288b783a1SJames Wright     case STAB_NONE:        // Galerkin
69388b783a1SJames Wright       break;
69488b783a1SJames Wright     case STAB_SU:        // SU
6955c677226SJed Brown       Tau_diagPrim(Tau_d, dXdx, s.Y.velocity, cv, context, mu, dt, s.U.density);
69688626eedSJames Wright       tau_strong_conv[0] = Tau_d[0] * strong_conv[0];
69788626eedSJames Wright       tau_strong_conv[1] = Tau_d[1] * strong_conv[1];
69888626eedSJames Wright       tau_strong_conv[2] = Tau_d[1] * strong_conv[2];
69988626eedSJames Wright       tau_strong_conv[3] = Tau_d[1] * strong_conv[3];
70088626eedSJames Wright       tau_strong_conv[4] = Tau_d[2] * strong_conv[4];
7015c677226SJed Brown       PrimitiveToConservative_fwd(s.U.density, s.Y.velocity, s.U.E_total, Rd, cv,
7025c677226SJed Brown                                   tau_strong_conv, tau_strong_conv_conservative);
70388b783a1SJames Wright       for (int j=0; j<3; j++)
70488b783a1SJames Wright         for (int k=0; k<5; k++)
70588b783a1SJames Wright           for (int l=0; l<5; l++)
70688626eedSJames Wright             stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l];
70788b783a1SJames Wright 
70888b783a1SJames Wright       for (int j=0; j<5; j++)
70988b783a1SJames Wright         for (int k=0; k<3; k++)
71088b783a1SJames Wright           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
71188b783a1SJames Wright                                 stab[j][1] * dXdx[k][1] +
71288b783a1SJames Wright                                 stab[j][2] * dXdx[k][2]);
71388b783a1SJames Wright       break;
71488b783a1SJames Wright     case STAB_SUPG:        // SUPG
7155c677226SJed Brown       Tau_diagPrim(Tau_d, dXdx, s.Y.velocity, cv, context, mu, dt, s.U.density);
71688626eedSJames Wright       tau_strong_res[0] = Tau_d[0] * strong_res[0];
71788626eedSJames Wright       tau_strong_res[1] = Tau_d[1] * strong_res[1];
71888626eedSJames Wright       tau_strong_res[2] = Tau_d[1] * strong_res[2];
71988626eedSJames Wright       tau_strong_res[3] = Tau_d[1] * strong_res[3];
72088626eedSJames Wright       tau_strong_res[4] = Tau_d[2] * strong_res[4];
72188626eedSJames Wright // Alternate route (useful later with primitive variable code)
72288626eedSJames Wright // this function was verified against PHASTA for as IC that was as close as possible
72388626eedSJames Wright //    computeFluxJacobian_NSp(jacob_F_conv_p, rho, u, E, Rd, cv);
72488626eedSJames Wright // it has also been verified to compute a correct through the following
72588626eedSJames Wright //   stab[k][j] += jacob_F_conv_p[j][k][l] * tau_strong_res[l] // flux Jacobian wrt primitive
72688626eedSJames Wright // applied in the triple loop below
72788626eedSJames Wright //  However, it is more flops than using the existing Jacobian wrt q after q_{,Y} viz
7285c677226SJed Brown       PrimitiveToConservative_fwd(s.U.density, s.Y.velocity, s.U.E_total, Rd, cv,
7295c677226SJed Brown                                   tau_strong_res, tau_strong_res_conservative);
73088b783a1SJames Wright       for (int j=0; j<3; j++)
73188b783a1SJames Wright         for (int k=0; k<5; k++)
73288b783a1SJames Wright           for (int l=0; l<5; l++)
73388626eedSJames Wright             stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_res_conservative[l];
73488b783a1SJames Wright 
73588b783a1SJames Wright       for (int j=0; j<5; j++)
73688b783a1SJames Wright         for (int k=0; k<3; k++)
73788b783a1SJames Wright           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
73888b783a1SJames Wright                                 stab[j][1] * dXdx[k][1] +
73988b783a1SJames Wright                                 stab[j][2] * dXdx[k][2]);
74088b783a1SJames Wright       break;
74188b783a1SJames Wright     }
74288b783a1SJames Wright 
74388b783a1SJames Wright   } // End Quadrature Point Loop
74488b783a1SJames Wright 
74588b783a1SJames Wright   // Return
74688b783a1SJames Wright   return 0;
74788b783a1SJames Wright }
74888b783a1SJames Wright // *****************************************************************************
74988b783a1SJames Wright #endif // newtonian_h
750