xref: /honee/qfunctions/newtonian.h (revision 752f40e3cf4e9e12ff9fff9133a2bd9a934f2aaf)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
33a8779fbSJames Wright //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
53a8779fbSJames Wright //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73a8779fbSJames Wright 
83a8779fbSJames Wright /// @file
93a8779fbSJames Wright /// Operator for Navier-Stokes example using PETSc
103a8779fbSJames Wright 
113a8779fbSJames Wright 
123a8779fbSJames Wright #ifndef newtonian_h
133a8779fbSJames Wright #define newtonian_h
143a8779fbSJames Wright 
153a8779fbSJames Wright #include <math.h>
163a8779fbSJames Wright #include <ceed.h>
1715a3537eSJed Brown #include "newtonian_types.h"
183a8779fbSJames Wright 
193a8779fbSJames Wright #ifndef M_PI
203a8779fbSJames Wright #define M_PI    3.14159265358979323846
213a8779fbSJames Wright #endif
223a8779fbSJames Wright 
23c1a52365SJed Brown typedef struct {
24c1a52365SJed Brown   CeedScalar pressure;
25c1a52365SJed Brown   CeedScalar velocity[3];
26c1a52365SJed Brown   CeedScalar temperature;
27c1a52365SJed Brown } StatePrimitive;
28c1a52365SJed Brown 
29c1a52365SJed Brown typedef struct {
30c1a52365SJed Brown   CeedScalar density;
31c1a52365SJed Brown   CeedScalar momentum[3];
32c1a52365SJed Brown   CeedScalar E_total;
33c1a52365SJed Brown } StateConservative;
34c1a52365SJed Brown 
35c1a52365SJed Brown typedef struct {
36c1a52365SJed Brown   StateConservative U;
37c1a52365SJed Brown   StatePrimitive Y;
38c1a52365SJed Brown } State;
39c1a52365SJed Brown 
40c1a52365SJed Brown CEED_QFUNCTION_HELPER CeedScalar Dot3(const CeedScalar u[3],
41c1a52365SJed Brown                                       const CeedScalar v[3]) {
42c1a52365SJed Brown   return u[0]*v[0] + u[1]*v[1] + u[2]*v[2];
43c1a52365SJed Brown }
44c1a52365SJed Brown 
45c1a52365SJed Brown CEED_QFUNCTION_HELPER StatePrimitive StatePrimitiveFromConservative(
46c1a52365SJed Brown   NewtonianIdealGasContext gas, StateConservative U, const CeedScalar x[3]) {
47c1a52365SJed Brown   StatePrimitive Y;
48c1a52365SJed Brown   for (int i=0; i<3; i++) Y.velocity[i] = U.momentum[i] / U.density;
49c1a52365SJed Brown   CeedScalar e_kinetic = .5 * Dot3(Y.velocity, Y.velocity);
50c1a52365SJed Brown   CeedScalar e_potential = -Dot3(gas->g, x);
51c1a52365SJed Brown   CeedScalar e_total = U.E_total / U.density;
52c1a52365SJed Brown   CeedScalar e_internal = e_total - e_kinetic - e_potential;
53c1a52365SJed Brown   Y.temperature = e_internal / gas->cv;
54c1a52365SJed Brown   Y.pressure = (gas->cp / gas->cv - 1) * U.density * e_internal;
55c1a52365SJed Brown   return Y;
56c1a52365SJed Brown }
57c1a52365SJed Brown 
58c1a52365SJed Brown CEED_QFUNCTION_HELPER StatePrimitive StatePrimitiveFromConservative_fwd(
59c1a52365SJed Brown   NewtonianIdealGasContext gas, State s, StateConservative dU,
60c1a52365SJed Brown   const CeedScalar x[3], const CeedScalar dx[3]) {
61c1a52365SJed Brown   StatePrimitive dY;
62c1a52365SJed Brown   for (int i=0; i<3; i++) {
63c1a52365SJed Brown     dY.velocity[i] = (dU.momentum[i] - s.Y.velocity[i] * dU.density) / s.U.density;
64c1a52365SJed Brown   }
65c1a52365SJed Brown   CeedScalar e_kinetic = .5 * Dot3(s.Y.velocity, s.Y.velocity);
66c1a52365SJed Brown   CeedScalar de_kinetic = Dot3(dY.velocity, s.Y.velocity);
67c1a52365SJed Brown   CeedScalar e_potential = -Dot3(gas->g, x);
68c1a52365SJed Brown   CeedScalar de_potential = -Dot3(gas->g, dx);
69c1a52365SJed Brown   CeedScalar e_total = s.U.E_total / s.U.density;
70c1a52365SJed Brown   CeedScalar de_total = (dU.E_total - e_total * dU.density) / s.U.density;
71c1a52365SJed Brown   CeedScalar e_internal = e_total - e_kinetic - e_potential;
72c1a52365SJed Brown   CeedScalar de_internal = de_total - de_kinetic - de_potential;
73c1a52365SJed Brown   dY.temperature = de_internal / gas->cv;
74c1a52365SJed Brown   dY.pressure = (gas->cp / gas->cv - 1)
75c1a52365SJed Brown                 * (dU.density * e_internal + s.U.density * de_internal);
76c1a52365SJed Brown   return dY;
77c1a52365SJed Brown }
78c1a52365SJed Brown 
79c1a52365SJed Brown CEED_QFUNCTION_HELPER State StateFromU(NewtonianIdealGasContext gas,
80c1a52365SJed Brown                                        const CeedScalar U[5], const CeedScalar x[3]) {
81c1a52365SJed Brown   State s;
82c1a52365SJed Brown   s.U.density = U[0];
83c1a52365SJed Brown   s.U.momentum[0] = U[1];
84c1a52365SJed Brown   s.U.momentum[1] = U[2];
85c1a52365SJed Brown   s.U.momentum[2] = U[3];
86c1a52365SJed Brown   s.U.E_total = U[4];
87c1a52365SJed Brown   s.Y = StatePrimitiveFromConservative(gas, s.U, x);
88c1a52365SJed Brown   return s;
89c1a52365SJed Brown }
90c1a52365SJed Brown 
912f7ce6c1SJed Brown CEED_QFUNCTION_HELPER State StateFromU_fwd(NewtonianIdealGasContext gas,
922f7ce6c1SJed Brown     State s, const CeedScalar dU[5],
932f7ce6c1SJed Brown     const CeedScalar x[3], const CeedScalar dx[3]) {
942f7ce6c1SJed Brown   State ds;
952f7ce6c1SJed Brown   ds.U.density = dU[0];
962f7ce6c1SJed Brown   ds.U.momentum[0] = dU[1];
972f7ce6c1SJed Brown   ds.U.momentum[1] = dU[2];
982f7ce6c1SJed Brown   ds.U.momentum[2] = dU[3];
992f7ce6c1SJed Brown   ds.U.E_total = dU[4];
1002f7ce6c1SJed Brown   ds.Y = StatePrimitiveFromConservative_fwd(gas, s, ds.U, x, dx);
1012f7ce6c1SJed Brown   return ds;
1022f7ce6c1SJed Brown }
1032f7ce6c1SJed Brown 
104c1a52365SJed Brown CEED_QFUNCTION_HELPER void FluxInviscid(NewtonianIdealGasContext gas, State s,
105c1a52365SJed Brown                                         StateConservative Flux[3]) {
106c1a52365SJed Brown   for (int i=0; i<3; i++) {
107c1a52365SJed Brown     Flux[i].density = s.U.momentum[i];
108c1a52365SJed Brown     for (int j=0; j<3; j++)
109c1a52365SJed Brown       Flux[i].momentum[j] = s.U.momentum[i] * s.Y.velocity[j]
110c1a52365SJed Brown                             + s.Y.pressure * (i == j);
111c1a52365SJed Brown     Flux[i].E_total = (s.U.E_total + s.Y.pressure) * s.Y.velocity[i];
112c1a52365SJed Brown   }
113c1a52365SJed Brown }
114c1a52365SJed Brown 
115c1a52365SJed Brown CEED_QFUNCTION_HELPER void FluxInviscid_fwd(NewtonianIdealGasContext gas,
116c1a52365SJed Brown     State s, State ds, StateConservative dFlux[3]) {
117c1a52365SJed Brown   for (int i=0; i<3; i++) {
118c1a52365SJed Brown     dFlux[i].density = ds.U.momentum[i];
119c1a52365SJed Brown     for (int j=0; j<3; j++)
120c1a52365SJed Brown       dFlux[i].momentum[j] = ds.U.momentum[i] * s.Y.velocity[j] +
121c1a52365SJed Brown                              s.U.momentum[i] * ds.Y.velocity[j] + ds.Y.pressure * (i == j);
122c1a52365SJed Brown     dFlux[i].E_total = (ds.U.E_total + ds.Y.pressure) * s.Y.velocity[i] +
123c1a52365SJed Brown                        (s.U.E_total + s.Y.pressure) * ds.Y.velocity[i];
124c1a52365SJed Brown   }
125c1a52365SJed Brown }
126c1a52365SJed Brown 
127c1a52365SJed Brown // Kelvin-Mandel notation
128c1a52365SJed Brown CEED_QFUNCTION_HELPER void KMStrainRate(const State grad_s[3],
129c1a52365SJed Brown                                         CeedScalar strain_rate[6]) {
130c1a52365SJed Brown   const CeedScalar weight = 1 / sqrt(2.);
131c1a52365SJed Brown   strain_rate[0] = grad_s[0].Y.velocity[0];
132c1a52365SJed Brown   strain_rate[1] = grad_s[1].Y.velocity[1];
133c1a52365SJed Brown   strain_rate[2] = grad_s[2].Y.velocity[2];
134c1a52365SJed Brown   strain_rate[3] = weight * (grad_s[2].Y.velocity[1] + grad_s[1].Y.velocity[2]);
135c1a52365SJed Brown   strain_rate[4] = weight * (grad_s[2].Y.velocity[0] + grad_s[0].Y.velocity[2]);
136c1a52365SJed Brown   strain_rate[5] = weight * (grad_s[1].Y.velocity[0] + grad_s[0].Y.velocity[1]);
137c1a52365SJed Brown }
138c1a52365SJed Brown 
139c1a52365SJed Brown CEED_QFUNCTION_HELPER void KMUnpack(const CeedScalar v[6], CeedScalar A[3][3]) {
140c1a52365SJed Brown   const CeedScalar weight = 1 / sqrt(2.);
141c1a52365SJed Brown   A[0][0] = v[0];
142c1a52365SJed Brown   A[1][1] = v[1];
143c1a52365SJed Brown   A[2][2] = v[2];
144c1a52365SJed Brown   A[2][1] = A[1][2] = weight * v[3];
145c1a52365SJed Brown   A[2][0] = A[0][2] = weight * v[4];
146c1a52365SJed Brown   A[1][0] = A[0][1] = weight * v[5];
147c1a52365SJed Brown }
148c1a52365SJed Brown 
149c1a52365SJed Brown CEED_QFUNCTION_HELPER void NewtonianStress(NewtonianIdealGasContext gas,
150c1a52365SJed Brown     const CeedScalar strain_rate[6], CeedScalar stress[6]) {
151c1a52365SJed Brown   CeedScalar div_u = strain_rate[0] + strain_rate[1] + strain_rate[2];
152c1a52365SJed Brown   for (int i=0; i<6; i++) {
153c1a52365SJed Brown     stress[i] = gas->mu * (2 * strain_rate[i] + gas->lambda * div_u * (i < 3));
154c1a52365SJed Brown   }
155c1a52365SJed Brown }
156c1a52365SJed Brown 
157c1a52365SJed Brown CEED_QFUNCTION_HELPER void ViscousEnergyFlux(NewtonianIdealGasContext gas,
158c1a52365SJed Brown     StatePrimitive Y, const State grad_s[3], const CeedScalar stress[3][3],
159c1a52365SJed Brown     CeedScalar Fe[3]) {
160c1a52365SJed Brown   for (int i=0; i<3; i++) {
161c1a52365SJed Brown     Fe[i] = - Y.velocity[0] * stress[0][i]
162c1a52365SJed Brown             - Y.velocity[1] * stress[1][i]
163c1a52365SJed Brown             - Y.velocity[2] * stress[2][i]
164c1a52365SJed Brown             - gas->k * grad_s[i].Y.temperature;
165c1a52365SJed Brown   }
166c1a52365SJed Brown }
167c1a52365SJed Brown 
1683a8779fbSJames Wright // *****************************************************************************
1693a8779fbSJames Wright // Helper function for computing flux Jacobian
1703a8779fbSJames Wright // *****************************************************************************
1713a8779fbSJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NS(CeedScalar dF[3][5][5],
1723a8779fbSJames Wright     const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
173bb8a0c61SJames Wright     const CeedScalar gamma, const CeedScalar g[3], const CeedScalar x[3]) {
1743a8779fbSJames Wright   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square
175bb8a0c61SJames Wright   CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]);
1763a8779fbSJames Wright   for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions
1773a8779fbSJames Wright     for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix
178bb8a0c61SJames Wright       dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2. - e_potential)) : 0.) -
179bb8a0c61SJames Wright                       u[i]*u[j];
1803a8779fbSJames Wright       for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix
1813a8779fbSJames Wright         dF[i][0][k+1]   = ((i==k) ? 1. : 0.);
1823a8779fbSJames Wright         dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) +
1833a8779fbSJames Wright                           ((i==k) ? u[j] : 0.) -
1843a8779fbSJames Wright                           ((i==j) ? u[k] : 0.) * (gamma-1.);
1853a8779fbSJames Wright         dF[i][4][k+1]   = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) -
1863a8779fbSJames Wright                           (gamma-1.)*u[i]*u[k];
1873a8779fbSJames Wright       }
1883a8779fbSJames Wright       dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.);
1893a8779fbSJames Wright     }
1903a8779fbSJames Wright     dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho);
1913a8779fbSJames Wright     dF[i][4][4] = u[i] * gamma;
1923a8779fbSJames Wright   }
1933a8779fbSJames Wright }
1943a8779fbSJames Wright 
1953a8779fbSJames Wright // *****************************************************************************
196bb8a0c61SJames Wright // Helper function for computing flux Jacobian of Primitive variables
197bb8a0c61SJames Wright // *****************************************************************************
198bb8a0c61SJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NSp(CeedScalar dF[3][5][5],
199bb8a0c61SJames Wright     const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
200bb8a0c61SJames Wright     const CeedScalar Rd, const CeedScalar cv) {
201bb8a0c61SJames Wright   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square
202bb8a0c61SJames Wright   // TODO Add in gravity's contribution
203bb8a0c61SJames Wright 
204bb8a0c61SJames Wright   CeedScalar T    = ( E / rho - u_sq / 2. ) / cv;
205bb8a0c61SJames Wright   CeedScalar drdT = -rho / T;
206bb8a0c61SJames Wright   CeedScalar drdP = 1. / ( Rd * T);
207bb8a0c61SJames Wright   CeedScalar etot =  E / rho ;
208bb8a0c61SJames Wright   CeedScalar e2p  = drdP * etot + 1. ;
209bb8a0c61SJames Wright   CeedScalar e3p  = ( E  + rho * Rd * T );
210bb8a0c61SJames Wright   CeedScalar e4p  = drdT * etot + rho * cv ;
211bb8a0c61SJames Wright 
212bb8a0c61SJames Wright   for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions
213bb8a0c61SJames Wright     for (CeedInt j=0; j<3; j++) { // j counts F^{m_j}
214bb8a0c61SJames Wright //        [row][col] of A_i
215bb8a0c61SJames Wright       dF[i][j+1][0] = drdP * u[i] * u[j] + ((i==j) ? 1. : 0.); // F^{{m_j} wrt p
216bb8a0c61SJames Wright       for (CeedInt k=0; k<3; k++) { // k counts the wrt vel_k
2172acc7cbcSKenneth E. Jansen         dF[i][0][k+1]   =  ((i==k) ? rho  : 0.);   // F^c wrt u_k
218bb8a0c61SJames Wright         dF[i][j+1][k+1] = (((j==k) ? u[i] : 0.) +  // F^m_j wrt u_k
219bb8a0c61SJames Wright                            ((i==k) ? u[j] : 0.) ) * rho;
220bb8a0c61SJames Wright         dF[i][4][k+1]   = rho * u[i] * u[k]
221bb8a0c61SJames Wright                           + ((i==k) ? e3p  : 0.) ; // F^e wrt u_k
222bb8a0c61SJames Wright       }
223bb8a0c61SJames Wright       dF[i][j+1][4] = drdT * u[i] * u[j]; // F^{m_j} wrt T
224bb8a0c61SJames Wright     }
225bb8a0c61SJames Wright     dF[i][4][0] = u[i] * e2p; // F^e wrt p
226bb8a0c61SJames Wright     dF[i][4][4] = u[i] * e4p; // F^e wrt T
227bb8a0c61SJames Wright     dF[i][0][0] = u[i] * drdP; // F^c wrt p
228bb8a0c61SJames Wright     dF[i][0][4] = u[i] * drdT; // F^c wrt T
229bb8a0c61SJames Wright   }
230bb8a0c61SJames Wright }
231bb8a0c61SJames Wright 
232bb8a0c61SJames Wright CEED_QFUNCTION_HELPER void PrimitiveToConservative_fwd(const CeedScalar rho,
233bb8a0c61SJames Wright     const CeedScalar u[3], const CeedScalar E, const CeedScalar Rd,
234bb8a0c61SJames Wright     const CeedScalar cv, const CeedScalar dY[5], CeedScalar dU[5]) {
235bb8a0c61SJames Wright   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2];
236bb8a0c61SJames Wright   CeedScalar T    = ( E / rho - u_sq / 2. ) / cv;
237bb8a0c61SJames Wright   CeedScalar drdT = -rho / T;
238bb8a0c61SJames Wright   CeedScalar drdP = 1. / ( Rd * T);
239bb8a0c61SJames Wright   dU[0] = drdP * dY[0] + drdT * dY[4];
240bb8a0c61SJames Wright   CeedScalar de_kinetic = 0;
241bb8a0c61SJames Wright   for (int i=0; i<3; i++) {
242bb8a0c61SJames Wright     dU[1+i] = dU[0] * u[i] + rho * dY[1+i];
243bb8a0c61SJames Wright     de_kinetic += u[i] * dY[1+i];
244bb8a0c61SJames Wright   }
245bb8a0c61SJames Wright   dU[4] = rho * cv * dY[4] + dU[0] * cv * T // internal energy: rho * e
246bb8a0c61SJames Wright           + rho * de_kinetic + .5 * dU[0] * u_sq; // kinetic energy: .5 * rho * |u|^2
247bb8a0c61SJames Wright }
248bb8a0c61SJames Wright 
249bb8a0c61SJames Wright // *****************************************************************************
250bb8a0c61SJames Wright // Helper function for computing Tau elements (stabilization constant)
251bb8a0c61SJames Wright //   Model from:
252bb8a0c61SJames Wright //     PHASTA
253bb8a0c61SJames Wright //
254bb8a0c61SJames Wright //   Tau[i] = itau=0 which is diagonal-Shakib (3 values still but not spatial)
255bb8a0c61SJames Wright //
256bb8a0c61SJames Wright // Where NOT UPDATED YET
257bb8a0c61SJames Wright // *****************************************************************************
258bb8a0c61SJames Wright CEED_QFUNCTION_HELPER void Tau_diagPrim(CeedScalar Tau_d[3],
259bb8a0c61SJames Wright                                         const CeedScalar dXdx[3][3], const CeedScalar u[3],
260bb8a0c61SJames Wright                                         const CeedScalar cv, const NewtonianIdealGasContext newt_ctx,
261bb8a0c61SJames Wright                                         const CeedScalar mu, const CeedScalar dt,
262bb8a0c61SJames Wright                                         const CeedScalar rho) {
263bb8a0c61SJames Wright   // Context
264bb8a0c61SJames Wright   const CeedScalar Ctau_t = newt_ctx->Ctau_t;
265bb8a0c61SJames Wright   const CeedScalar Ctau_v = newt_ctx->Ctau_v;
266bb8a0c61SJames Wright   const CeedScalar Ctau_C = newt_ctx->Ctau_C;
267bb8a0c61SJames Wright   const CeedScalar Ctau_M = newt_ctx->Ctau_M;
268bb8a0c61SJames Wright   const CeedScalar Ctau_E = newt_ctx->Ctau_E;
269bb8a0c61SJames Wright   CeedScalar gijd[6];
270bb8a0c61SJames Wright   CeedScalar tau;
271bb8a0c61SJames Wright   CeedScalar dts;
272bb8a0c61SJames Wright   CeedScalar fact;
273bb8a0c61SJames Wright 
274bb8a0c61SJames Wright   //*INDENT-OFF*
275bb8a0c61SJames Wright   gijd[0] =   dXdx[0][0] * dXdx[0][0]
276bb8a0c61SJames Wright             + dXdx[1][0] * dXdx[1][0]
277bb8a0c61SJames Wright             + dXdx[2][0] * dXdx[2][0];
278bb8a0c61SJames Wright 
279bb8a0c61SJames Wright   gijd[1] =   dXdx[0][0] * dXdx[0][1]
280bb8a0c61SJames Wright             + dXdx[1][0] * dXdx[1][1]
281bb8a0c61SJames Wright             + dXdx[2][0] * dXdx[2][1];
282bb8a0c61SJames Wright 
283bb8a0c61SJames Wright   gijd[2] =   dXdx[0][1] * dXdx[0][1]
284bb8a0c61SJames Wright             + dXdx[1][1] * dXdx[1][1]
285bb8a0c61SJames Wright             + dXdx[2][1] * dXdx[2][1];
286bb8a0c61SJames Wright 
287bb8a0c61SJames Wright   gijd[3] =   dXdx[0][0] * dXdx[0][2]
288bb8a0c61SJames Wright             + dXdx[1][0] * dXdx[1][2]
289bb8a0c61SJames Wright             + dXdx[2][0] * dXdx[2][2];
290bb8a0c61SJames Wright 
291bb8a0c61SJames Wright   gijd[4] =   dXdx[0][1] * dXdx[0][2]
292bb8a0c61SJames Wright             + dXdx[1][1] * dXdx[1][2]
293bb8a0c61SJames Wright             + dXdx[2][1] * dXdx[2][2];
294bb8a0c61SJames Wright 
295bb8a0c61SJames Wright   gijd[5] =   dXdx[0][2] * dXdx[0][2]
296bb8a0c61SJames Wright             + dXdx[1][2] * dXdx[1][2]
297bb8a0c61SJames Wright             + dXdx[2][2] * dXdx[2][2];
298bb8a0c61SJames Wright   //*INDENT-ON*
299bb8a0c61SJames Wright 
300bb8a0c61SJames Wright   dts = Ctau_t / dt ;
301bb8a0c61SJames Wright 
302bb8a0c61SJames Wright   tau = rho*rho*((4. * dts * dts)
303bb8a0c61SJames Wright                  + u[0] * ( u[0] * gijd[0] + 2. * ( u[1] * gijd[1] + u[2] * gijd[3]))
304bb8a0c61SJames Wright                  + u[1] * ( u[1] * gijd[2] + 2. *   u[2] * gijd[4])
305bb8a0c61SJames Wright                  + u[2] *   u[2] * gijd[5])
306bb8a0c61SJames Wright         + Ctau_v* mu * mu *
307bb8a0c61SJames Wright         (gijd[0]*gijd[0] + gijd[2]*gijd[2] + gijd[5]*gijd[5] +
308bb8a0c61SJames Wright          + 2. * (gijd[1]*gijd[1] + gijd[3]*gijd[3] + gijd[4]*gijd[4]));
309bb8a0c61SJames Wright 
310bb8a0c61SJames Wright   fact=sqrt(tau);
311bb8a0c61SJames Wright 
312bb8a0c61SJames Wright   Tau_d[0] = Ctau_C * fact / (rho*(gijd[0] + gijd[2] + gijd[5]))*0.125;
313bb8a0c61SJames Wright 
314bb8a0c61SJames Wright   Tau_d[1] = Ctau_M / fact;
315bb8a0c61SJames Wright   Tau_d[2] = Ctau_E / ( fact * cv );
316bb8a0c61SJames Wright 
317bb8a0c61SJames Wright // consider putting back the way I initially had it  Ctau_E * Tau_d[1] /cv
318bb8a0c61SJames Wright //  to avoid a division if the compiler is smart enough to see that cv IS
319bb8a0c61SJames Wright // a constant that it could invert once for all elements
320bb8a0c61SJames Wright // but in that case energy tau is scaled by the product of Ctau_E * Ctau_M
321bb8a0c61SJames Wright // OR we could absorb cv into Ctau_E but this puts more burden on user to
322bb8a0c61SJames Wright // know how to change constants with a change of fluid or units.  Same for
323bb8a0c61SJames Wright // Ctau_v * mu * mu IF AND ONLY IF we don't add viscosity law =f(T)
324bb8a0c61SJames Wright }
325bb8a0c61SJames Wright 
326bb8a0c61SJames Wright // *****************************************************************************
3273a8779fbSJames Wright // This QFunction sets a "still" initial condition for generic Newtonian IG problems
3283a8779fbSJames Wright // *****************************************************************************
3293a8779fbSJames Wright CEED_QFUNCTION(ICsNewtonianIG)(void *ctx, CeedInt Q,
3303a8779fbSJames Wright                                const CeedScalar *const *in, CeedScalar *const *out) {
3313a8779fbSJames Wright   // Inputs
3323a8779fbSJames Wright   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
3333a8779fbSJames Wright 
3343a8779fbSJames Wright   // Outputs
3353a8779fbSJames Wright   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
3363a8779fbSJames Wright 
337bb8a0c61SJames Wright   // Context
338bb8a0c61SJames Wright   const SetupContext context = (SetupContext)ctx;
339bb8a0c61SJames Wright   const CeedScalar theta0    = context->theta0;
340bb8a0c61SJames Wright   const CeedScalar P0        = context->P0;
341bb8a0c61SJames Wright   const CeedScalar cv        = context->cv;
342bb8a0c61SJames Wright   const CeedScalar cp        = context->cp;
343bb8a0c61SJames Wright   const CeedScalar *g        = context->g;
344bb8a0c61SJames Wright   const CeedScalar Rd        = cp - cv;
345bb8a0c61SJames Wright 
3463a8779fbSJames Wright   // Quadrature Point Loop
3473a8779fbSJames Wright   CeedPragmaSIMD
3483a8779fbSJames Wright   for (CeedInt i=0; i<Q; i++) {
3493a8779fbSJames Wright     CeedScalar q[5] = {0.};
3503a8779fbSJames Wright 
3513a8779fbSJames Wright     // Setup
3523a8779fbSJames Wright     // -- Coordinates
353bb8a0c61SJames Wright     const CeedScalar x[3] = {X[0][i], X[1][i], X[2][i]};
354bb8a0c61SJames Wright     const CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]);
3553a8779fbSJames Wright 
3563a8779fbSJames Wright     // -- Density
357bb8a0c61SJames Wright     const CeedScalar rho = P0 / (Rd*theta0);
3583a8779fbSJames Wright 
3593a8779fbSJames Wright     // Initial Conditions
3603a8779fbSJames Wright     q[0] = rho;
3613a8779fbSJames Wright     q[1] = 0.0;
3623a8779fbSJames Wright     q[2] = 0.0;
3633a8779fbSJames Wright     q[3] = 0.0;
364bb8a0c61SJames Wright     q[4] = rho * (cv*theta0 + e_potential);
3653a8779fbSJames Wright 
3663a8779fbSJames Wright     for (CeedInt j=0; j<5; j++)
3673a8779fbSJames Wright       q0[j][i] = q[j];
3683a8779fbSJames Wright   } // End of Quadrature Point Loop
3693a8779fbSJames Wright   return 0;
3703a8779fbSJames Wright }
3713a8779fbSJames Wright 
3723a8779fbSJames Wright // *****************************************************************************
3733a8779fbSJames Wright // This QFunction implements the following formulation of Navier-Stokes with
3743a8779fbSJames Wright //   explicit time stepping method
3753a8779fbSJames Wright //
3763a8779fbSJames Wright // This is 3D compressible Navier-Stokes in conservation form with state
3773a8779fbSJames Wright //   variables of density, momentum density, and total energy density.
3783a8779fbSJames Wright //
3793a8779fbSJames Wright // State Variables: q = ( rho, U1, U2, U3, E )
3803a8779fbSJames Wright //   rho - Mass Density
3813a8779fbSJames Wright //   Ui  - Momentum Density,      Ui = rho ui
3823a8779fbSJames Wright //   E   - Total Energy Density,  E  = rho (cv T + (u u)/2 + g z)
3833a8779fbSJames Wright //
3843a8779fbSJames Wright // Navier-Stokes Equations:
3853a8779fbSJames Wright //   drho/dt + div( U )                               = 0
3863a8779fbSJames Wright //   dU/dt   + div( rho (u x u) + P I3 ) + rho g khat = div( Fu )
3873a8779fbSJames Wright //   dE/dt   + div( (E + P) u )                       = div( Fe )
3883a8779fbSJames Wright //
3893a8779fbSJames Wright // Viscous Stress:
3903a8779fbSJames Wright //   Fu = mu (grad( u ) + grad( u )^T + lambda div ( u ) I3)
3913a8779fbSJames Wright //
3923a8779fbSJames Wright // Thermal Stress:
3933a8779fbSJames Wright //   Fe = u Fu + k grad( T )
394bb8a0c61SJames Wright // Equation of State
3953a8779fbSJames Wright //   P = (gamma - 1) (E - rho (u u) / 2 - rho g z)
3963a8779fbSJames Wright //
3973a8779fbSJames Wright // Stabilization:
3983a8779fbSJames Wright //   Tau = diag(TauC, TauM, TauM, TauM, TauE)
3993a8779fbSJames Wright //     f1 = rho  sqrt(ui uj gij)
4003a8779fbSJames Wright //     gij = dXi/dX * dXi/dX
4013a8779fbSJames Wright //     TauC = Cc f1 / (8 gii)
4023a8779fbSJames Wright //     TauM = min( 1 , 1 / f1 )
4033a8779fbSJames Wright //     TauE = TauM / (Ce cv)
4043a8779fbSJames Wright //
4053a8779fbSJames Wright //  SU   = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) )
4063a8779fbSJames Wright //
4073a8779fbSJames Wright // Constants:
4083a8779fbSJames Wright //   lambda = - 2 / 3,  From Stokes hypothesis
4093a8779fbSJames Wright //   mu              ,  Dynamic viscosity
4103a8779fbSJames Wright //   k               ,  Thermal conductivity
4113a8779fbSJames Wright //   cv              ,  Specific heat, constant volume
4123a8779fbSJames Wright //   cp              ,  Specific heat, constant pressure
4133a8779fbSJames Wright //   g               ,  Gravity
4143a8779fbSJames Wright //   gamma  = cp / cv,  Specific heat ratio
4153a8779fbSJames Wright //
4163a8779fbSJames Wright // We require the product of the inverse of the Jacobian (dXdx_j,k) and
4173a8779fbSJames Wright // its transpose (dXdx_k,j) to properly compute integrals of the form:
4183a8779fbSJames Wright // int( gradv gradu )
4193a8779fbSJames Wright //
4203a8779fbSJames Wright // *****************************************************************************
421c1a52365SJed Brown CEED_QFUNCTION(RHSFunction_Newtonian)(void *ctx, CeedInt Q,
4223a8779fbSJames Wright                                       const CeedScalar *const *in, CeedScalar *const *out) {
4233a8779fbSJames Wright   // *INDENT-OFF*
4243a8779fbSJames Wright   // Inputs
4253a8779fbSJames Wright   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
426*752f40e3SJed Brown                    (*Grad_q)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
4273a8779fbSJames Wright                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
4283a8779fbSJames Wright                    (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
4293a8779fbSJames Wright   // Outputs
4303a8779fbSJames Wright   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
431*752f40e3SJed Brown              (*Grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
4323a8779fbSJames Wright   // *INDENT-ON*
4333a8779fbSJames Wright 
4343a8779fbSJames Wright   // Context
4353a8779fbSJames Wright   NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx;
4363a8779fbSJames Wright   const CeedScalar mu     = context->mu;
4373a8779fbSJames Wright   const CeedScalar cv     = context->cv;
4383a8779fbSJames Wright   const CeedScalar cp     = context->cp;
439bb8a0c61SJames Wright   const CeedScalar *g     = context->g;
440bb8a0c61SJames Wright   const CeedScalar dt     = context->dt;
4413a8779fbSJames Wright   const CeedScalar gamma  = cp / cv;
442bb8a0c61SJames Wright   const CeedScalar Rd     = cp - cv;
4433a8779fbSJames Wright 
4443a8779fbSJames Wright   CeedPragmaSIMD
4453a8779fbSJames Wright   // Quadrature Point Loop
4463a8779fbSJames Wright   for (CeedInt i=0; i<Q; i++) {
447c1a52365SJed Brown     CeedScalar U[5];
448c1a52365SJed Brown     for (int j=0; j<5; j++) U[j] = q[j][i];
449c1a52365SJed Brown     const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]};
450c1a52365SJed Brown     State s = StateFromU(context, U, x_i);
451c1a52365SJed Brown 
4523a8779fbSJames Wright     // -- Interp-to-Interp q_data
4533a8779fbSJames Wright     const CeedScalar wdetJ      =   q_data[0][i];
4543a8779fbSJames Wright     // -- Interp-to-Grad q_data
4553a8779fbSJames Wright     // ---- Inverse of change of coordinate matrix: X_i,j
4563a8779fbSJames Wright     // *INDENT-OFF*
4573a8779fbSJames Wright     const CeedScalar dXdx[3][3] = {{q_data[1][i],
4583a8779fbSJames Wright                                     q_data[2][i],
4593a8779fbSJames Wright                                     q_data[3][i]},
4603a8779fbSJames Wright                                    {q_data[4][i],
4613a8779fbSJames Wright                                     q_data[5][i],
4623a8779fbSJames Wright                                     q_data[6][i]},
4633a8779fbSJames Wright                                    {q_data[7][i],
4643a8779fbSJames Wright                                     q_data[8][i],
4653a8779fbSJames Wright                                     q_data[9][i]}
4663a8779fbSJames Wright                                   };
4673a8779fbSJames Wright     // *INDENT-ON*
4683a8779fbSJames Wright 
469c1a52365SJed Brown     State grad_s[3];
470c1a52365SJed Brown     for (int j=0; j<3; j++) {
4712f7ce6c1SJed Brown       CeedScalar dx_i[3] = {0}, dU[5];
472*752f40e3SJed Brown       for (int k=0; k<5; k++) dU[k] = Grad_q[0][k][i] * dXdx[0][j]
473*752f40e3SJed Brown                                         + Grad_q[1][k][i] * dXdx[1][j]
474*752f40e3SJed Brown                                         + Grad_q[2][k][i] * dXdx[2][j];
475c1a52365SJed Brown       dx_i[j] = 1.;
4762f7ce6c1SJed Brown       grad_s[j] = StateFromU_fwd(context, s, dU, x_i, dx_i);
477c1a52365SJed Brown     }
478c1a52365SJed Brown 
479c1a52365SJed Brown     CeedScalar strain_rate[6], kmstress[6], stress[3][3], Fe[3];
480c1a52365SJed Brown     KMStrainRate(grad_s, strain_rate);
481c1a52365SJed Brown     NewtonianStress(context, strain_rate, kmstress);
482c1a52365SJed Brown     KMUnpack(kmstress, stress);
483c1a52365SJed Brown     ViscousEnergyFlux(context, s.Y, grad_s, stress, Fe);
484c1a52365SJed Brown 
485c1a52365SJed Brown     StateConservative F_inviscid[3];
486c1a52365SJed Brown     FluxInviscid(context, s, F_inviscid);
487c1a52365SJed Brown 
488c1a52365SJed Brown     // Total flux
489c1a52365SJed Brown     CeedScalar Flux[5][3];
490c1a52365SJed Brown     for (int j=0; j<3; j++) {
491c1a52365SJed Brown       Flux[0][j] = F_inviscid[j].density;
492c1a52365SJed Brown       for (int k=0; k<3; k++)
493c1a52365SJed Brown         Flux[k+1][j] = F_inviscid[j].momentum[k] - stress[k][j];
494c1a52365SJed Brown       Flux[4][j] = F_inviscid[j].E_total + Fe[j];
495c1a52365SJed Brown     }
496c1a52365SJed Brown 
497c1a52365SJed Brown     for (int j=0; j<3; j++) {
498c1a52365SJed Brown       for (int k=0; k<5; k++) {
499*752f40e3SJed Brown         Grad_v[j][k][i] = wdetJ * (dXdx[j][0] * Flux[k][0] +
500c1a52365SJed Brown                                    dXdx[j][1] * Flux[k][1] +
501c1a52365SJed Brown                                    dXdx[j][2] * Flux[k][2]);
502c1a52365SJed Brown       }
503c1a52365SJed Brown     }
504c1a52365SJed Brown 
505c1a52365SJed Brown     const CeedScalar body_force[5] = {0, s.U.density *g[0], s.U.density *g[1], s.U.density *g[2], 0};
506c1a52365SJed Brown     for (int j=0; j<5; j++)
507c1a52365SJed Brown       v[j][i] = wdetJ * body_force[j];
5083a8779fbSJames Wright 
5093a8779fbSJames Wright     // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
510c1a52365SJed Brown     CeedScalar jacob_F_conv[3][5][5] = {0};
511c1a52365SJed Brown     computeFluxJacobian_NS(jacob_F_conv, s.U.density, s.Y.velocity, s.U.E_total,
512c1a52365SJed Brown                            gamma, g, x_i);
513c1a52365SJed Brown     CeedScalar grad_U[5][3];
5143a8779fbSJames Wright     for (int j=0; j<3; j++) {
515c1a52365SJed Brown       grad_U[0][j] = grad_s[j].U.density;
516c1a52365SJed Brown       for (int k=0; k<3; k++) grad_U[k+1][j] = grad_s[j].U.momentum[k];
517c1a52365SJed Brown       grad_U[4][j] = grad_s[j].U.E_total;
5183a8779fbSJames Wright     }
5193a8779fbSJames Wright 
5203a8779fbSJames Wright     // strong_conv = dF/dq * dq/dx    (Strong convection)
5213a8779fbSJames Wright     CeedScalar strong_conv[5] = {0};
5223a8779fbSJames Wright     for (int j=0; j<3; j++)
5233a8779fbSJames Wright       for (int k=0; k<5; k++)
5243a8779fbSJames Wright         for (int l=0; l<5; l++)
525c1a52365SJed Brown           strong_conv[k] += jacob_F_conv[j][k][l] * grad_U[l][j];
5263a8779fbSJames Wright 
527bb8a0c61SJames Wright     // -- Stabilization method: none, SU, or SUPG
528bb8a0c61SJames Wright     CeedScalar stab[5][3] = {{0.}};
529bb8a0c61SJames Wright     CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0};
530bb8a0c61SJames Wright     CeedScalar Tau_d[3] = {0.};
5313a8779fbSJames Wright     switch (context->stabilization) {
5323a8779fbSJames Wright     case STAB_NONE:        // Galerkin
5333a8779fbSJames Wright       break;
5343a8779fbSJames Wright     case STAB_SU:        // SU
535c1a52365SJed Brown       Tau_diagPrim(Tau_d, dXdx, s.Y.velocity, cv, context, mu, dt, s.U.density);
536bb8a0c61SJames Wright       tau_strong_conv[0] = Tau_d[0] * strong_conv[0];
537bb8a0c61SJames Wright       tau_strong_conv[1] = Tau_d[1] * strong_conv[1];
538bb8a0c61SJames Wright       tau_strong_conv[2] = Tau_d[1] * strong_conv[2];
539bb8a0c61SJames Wright       tau_strong_conv[3] = Tau_d[1] * strong_conv[3];
540bb8a0c61SJames Wright       tau_strong_conv[4] = Tau_d[2] * strong_conv[4];
541c1a52365SJed Brown       PrimitiveToConservative_fwd(s.U.density, s.Y.velocity, s.U.E_total, Rd, cv,
542c1a52365SJed Brown                                   tau_strong_conv,
543bb8a0c61SJames Wright                                   tau_strong_conv_conservative);
5443a8779fbSJames Wright       for (int j=0; j<3; j++)
5453a8779fbSJames Wright         for (int k=0; k<5; k++)
5463a8779fbSJames Wright           for (int l=0; l<5; l++)
547bb8a0c61SJames Wright             stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l];
5483a8779fbSJames Wright 
5493a8779fbSJames Wright       for (int j=0; j<5; j++)
5503a8779fbSJames Wright         for (int k=0; k<3; k++)
551*752f40e3SJed Brown           Grad_v[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] +
5523a8779fbSJames Wright                                     stab[j][1] * dXdx[k][1] +
5533a8779fbSJames Wright                                     stab[j][2] * dXdx[k][2]);
5543a8779fbSJames Wright       break;
5553a8779fbSJames Wright     case STAB_SUPG:        // SUPG is not implemented for explicit scheme
5563a8779fbSJames Wright       break;
5573a8779fbSJames Wright     }
5583a8779fbSJames Wright 
5593a8779fbSJames Wright   } // End Quadrature Point Loop
5603a8779fbSJames Wright 
5613a8779fbSJames Wright   // Return
5623a8779fbSJames Wright   return 0;
5633a8779fbSJames Wright }
5643a8779fbSJames Wright 
5653a8779fbSJames Wright // *****************************************************************************
5663a8779fbSJames Wright // This QFunction implements the Navier-Stokes equations (mentioned above) with
5673a8779fbSJames Wright //   implicit time stepping method
5683a8779fbSJames Wright //
5693a8779fbSJames Wright //  SU   = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) )
5703a8779fbSJames Wright //  SUPG = Galerkin + grad(v) . ( Ai^T * Tau * (q_dot + Aj q,j - body force) )
5713a8779fbSJames Wright //                                       (diffussive terms will be added later)
5723a8779fbSJames Wright //
5733a8779fbSJames Wright // *****************************************************************************
5743a8779fbSJames Wright CEED_QFUNCTION(IFunction_Newtonian)(void *ctx, CeedInt Q,
5753a8779fbSJames Wright                                     const CeedScalar *const *in,
5763a8779fbSJames Wright                                     CeedScalar *const *out) {
5773a8779fbSJames Wright   // *INDENT-OFF*
5783a8779fbSJames Wright   // Inputs
5793a8779fbSJames Wright   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
580*752f40e3SJed Brown                    (*Grad_q)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
5813a8779fbSJames Wright                    (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
5823a8779fbSJames Wright                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3],
5833a8779fbSJames Wright                    (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4];
5843a8779fbSJames Wright   // Outputs
5853a8779fbSJames Wright   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
586*752f40e3SJed Brown              (*Grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1],
587*752f40e3SJed Brown              (*jac_data)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[2];
5883a8779fbSJames Wright   // *INDENT-ON*
5893a8779fbSJames Wright   // Context
5903a8779fbSJames Wright   NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx;
5913a8779fbSJames Wright   const CeedScalar mu     = context->mu;
5923a8779fbSJames Wright   const CeedScalar cv     = context->cv;
5933a8779fbSJames Wright   const CeedScalar cp     = context->cp;
594bb8a0c61SJames Wright   const CeedScalar *g     = context->g;
595bb8a0c61SJames Wright   const CeedScalar dt     = context->dt;
5963a8779fbSJames Wright   const CeedScalar gamma  = cp / cv;
597bb8a0c61SJames Wright   const CeedScalar Rd     = cp-cv;
5983a8779fbSJames Wright 
5993a8779fbSJames Wright   CeedPragmaSIMD
6003a8779fbSJames Wright   // Quadrature Point Loop
6013a8779fbSJames Wright   for (CeedInt i=0; i<Q; i++) {
602c1a52365SJed Brown     CeedScalar U[5];
603c1a52365SJed Brown     for (int j=0; j<5; j++) U[j] = q[j][i];
604c1a52365SJed Brown     const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]};
605c1a52365SJed Brown     State s = StateFromU(context, U, x_i);
606c1a52365SJed Brown 
6073a8779fbSJames Wright     // -- Interp-to-Interp q_data
6083a8779fbSJames Wright     const CeedScalar wdetJ      =   q_data[0][i];
6093a8779fbSJames Wright     // -- Interp-to-Grad q_data
6103a8779fbSJames Wright     // ---- Inverse of change of coordinate matrix: X_i,j
6113a8779fbSJames Wright     // *INDENT-OFF*
6123a8779fbSJames Wright     const CeedScalar dXdx[3][3] = {{q_data[1][i],
6133a8779fbSJames Wright                                     q_data[2][i],
6143a8779fbSJames Wright                                     q_data[3][i]},
6153a8779fbSJames Wright                                    {q_data[4][i],
6163a8779fbSJames Wright                                     q_data[5][i],
6173a8779fbSJames Wright                                     q_data[6][i]},
6183a8779fbSJames Wright                                    {q_data[7][i],
6193a8779fbSJames Wright                                     q_data[8][i],
6203a8779fbSJames Wright                                     q_data[9][i]}
6213a8779fbSJames Wright                                   };
6223a8779fbSJames Wright     // *INDENT-ON*
623c1a52365SJed Brown     State grad_s[3];
6243a8779fbSJames Wright     for (int j=0; j<3; j++) {
6252f7ce6c1SJed Brown       CeedScalar dx_i[3] = {0}, dU[5];
626*752f40e3SJed Brown       for (int k=0; k<5; k++) dU[k] = Grad_q[0][k][i] * dXdx[0][j]
627*752f40e3SJed Brown                                         + Grad_q[1][k][i] * dXdx[1][j]
628*752f40e3SJed Brown                                         + Grad_q[2][k][i] * dXdx[2][j];
629c1a52365SJed Brown       dx_i[j] = 1.;
6302f7ce6c1SJed Brown       grad_s[j] = StateFromU_fwd(context, s, dU, x_i, dx_i);
6313a8779fbSJames Wright     }
632c1a52365SJed Brown 
633c1a52365SJed Brown     CeedScalar strain_rate[6], kmstress[6], stress[3][3], Fe[3];
634c1a52365SJed Brown     KMStrainRate(grad_s, strain_rate);
635c1a52365SJed Brown     NewtonianStress(context, strain_rate, kmstress);
636c1a52365SJed Brown     KMUnpack(kmstress, stress);
637c1a52365SJed Brown     ViscousEnergyFlux(context, s.Y, grad_s, stress, Fe);
638c1a52365SJed Brown 
639c1a52365SJed Brown     StateConservative F_inviscid[3];
640c1a52365SJed Brown     FluxInviscid(context, s, F_inviscid);
641c1a52365SJed Brown 
642c1a52365SJed Brown 
643c1a52365SJed Brown     // Total flux
644c1a52365SJed Brown     CeedScalar Flux[5][3];
645c1a52365SJed Brown     for (int j=0; j<3; j++) {
646c1a52365SJed Brown       Flux[0][j] = F_inviscid[j].density;
6473a8779fbSJames Wright       for (int k=0; k<3; k++)
648c1a52365SJed Brown         Flux[k+1][j] = F_inviscid[j].momentum[k] - stress[k][j];
649c1a52365SJed Brown       Flux[4][j] = F_inviscid[j].E_total + Fe[j];
650c1a52365SJed Brown     }
651c1a52365SJed Brown 
652c1a52365SJed Brown     for (int j=0; j<3; j++) {
653c1a52365SJed Brown       for (int k=0; k<5; k++) {
654*752f40e3SJed Brown         Grad_v[j][k][i] = -wdetJ * (dXdx[j][0] * Flux[k][0] +
655c1a52365SJed Brown                                     dXdx[j][1] * Flux[k][1] +
656c1a52365SJed Brown                                     dXdx[j][2] * Flux[k][2]);
657c1a52365SJed Brown       }
658c1a52365SJed Brown     }
659c1a52365SJed Brown 
660c1a52365SJed Brown     const CeedScalar body_force[5] = {0, s.U.density *g[0], s.U.density *g[1], s.U.density *g[2], 0};
661c1a52365SJed Brown     for (int j=0; j<5; j++)
662c1a52365SJed Brown       v[j][i] = wdetJ * (q_dot[j][i] - body_force[j]);
6633a8779fbSJames Wright 
6643a8779fbSJames Wright     // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
665c1a52365SJed Brown     CeedScalar jacob_F_conv[3][5][5] = {0};
666c1a52365SJed Brown     computeFluxJacobian_NS(jacob_F_conv, s.U.density, s.Y.velocity, s.U.E_total,
667c1a52365SJed Brown                            gamma, g, x_i);
668c1a52365SJed Brown     CeedScalar grad_U[5][3];
6693a8779fbSJames Wright     for (int j=0; j<3; j++) {
670c1a52365SJed Brown       grad_U[0][j] = grad_s[j].U.density;
671c1a52365SJed Brown       for (int k=0; k<3; k++) grad_U[k+1][j] = grad_s[j].U.momentum[k];
672c1a52365SJed Brown       grad_U[4][j] = grad_s[j].U.E_total;
6733a8779fbSJames Wright     }
674c1a52365SJed Brown 
6753a8779fbSJames Wright     // strong_conv = dF/dq * dq/dx    (Strong convection)
6763a8779fbSJames Wright     CeedScalar strong_conv[5] = {0};
6773a8779fbSJames Wright     for (int j=0; j<3; j++)
6783a8779fbSJames Wright       for (int k=0; k<5; k++)
6793a8779fbSJames Wright         for (int l=0; l<5; l++)
680c1a52365SJed Brown           strong_conv[k] += jacob_F_conv[j][k][l] * grad_U[l][j];
6813a8779fbSJames Wright 
6823a8779fbSJames Wright     // Strong residual
6833a8779fbSJames Wright     CeedScalar strong_res[5];
6843a8779fbSJames Wright     for (int j=0; j<5; j++)
6853a8779fbSJames Wright       strong_res[j] = q_dot[j][i] + strong_conv[j] - body_force[j];
6863a8779fbSJames Wright 
6873a8779fbSJames Wright     // -- Stabilization method: none, SU, or SUPG
688bb8a0c61SJames Wright     CeedScalar stab[5][3] = {{0.}};
689bb8a0c61SJames Wright     CeedScalar tau_strong_res[5] = {0.}, tau_strong_res_conservative[5] = {0};
690bb8a0c61SJames Wright     CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0};
691bb8a0c61SJames Wright     CeedScalar Tau_d[3] = {0.};
6923a8779fbSJames Wright     switch (context->stabilization) {
6933a8779fbSJames Wright     case STAB_NONE:        // Galerkin
6943a8779fbSJames Wright       break;
6953a8779fbSJames Wright     case STAB_SU:        // SU
696c1a52365SJed Brown       Tau_diagPrim(Tau_d, dXdx, s.Y.velocity, cv, context, mu, dt, s.U.density);
697bb8a0c61SJames Wright       tau_strong_conv[0] = Tau_d[0] * strong_conv[0];
698bb8a0c61SJames Wright       tau_strong_conv[1] = Tau_d[1] * strong_conv[1];
699bb8a0c61SJames Wright       tau_strong_conv[2] = Tau_d[1] * strong_conv[2];
700bb8a0c61SJames Wright       tau_strong_conv[3] = Tau_d[1] * strong_conv[3];
701bb8a0c61SJames Wright       tau_strong_conv[4] = Tau_d[2] * strong_conv[4];
702c1a52365SJed Brown       PrimitiveToConservative_fwd(s.U.density, s.Y.velocity, s.U.E_total, Rd, cv,
703c1a52365SJed Brown                                   tau_strong_conv, tau_strong_conv_conservative);
7043a8779fbSJames Wright       for (int j=0; j<3; j++)
7053a8779fbSJames Wright         for (int k=0; k<5; k++)
7063a8779fbSJames Wright           for (int l=0; l<5; l++)
707bb8a0c61SJames Wright             stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l];
7083a8779fbSJames Wright 
7093a8779fbSJames Wright       for (int j=0; j<5; j++)
7103a8779fbSJames Wright         for (int k=0; k<3; k++)
711*752f40e3SJed Brown           Grad_v[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
7123a8779fbSJames Wright                                     stab[j][1] * dXdx[k][1] +
7133a8779fbSJames Wright                                     stab[j][2] * dXdx[k][2]);
7143a8779fbSJames Wright       break;
7153a8779fbSJames Wright     case STAB_SUPG:        // SUPG
716c1a52365SJed Brown       Tau_diagPrim(Tau_d, dXdx, s.Y.velocity, cv, context, mu, dt, s.U.density);
717bb8a0c61SJames Wright       tau_strong_res[0] = Tau_d[0] * strong_res[0];
718bb8a0c61SJames Wright       tau_strong_res[1] = Tau_d[1] * strong_res[1];
719bb8a0c61SJames Wright       tau_strong_res[2] = Tau_d[1] * strong_res[2];
720bb8a0c61SJames Wright       tau_strong_res[3] = Tau_d[1] * strong_res[3];
721bb8a0c61SJames Wright       tau_strong_res[4] = Tau_d[2] * strong_res[4];
722bb8a0c61SJames Wright // Alternate route (useful later with primitive variable code)
723bb8a0c61SJames Wright // this function was verified against PHASTA for as IC that was as close as possible
724bb8a0c61SJames Wright //    computeFluxJacobian_NSp(jacob_F_conv_p, rho, u, E, Rd, cv);
725bb8a0c61SJames Wright // it has also been verified to compute a correct through the following
726bb8a0c61SJames Wright //   stab[k][j] += jacob_F_conv_p[j][k][l] * tau_strong_res[l] // flux Jacobian wrt primitive
727bb8a0c61SJames Wright // applied in the triple loop below
728bb8a0c61SJames Wright //  However, it is more flops than using the existing Jacobian wrt q after q_{,Y} viz
729c1a52365SJed Brown       PrimitiveToConservative_fwd(s.U.density, s.Y.velocity, s.U.E_total, Rd, cv,
730c1a52365SJed Brown                                   tau_strong_res, tau_strong_res_conservative);
7313a8779fbSJames Wright       for (int j=0; j<3; j++)
7323a8779fbSJames Wright         for (int k=0; k<5; k++)
7333a8779fbSJames Wright           for (int l=0; l<5; l++)
734bb8a0c61SJames Wright             stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_res_conservative[l];
7353a8779fbSJames Wright 
7363a8779fbSJames Wright       for (int j=0; j<5; j++)
7373a8779fbSJames Wright         for (int k=0; k<3; k++)
738*752f40e3SJed Brown           Grad_v[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
7393a8779fbSJames Wright                                     stab[j][1] * dXdx[k][1] +
7403a8779fbSJames Wright                                     stab[j][2] * dXdx[k][2]);
7413a8779fbSJames Wright       break;
7423a8779fbSJames Wright     }
743*752f40e3SJed Brown     for (int j=0; j<5; j++) jac_data[j][i] = U[j];
744*752f40e3SJed Brown     for (int j=0; j<3; j++) jac_data[5+j][i] = Tau_d[j];
7453a8779fbSJames Wright 
7463a8779fbSJames Wright   } // End Quadrature Point Loop
7473a8779fbSJames Wright 
7483a8779fbSJames Wright   // Return
7493a8779fbSJames Wright   return 0;
7503a8779fbSJames Wright }
7513a8779fbSJames Wright // *****************************************************************************
7523a8779fbSJames Wright #endif // newtonian_h
753