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" 18475b2820SJames Wright #include "newtonian_state.h" 193a8779fbSJames Wright 203a8779fbSJames Wright #ifndef M_PI 213a8779fbSJames Wright #define M_PI 3.14159265358979323846 223a8779fbSJames Wright #endif 233a8779fbSJames Wright 243a8779fbSJames Wright // ***************************************************************************** 253a8779fbSJames Wright // Helper function for computing flux Jacobian 263a8779fbSJames Wright // ***************************************************************************** 273a8779fbSJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NS(CeedScalar dF[3][5][5], 283a8779fbSJames Wright const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 29bb8a0c61SJames Wright const CeedScalar gamma, const CeedScalar g[3], const CeedScalar x[3]) { 303a8779fbSJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 31bb8a0c61SJames Wright CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]); 323a8779fbSJames Wright for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 333a8779fbSJames Wright for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix 34bb8a0c61SJames Wright dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2. - e_potential)) : 0.) - 35bb8a0c61SJames Wright u[i]*u[j]; 363a8779fbSJames Wright for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix 373a8779fbSJames Wright dF[i][0][k+1] = ((i==k) ? 1. : 0.); 383a8779fbSJames Wright dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) + 393a8779fbSJames Wright ((i==k) ? u[j] : 0.) - 403a8779fbSJames Wright ((i==j) ? u[k] : 0.) * (gamma-1.); 413a8779fbSJames Wright dF[i][4][k+1] = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) - 423a8779fbSJames Wright (gamma-1.)*u[i]*u[k]; 433a8779fbSJames Wright } 443a8779fbSJames Wright dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.); 453a8779fbSJames Wright } 463a8779fbSJames Wright dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho); 473a8779fbSJames Wright dF[i][4][4] = u[i] * gamma; 483a8779fbSJames Wright } 493a8779fbSJames Wright } 503a8779fbSJames Wright 513a8779fbSJames Wright // ***************************************************************************** 52bb8a0c61SJames Wright // Helper function for computing flux Jacobian of Primitive variables 53bb8a0c61SJames Wright // ***************************************************************************** 54bb8a0c61SJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NSp(CeedScalar dF[3][5][5], 55bb8a0c61SJames Wright const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 56bb8a0c61SJames Wright const CeedScalar Rd, const CeedScalar cv) { 57bb8a0c61SJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 58bb8a0c61SJames Wright // TODO Add in gravity's contribution 59bb8a0c61SJames Wright 60bb8a0c61SJames Wright CeedScalar T = ( E / rho - u_sq / 2. ) / cv; 61bb8a0c61SJames Wright CeedScalar drdT = -rho / T; 62bb8a0c61SJames Wright CeedScalar drdP = 1. / ( Rd * T); 63bb8a0c61SJames Wright CeedScalar etot = E / rho ; 64bb8a0c61SJames Wright CeedScalar e2p = drdP * etot + 1. ; 65bb8a0c61SJames Wright CeedScalar e3p = ( E + rho * Rd * T ); 66bb8a0c61SJames Wright CeedScalar e4p = drdT * etot + rho * cv ; 67bb8a0c61SJames Wright 68bb8a0c61SJames Wright for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 69bb8a0c61SJames Wright for (CeedInt j=0; j<3; j++) { // j counts F^{m_j} 70bb8a0c61SJames Wright // [row][col] of A_i 71bb8a0c61SJames Wright dF[i][j+1][0] = drdP * u[i] * u[j] + ((i==j) ? 1. : 0.); // F^{{m_j} wrt p 72bb8a0c61SJames Wright for (CeedInt k=0; k<3; k++) { // k counts the wrt vel_k 732acc7cbcSKenneth E. Jansen dF[i][0][k+1] = ((i==k) ? rho : 0.); // F^c wrt u_k 74bb8a0c61SJames Wright dF[i][j+1][k+1] = (((j==k) ? u[i] : 0.) + // F^m_j wrt u_k 75bb8a0c61SJames Wright ((i==k) ? u[j] : 0.) ) * rho; 76bb8a0c61SJames Wright dF[i][4][k+1] = rho * u[i] * u[k] 77bb8a0c61SJames Wright + ((i==k) ? e3p : 0.) ; // F^e wrt u_k 78bb8a0c61SJames Wright } 79bb8a0c61SJames Wright dF[i][j+1][4] = drdT * u[i] * u[j]; // F^{m_j} wrt T 80bb8a0c61SJames Wright } 81bb8a0c61SJames Wright dF[i][4][0] = u[i] * e2p; // F^e wrt p 82bb8a0c61SJames Wright dF[i][4][4] = u[i] * e4p; // F^e wrt T 83bb8a0c61SJames Wright dF[i][0][0] = u[i] * drdP; // F^c wrt p 84bb8a0c61SJames Wright dF[i][0][4] = u[i] * drdT; // F^c wrt T 85bb8a0c61SJames Wright } 86bb8a0c61SJames Wright } 87bb8a0c61SJames Wright 88bb8a0c61SJames Wright CEED_QFUNCTION_HELPER void PrimitiveToConservative_fwd(const CeedScalar rho, 89bb8a0c61SJames Wright const CeedScalar u[3], const CeedScalar E, const CeedScalar Rd, 90bb8a0c61SJames Wright const CeedScalar cv, const CeedScalar dY[5], CeedScalar dU[5]) { 91bb8a0c61SJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; 92bb8a0c61SJames Wright CeedScalar T = ( E / rho - u_sq / 2. ) / cv; 93bb8a0c61SJames Wright CeedScalar drdT = -rho / T; 94bb8a0c61SJames Wright CeedScalar drdP = 1. / ( Rd * T); 95bb8a0c61SJames Wright dU[0] = drdP * dY[0] + drdT * dY[4]; 96bb8a0c61SJames Wright CeedScalar de_kinetic = 0; 97493642f1SJames Wright for (CeedInt i=0; i<3; i++) { 98bb8a0c61SJames Wright dU[1+i] = dU[0] * u[i] + rho * dY[1+i]; 99bb8a0c61SJames Wright de_kinetic += u[i] * dY[1+i]; 100bb8a0c61SJames Wright } 101bb8a0c61SJames Wright dU[4] = rho * cv * dY[4] + dU[0] * cv * T // internal energy: rho * e 102bb8a0c61SJames Wright + rho * de_kinetic + .5 * dU[0] * u_sq; // kinetic energy: .5 * rho * |u|^2 103bb8a0c61SJames Wright } 104bb8a0c61SJames Wright 105bb8a0c61SJames Wright // ***************************************************************************** 106bb8a0c61SJames Wright // Helper function for computing Tau elements (stabilization constant) 107bb8a0c61SJames Wright // Model from: 108bb8a0c61SJames Wright // PHASTA 109bb8a0c61SJames Wright // 110bb8a0c61SJames Wright // Tau[i] = itau=0 which is diagonal-Shakib (3 values still but not spatial) 111bb8a0c61SJames Wright // 112bb8a0c61SJames Wright // Where NOT UPDATED YET 113bb8a0c61SJames Wright // ***************************************************************************** 114bb8a0c61SJames Wright CEED_QFUNCTION_HELPER void Tau_diagPrim(CeedScalar Tau_d[3], 115bb8a0c61SJames Wright const CeedScalar dXdx[3][3], const CeedScalar u[3], 116bb8a0c61SJames Wright const CeedScalar cv, const NewtonianIdealGasContext newt_ctx, 117bb8a0c61SJames Wright const CeedScalar mu, const CeedScalar dt, 118bb8a0c61SJames Wright const CeedScalar rho) { 119bb8a0c61SJames Wright // Context 120bb8a0c61SJames Wright const CeedScalar Ctau_t = newt_ctx->Ctau_t; 121bb8a0c61SJames Wright const CeedScalar Ctau_v = newt_ctx->Ctau_v; 122bb8a0c61SJames Wright const CeedScalar Ctau_C = newt_ctx->Ctau_C; 123bb8a0c61SJames Wright const CeedScalar Ctau_M = newt_ctx->Ctau_M; 124bb8a0c61SJames Wright const CeedScalar Ctau_E = newt_ctx->Ctau_E; 125bb8a0c61SJames Wright CeedScalar gijd[6]; 126bb8a0c61SJames Wright CeedScalar tau; 127bb8a0c61SJames Wright CeedScalar dts; 128bb8a0c61SJames Wright CeedScalar fact; 129bb8a0c61SJames Wright 130bb8a0c61SJames Wright //*INDENT-OFF* 131bb8a0c61SJames Wright gijd[0] = dXdx[0][0] * dXdx[0][0] 132bb8a0c61SJames Wright + dXdx[1][0] * dXdx[1][0] 133bb8a0c61SJames Wright + dXdx[2][0] * dXdx[2][0]; 134bb8a0c61SJames Wright 135bb8a0c61SJames Wright gijd[1] = dXdx[0][0] * dXdx[0][1] 136bb8a0c61SJames Wright + dXdx[1][0] * dXdx[1][1] 137bb8a0c61SJames Wright + dXdx[2][0] * dXdx[2][1]; 138bb8a0c61SJames Wright 139bb8a0c61SJames Wright gijd[2] = dXdx[0][1] * dXdx[0][1] 140bb8a0c61SJames Wright + dXdx[1][1] * dXdx[1][1] 141bb8a0c61SJames Wright + dXdx[2][1] * dXdx[2][1]; 142bb8a0c61SJames Wright 143bb8a0c61SJames Wright gijd[3] = dXdx[0][0] * dXdx[0][2] 144bb8a0c61SJames Wright + dXdx[1][0] * dXdx[1][2] 145bb8a0c61SJames Wright + dXdx[2][0] * dXdx[2][2]; 146bb8a0c61SJames Wright 147bb8a0c61SJames Wright gijd[4] = dXdx[0][1] * dXdx[0][2] 148bb8a0c61SJames Wright + dXdx[1][1] * dXdx[1][2] 149bb8a0c61SJames Wright + dXdx[2][1] * dXdx[2][2]; 150bb8a0c61SJames Wright 151bb8a0c61SJames Wright gijd[5] = dXdx[0][2] * dXdx[0][2] 152bb8a0c61SJames Wright + dXdx[1][2] * dXdx[1][2] 153bb8a0c61SJames Wright + dXdx[2][2] * dXdx[2][2]; 154bb8a0c61SJames Wright //*INDENT-ON* 155bb8a0c61SJames Wright 156bb8a0c61SJames Wright dts = Ctau_t / dt ; 157bb8a0c61SJames Wright 158bb8a0c61SJames Wright tau = rho*rho*((4. * dts * dts) 159bb8a0c61SJames Wright + u[0] * ( u[0] * gijd[0] + 2. * ( u[1] * gijd[1] + u[2] * gijd[3])) 160bb8a0c61SJames Wright + u[1] * ( u[1] * gijd[2] + 2. * u[2] * gijd[4]) 161bb8a0c61SJames Wright + u[2] * u[2] * gijd[5]) 162bb8a0c61SJames Wright + Ctau_v* mu * mu * 163bb8a0c61SJames Wright (gijd[0]*gijd[0] + gijd[2]*gijd[2] + gijd[5]*gijd[5] + 164bb8a0c61SJames Wright + 2. * (gijd[1]*gijd[1] + gijd[3]*gijd[3] + gijd[4]*gijd[4])); 165bb8a0c61SJames Wright 166bb8a0c61SJames Wright fact=sqrt(tau); 167bb8a0c61SJames Wright 168bb8a0c61SJames Wright Tau_d[0] = Ctau_C * fact / (rho*(gijd[0] + gijd[2] + gijd[5]))*0.125; 169bb8a0c61SJames Wright 170bb8a0c61SJames Wright Tau_d[1] = Ctau_M / fact; 171bb8a0c61SJames Wright Tau_d[2] = Ctau_E / ( fact * cv ); 172bb8a0c61SJames Wright 173bb8a0c61SJames Wright // consider putting back the way I initially had it Ctau_E * Tau_d[1] /cv 174bb8a0c61SJames Wright // to avoid a division if the compiler is smart enough to see that cv IS 175bb8a0c61SJames Wright // a constant that it could invert once for all elements 176bb8a0c61SJames Wright // but in that case energy tau is scaled by the product of Ctau_E * Ctau_M 177bb8a0c61SJames Wright // OR we could absorb cv into Ctau_E but this puts more burden on user to 178bb8a0c61SJames Wright // know how to change constants with a change of fluid or units. Same for 179bb8a0c61SJames Wright // Ctau_v * mu * mu IF AND ONLY IF we don't add viscosity law =f(T) 180bb8a0c61SJames Wright } 181bb8a0c61SJames Wright 182bb8a0c61SJames Wright // ***************************************************************************** 1833a8779fbSJames Wright // This QFunction sets a "still" initial condition for generic Newtonian IG problems 1843a8779fbSJames Wright // ***************************************************************************** 1853a8779fbSJames Wright CEED_QFUNCTION(ICsNewtonianIG)(void *ctx, CeedInt Q, 1863a8779fbSJames Wright const CeedScalar *const *in, CeedScalar *const *out) { 1873a8779fbSJames Wright // Inputs 1883a8779fbSJames Wright const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 1893a8779fbSJames Wright 1903a8779fbSJames Wright // Outputs 1913a8779fbSJames Wright CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 1923a8779fbSJames Wright 193bb8a0c61SJames Wright // Context 194bb8a0c61SJames Wright const SetupContext context = (SetupContext)ctx; 195bb8a0c61SJames Wright const CeedScalar theta0 = context->theta0; 196bb8a0c61SJames Wright const CeedScalar P0 = context->P0; 197bb8a0c61SJames Wright const CeedScalar cv = context->cv; 198bb8a0c61SJames Wright const CeedScalar cp = context->cp; 199bb8a0c61SJames Wright const CeedScalar *g = context->g; 200bb8a0c61SJames Wright const CeedScalar Rd = cp - cv; 201bb8a0c61SJames Wright 2023a8779fbSJames Wright // Quadrature Point Loop 2033a8779fbSJames Wright CeedPragmaSIMD 2043a8779fbSJames Wright for (CeedInt i=0; i<Q; i++) { 2053a8779fbSJames Wright CeedScalar q[5] = {0.}; 2063a8779fbSJames Wright 2073a8779fbSJames Wright // Setup 2083a8779fbSJames Wright // -- Coordinates 209bb8a0c61SJames Wright const CeedScalar x[3] = {X[0][i], X[1][i], X[2][i]}; 210bb8a0c61SJames Wright const CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]); 2113a8779fbSJames Wright 2123a8779fbSJames Wright // -- Density 213bb8a0c61SJames Wright const CeedScalar rho = P0 / (Rd*theta0); 2143a8779fbSJames Wright 2153a8779fbSJames Wright // Initial Conditions 2163a8779fbSJames Wright q[0] = rho; 2173a8779fbSJames Wright q[1] = 0.0; 2183a8779fbSJames Wright q[2] = 0.0; 2193a8779fbSJames Wright q[3] = 0.0; 220bb8a0c61SJames Wright q[4] = rho * (cv*theta0 + e_potential); 2213a8779fbSJames Wright 2223a8779fbSJames Wright for (CeedInt j=0; j<5; j++) 2233a8779fbSJames Wright q0[j][i] = q[j]; 2243a8779fbSJames Wright } // End of Quadrature Point Loop 2253a8779fbSJames Wright return 0; 2263a8779fbSJames Wright } 2273a8779fbSJames Wright 2283a8779fbSJames Wright // ***************************************************************************** 2293a8779fbSJames Wright // This QFunction implements the following formulation of Navier-Stokes with 2303a8779fbSJames Wright // explicit time stepping method 2313a8779fbSJames Wright // 2323a8779fbSJames Wright // This is 3D compressible Navier-Stokes in conservation form with state 2333a8779fbSJames Wright // variables of density, momentum density, and total energy density. 2343a8779fbSJames Wright // 2353a8779fbSJames Wright // State Variables: q = ( rho, U1, U2, U3, E ) 2363a8779fbSJames Wright // rho - Mass Density 2373a8779fbSJames Wright // Ui - Momentum Density, Ui = rho ui 2383a8779fbSJames Wright // E - Total Energy Density, E = rho (cv T + (u u)/2 + g z) 2393a8779fbSJames Wright // 2403a8779fbSJames Wright // Navier-Stokes Equations: 2413a8779fbSJames Wright // drho/dt + div( U ) = 0 2423a8779fbSJames Wright // dU/dt + div( rho (u x u) + P I3 ) + rho g khat = div( Fu ) 2433a8779fbSJames Wright // dE/dt + div( (E + P) u ) = div( Fe ) 2443a8779fbSJames Wright // 2453a8779fbSJames Wright // Viscous Stress: 2463a8779fbSJames Wright // Fu = mu (grad( u ) + grad( u )^T + lambda div ( u ) I3) 2473a8779fbSJames Wright // 2483a8779fbSJames Wright // Thermal Stress: 2493a8779fbSJames Wright // Fe = u Fu + k grad( T ) 250bb8a0c61SJames Wright // Equation of State 2513a8779fbSJames Wright // P = (gamma - 1) (E - rho (u u) / 2 - rho g z) 2523a8779fbSJames Wright // 2533a8779fbSJames Wright // Stabilization: 2543a8779fbSJames Wright // Tau = diag(TauC, TauM, TauM, TauM, TauE) 2553a8779fbSJames Wright // f1 = rho sqrt(ui uj gij) 2563a8779fbSJames Wright // gij = dXi/dX * dXi/dX 2573a8779fbSJames Wright // TauC = Cc f1 / (8 gii) 2583a8779fbSJames Wright // TauM = min( 1 , 1 / f1 ) 2593a8779fbSJames Wright // TauE = TauM / (Ce cv) 2603a8779fbSJames Wright // 2613a8779fbSJames Wright // SU = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) ) 2623a8779fbSJames Wright // 2633a8779fbSJames Wright // Constants: 2643a8779fbSJames Wright // lambda = - 2 / 3, From Stokes hypothesis 2653a8779fbSJames Wright // mu , Dynamic viscosity 2663a8779fbSJames Wright // k , Thermal conductivity 2673a8779fbSJames Wright // cv , Specific heat, constant volume 2683a8779fbSJames Wright // cp , Specific heat, constant pressure 2693a8779fbSJames Wright // g , Gravity 2703a8779fbSJames Wright // gamma = cp / cv, Specific heat ratio 2713a8779fbSJames Wright // 2723a8779fbSJames Wright // We require the product of the inverse of the Jacobian (dXdx_j,k) and 2733a8779fbSJames Wright // its transpose (dXdx_k,j) to properly compute integrals of the form: 2743a8779fbSJames Wright // int( gradv gradu ) 2753a8779fbSJames Wright // 2763a8779fbSJames Wright // ***************************************************************************** 277c1a52365SJed Brown CEED_QFUNCTION(RHSFunction_Newtonian)(void *ctx, CeedInt Q, 2783a8779fbSJames Wright const CeedScalar *const *in, CeedScalar *const *out) { 2793a8779fbSJames Wright // *INDENT-OFF* 2803a8779fbSJames Wright // Inputs 2813a8779fbSJames Wright const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 282752f40e3SJed Brown (*Grad_q)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 2833a8779fbSJames Wright (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 2843a8779fbSJames Wright (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 2853a8779fbSJames Wright // Outputs 2863a8779fbSJames Wright CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 287752f40e3SJed Brown (*Grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 2883a8779fbSJames Wright // *INDENT-ON* 2893a8779fbSJames Wright 2903a8779fbSJames Wright // Context 2913a8779fbSJames Wright NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx; 2923a8779fbSJames Wright const CeedScalar mu = context->mu; 2933a8779fbSJames Wright const CeedScalar cv = context->cv; 2943a8779fbSJames Wright const CeedScalar cp = context->cp; 295bb8a0c61SJames Wright const CeedScalar *g = context->g; 296bb8a0c61SJames Wright const CeedScalar dt = context->dt; 2973a8779fbSJames Wright const CeedScalar gamma = cp / cv; 298bb8a0c61SJames Wright const CeedScalar Rd = cp - cv; 2993a8779fbSJames Wright 3003a8779fbSJames Wright CeedPragmaSIMD 3013a8779fbSJames Wright // Quadrature Point Loop 3023a8779fbSJames Wright for (CeedInt i=0; i<Q; i++) { 303c1a52365SJed Brown CeedScalar U[5]; 304c1a52365SJed Brown for (int j=0; j<5; j++) U[j] = q[j][i]; 305c1a52365SJed Brown const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 306c1a52365SJed Brown State s = StateFromU(context, U, x_i); 307c1a52365SJed Brown 3083a8779fbSJames Wright // -- Interp-to-Interp q_data 3093a8779fbSJames Wright const CeedScalar wdetJ = q_data[0][i]; 3103a8779fbSJames Wright // -- Interp-to-Grad q_data 3113a8779fbSJames Wright // ---- Inverse of change of coordinate matrix: X_i,j 3123a8779fbSJames Wright // *INDENT-OFF* 3133a8779fbSJames Wright const CeedScalar dXdx[3][3] = {{q_data[1][i], 3143a8779fbSJames Wright q_data[2][i], 3153a8779fbSJames Wright q_data[3][i]}, 3163a8779fbSJames Wright {q_data[4][i], 3173a8779fbSJames Wright q_data[5][i], 3183a8779fbSJames Wright q_data[6][i]}, 3193a8779fbSJames Wright {q_data[7][i], 3203a8779fbSJames Wright q_data[8][i], 3213a8779fbSJames Wright q_data[9][i]} 3223a8779fbSJames Wright }; 3233a8779fbSJames Wright // *INDENT-ON* 3243a8779fbSJames Wright 325c1a52365SJed Brown State grad_s[3]; 326eef2387dSJed Brown for (CeedInt j=0; j<3; j++) { 3272f7ce6c1SJed Brown CeedScalar dx_i[3] = {0}, dU[5]; 3282556a851SJed Brown for (CeedInt k=0; k<5; k++) 3292556a851SJed Brown dU[k] = Grad_q[0][k][i] * dXdx[0][j] + 3302556a851SJed Brown Grad_q[1][k][i] * dXdx[1][j] + 3312556a851SJed Brown Grad_q[2][k][i] * dXdx[2][j]; 332c1a52365SJed Brown dx_i[j] = 1.; 3332f7ce6c1SJed Brown grad_s[j] = StateFromU_fwd(context, s, dU, x_i, dx_i); 334c1a52365SJed Brown } 335c1a52365SJed Brown 336c1a52365SJed Brown CeedScalar strain_rate[6], kmstress[6], stress[3][3], Fe[3]; 337c1a52365SJed Brown KMStrainRate(grad_s, strain_rate); 338c1a52365SJed Brown NewtonianStress(context, strain_rate, kmstress); 339c1a52365SJed Brown KMUnpack(kmstress, stress); 340c1a52365SJed Brown ViscousEnergyFlux(context, s.Y, grad_s, stress, Fe); 341c1a52365SJed Brown 342c1a52365SJed Brown StateConservative F_inviscid[3]; 343c1a52365SJed Brown FluxInviscid(context, s, F_inviscid); 344c1a52365SJed Brown 345c1a52365SJed Brown // Total flux 346c1a52365SJed Brown CeedScalar Flux[5][3]; 347eef2387dSJed Brown for (CeedInt j=0; j<3; j++) { 348c1a52365SJed Brown Flux[0][j] = F_inviscid[j].density; 349eef2387dSJed Brown for (CeedInt k=0; k<3; k++) 350c1a52365SJed Brown Flux[k+1][j] = F_inviscid[j].momentum[k] - stress[k][j]; 351c1a52365SJed Brown Flux[4][j] = F_inviscid[j].E_total + Fe[j]; 352c1a52365SJed Brown } 353c1a52365SJed Brown 354eef2387dSJed Brown for (CeedInt j=0; j<3; j++) { 355eef2387dSJed Brown for (CeedInt k=0; k<5; k++) { 356752f40e3SJed Brown Grad_v[j][k][i] = wdetJ * (dXdx[j][0] * Flux[k][0] + 357c1a52365SJed Brown dXdx[j][1] * Flux[k][1] + 358c1a52365SJed Brown dXdx[j][2] * Flux[k][2]); 359c1a52365SJed Brown } 360c1a52365SJed Brown } 361c1a52365SJed Brown 362c1a52365SJed Brown const CeedScalar body_force[5] = {0, s.U.density *g[0], s.U.density *g[1], s.U.density *g[2], 0}; 363c1a52365SJed Brown for (int j=0; j<5; j++) 364c1a52365SJed Brown v[j][i] = wdetJ * body_force[j]; 3653a8779fbSJames Wright 3663a8779fbSJames Wright // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 367c1a52365SJed Brown CeedScalar jacob_F_conv[3][5][5] = {0}; 368c1a52365SJed Brown computeFluxJacobian_NS(jacob_F_conv, s.U.density, s.Y.velocity, s.U.E_total, 369c1a52365SJed Brown gamma, g, x_i); 370c1a52365SJed Brown CeedScalar grad_U[5][3]; 371493642f1SJames Wright for (CeedInt j=0; j<3; j++) { 372c1a52365SJed Brown grad_U[0][j] = grad_s[j].U.density; 373eef2387dSJed Brown for (CeedInt k=0; k<3; k++) grad_U[k+1][j] = grad_s[j].U.momentum[k]; 374c1a52365SJed Brown grad_U[4][j] = grad_s[j].U.E_total; 3753a8779fbSJames Wright } 3763a8779fbSJames Wright 3773a8779fbSJames Wright // strong_conv = dF/dq * dq/dx (Strong convection) 3783a8779fbSJames Wright CeedScalar strong_conv[5] = {0}; 379493642f1SJames Wright for (CeedInt j=0; j<3; j++) 380493642f1SJames Wright for (CeedInt k=0; k<5; k++) 381493642f1SJames Wright for (CeedInt l=0; l<5; l++) 382c1a52365SJed Brown strong_conv[k] += jacob_F_conv[j][k][l] * grad_U[l][j]; 3833a8779fbSJames Wright 384bb8a0c61SJames Wright // -- Stabilization method: none, SU, or SUPG 385bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 386bb8a0c61SJames Wright CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0}; 387bb8a0c61SJames Wright CeedScalar Tau_d[3] = {0.}; 3883a8779fbSJames Wright switch (context->stabilization) { 3893a8779fbSJames Wright case STAB_NONE: // Galerkin 3903a8779fbSJames Wright break; 3913a8779fbSJames Wright case STAB_SU: // SU 392c1a52365SJed Brown Tau_diagPrim(Tau_d, dXdx, s.Y.velocity, cv, context, mu, dt, s.U.density); 393bb8a0c61SJames Wright tau_strong_conv[0] = Tau_d[0] * strong_conv[0]; 394bb8a0c61SJames Wright tau_strong_conv[1] = Tau_d[1] * strong_conv[1]; 395bb8a0c61SJames Wright tau_strong_conv[2] = Tau_d[1] * strong_conv[2]; 396bb8a0c61SJames Wright tau_strong_conv[3] = Tau_d[1] * strong_conv[3]; 397bb8a0c61SJames Wright tau_strong_conv[4] = Tau_d[2] * strong_conv[4]; 398c1a52365SJed Brown PrimitiveToConservative_fwd(s.U.density, s.Y.velocity, s.U.E_total, Rd, cv, 399c1a52365SJed Brown tau_strong_conv, 400bb8a0c61SJames Wright tau_strong_conv_conservative); 401493642f1SJames Wright for (CeedInt j=0; j<3; j++) 402493642f1SJames Wright for (CeedInt k=0; k<5; k++) 403493642f1SJames Wright for (CeedInt l=0; l<5; l++) 404bb8a0c61SJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l]; 4053a8779fbSJames Wright 406493642f1SJames Wright for (CeedInt j=0; j<5; j++) 407493642f1SJames Wright for (CeedInt k=0; k<3; k++) 408752f40e3SJed Brown Grad_v[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] + 4093a8779fbSJames Wright stab[j][1] * dXdx[k][1] + 4103a8779fbSJames Wright stab[j][2] * dXdx[k][2]); 4113a8779fbSJames Wright break; 4123a8779fbSJames Wright case STAB_SUPG: // SUPG is not implemented for explicit scheme 4133a8779fbSJames Wright break; 4143a8779fbSJames Wright } 4153a8779fbSJames Wright 4163a8779fbSJames Wright } // End Quadrature Point Loop 4173a8779fbSJames Wright 4183a8779fbSJames Wright // Return 4193a8779fbSJames Wright return 0; 4203a8779fbSJames Wright } 4213a8779fbSJames Wright 4223a8779fbSJames Wright // ***************************************************************************** 4233a8779fbSJames Wright // This QFunction implements the Navier-Stokes equations (mentioned above) with 4243a8779fbSJames Wright // implicit time stepping method 4253a8779fbSJames Wright // 4263a8779fbSJames Wright // SU = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) ) 4273a8779fbSJames Wright // SUPG = Galerkin + grad(v) . ( Ai^T * Tau * (q_dot + Aj q,j - body force) ) 4283a8779fbSJames Wright // (diffussive terms will be added later) 4293a8779fbSJames Wright // 4303a8779fbSJames Wright // ***************************************************************************** 4313a8779fbSJames Wright CEED_QFUNCTION(IFunction_Newtonian)(void *ctx, CeedInt Q, 4323a8779fbSJames Wright const CeedScalar *const *in, 4333a8779fbSJames Wright CeedScalar *const *out) { 4343a8779fbSJames Wright // *INDENT-OFF* 4353a8779fbSJames Wright // Inputs 4363a8779fbSJames Wright const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 437752f40e3SJed Brown (*Grad_q)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 4383a8779fbSJames Wright (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 4393a8779fbSJames Wright (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3], 4403a8779fbSJames Wright (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 4413a8779fbSJames Wright // Outputs 4423a8779fbSJames Wright CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 443752f40e3SJed Brown (*Grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1], 444752f40e3SJed Brown (*jac_data)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[2]; 4453a8779fbSJames Wright // *INDENT-ON* 4463a8779fbSJames Wright // Context 4473a8779fbSJames Wright NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx; 4483a8779fbSJames Wright const CeedScalar mu = context->mu; 4493a8779fbSJames Wright const CeedScalar cv = context->cv; 4503a8779fbSJames Wright const CeedScalar cp = context->cp; 451bb8a0c61SJames Wright const CeedScalar *g = context->g; 452bb8a0c61SJames Wright const CeedScalar dt = context->dt; 4533a8779fbSJames Wright const CeedScalar gamma = cp / cv; 454bb8a0c61SJames Wright const CeedScalar Rd = cp-cv; 4553a8779fbSJames Wright 4563a8779fbSJames Wright CeedPragmaSIMD 4573a8779fbSJames Wright // Quadrature Point Loop 4583a8779fbSJames Wright for (CeedInt i=0; i<Q; i++) { 459c1a52365SJed Brown CeedScalar U[5]; 460eef2387dSJed Brown for (CeedInt j=0; j<5; j++) U[j] = q[j][i]; 461c1a52365SJed Brown const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 462c1a52365SJed Brown State s = StateFromU(context, U, x_i); 463c1a52365SJed Brown 4643a8779fbSJames Wright // -- Interp-to-Interp q_data 4653a8779fbSJames Wright const CeedScalar wdetJ = q_data[0][i]; 4663a8779fbSJames Wright // -- Interp-to-Grad q_data 4673a8779fbSJames Wright // ---- Inverse of change of coordinate matrix: X_i,j 4683a8779fbSJames Wright // *INDENT-OFF* 4693a8779fbSJames Wright const CeedScalar dXdx[3][3] = {{q_data[1][i], 4703a8779fbSJames Wright q_data[2][i], 4713a8779fbSJames Wright q_data[3][i]}, 4723a8779fbSJames Wright {q_data[4][i], 4733a8779fbSJames Wright q_data[5][i], 4743a8779fbSJames Wright q_data[6][i]}, 4753a8779fbSJames Wright {q_data[7][i], 4763a8779fbSJames Wright q_data[8][i], 4773a8779fbSJames Wright q_data[9][i]} 4783a8779fbSJames Wright }; 4793a8779fbSJames Wright // *INDENT-ON* 480c1a52365SJed Brown State grad_s[3]; 481493642f1SJames Wright for (CeedInt j=0; j<3; j++) { 4822f7ce6c1SJed Brown CeedScalar dx_i[3] = {0}, dU[5]; 4832556a851SJed Brown for (CeedInt k=0; k<5; k++) 4842556a851SJed Brown dU[k] = Grad_q[0][k][i] * dXdx[0][j] + 4852556a851SJed Brown Grad_q[1][k][i] * dXdx[1][j] + 4862556a851SJed Brown Grad_q[2][k][i] * dXdx[2][j]; 487c1a52365SJed Brown dx_i[j] = 1.; 4882f7ce6c1SJed Brown grad_s[j] = StateFromU_fwd(context, s, dU, x_i, dx_i); 4893a8779fbSJames Wright } 490c1a52365SJed Brown 491c1a52365SJed Brown CeedScalar strain_rate[6], kmstress[6], stress[3][3], Fe[3]; 492c1a52365SJed Brown KMStrainRate(grad_s, strain_rate); 493c1a52365SJed Brown NewtonianStress(context, strain_rate, kmstress); 494c1a52365SJed Brown KMUnpack(kmstress, stress); 495c1a52365SJed Brown ViscousEnergyFlux(context, s.Y, grad_s, stress, Fe); 496c1a52365SJed Brown 497c1a52365SJed Brown StateConservative F_inviscid[3]; 498c1a52365SJed Brown FluxInviscid(context, s, F_inviscid); 499c1a52365SJed Brown 500c1a52365SJed Brown 501c1a52365SJed Brown // Total flux 502c1a52365SJed Brown CeedScalar Flux[5][3]; 503eef2387dSJed Brown for (CeedInt j=0; j<3; j++) { 504c1a52365SJed Brown Flux[0][j] = F_inviscid[j].density; 505493642f1SJames Wright for (CeedInt k=0; k<3; k++) 506c1a52365SJed Brown Flux[k+1][j] = F_inviscid[j].momentum[k] - stress[k][j]; 507c1a52365SJed Brown Flux[4][j] = F_inviscid[j].E_total + Fe[j]; 508c1a52365SJed Brown } 509c1a52365SJed Brown 510eef2387dSJed Brown for (CeedInt j=0; j<3; j++) { 511eef2387dSJed Brown for (CeedInt k=0; k<5; k++) { 512752f40e3SJed Brown Grad_v[j][k][i] = -wdetJ * (dXdx[j][0] * Flux[k][0] + 513c1a52365SJed Brown dXdx[j][1] * Flux[k][1] + 514c1a52365SJed Brown dXdx[j][2] * Flux[k][2]); 515c1a52365SJed Brown } 516c1a52365SJed Brown } 517c1a52365SJed Brown 518c1a52365SJed Brown const CeedScalar body_force[5] = {0, s.U.density *g[0], s.U.density *g[1], s.U.density *g[2], 0}; 519eef2387dSJed Brown for (CeedInt j=0; j<5; j++) 520c1a52365SJed Brown v[j][i] = wdetJ * (q_dot[j][i] - body_force[j]); 5213a8779fbSJames Wright 5223a8779fbSJames Wright // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 523c1a52365SJed Brown CeedScalar jacob_F_conv[3][5][5] = {0}; 524c1a52365SJed Brown computeFluxJacobian_NS(jacob_F_conv, s.U.density, s.Y.velocity, s.U.E_total, 525c1a52365SJed Brown gamma, g, x_i); 526c1a52365SJed Brown CeedScalar grad_U[5][3]; 527493642f1SJames Wright for (CeedInt j=0; j<3; j++) { 528c1a52365SJed Brown grad_U[0][j] = grad_s[j].U.density; 529eef2387dSJed Brown for (CeedInt k=0; k<3; k++) grad_U[k+1][j] = grad_s[j].U.momentum[k]; 530c1a52365SJed Brown grad_U[4][j] = grad_s[j].U.E_total; 5313a8779fbSJames Wright } 532c1a52365SJed Brown 5333a8779fbSJames Wright // strong_conv = dF/dq * dq/dx (Strong convection) 5343a8779fbSJames Wright CeedScalar strong_conv[5] = {0}; 535493642f1SJames Wright for (CeedInt j=0; j<3; j++) 536493642f1SJames Wright for (CeedInt k=0; k<5; k++) 537493642f1SJames Wright for (CeedInt l=0; l<5; l++) 538c1a52365SJed Brown strong_conv[k] += jacob_F_conv[j][k][l] * grad_U[l][j]; 5393a8779fbSJames Wright 5403a8779fbSJames Wright // Strong residual 5413a8779fbSJames Wright CeedScalar strong_res[5]; 542493642f1SJames Wright for (CeedInt j=0; j<5; j++) 5433a8779fbSJames Wright strong_res[j] = q_dot[j][i] + strong_conv[j] - body_force[j]; 5443a8779fbSJames Wright 5453a8779fbSJames Wright // -- Stabilization method: none, SU, or SUPG 546bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 547bb8a0c61SJames Wright CeedScalar tau_strong_res[5] = {0.}, tau_strong_res_conservative[5] = {0}; 548bb8a0c61SJames Wright CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0}; 549bb8a0c61SJames Wright CeedScalar Tau_d[3] = {0.}; 5503a8779fbSJames Wright switch (context->stabilization) { 5513a8779fbSJames Wright case STAB_NONE: // Galerkin 5523a8779fbSJames Wright break; 5533a8779fbSJames Wright case STAB_SU: // SU 554c1a52365SJed Brown Tau_diagPrim(Tau_d, dXdx, s.Y.velocity, cv, context, mu, dt, s.U.density); 555bb8a0c61SJames Wright tau_strong_conv[0] = Tau_d[0] * strong_conv[0]; 556bb8a0c61SJames Wright tau_strong_conv[1] = Tau_d[1] * strong_conv[1]; 557bb8a0c61SJames Wright tau_strong_conv[2] = Tau_d[1] * strong_conv[2]; 558bb8a0c61SJames Wright tau_strong_conv[3] = Tau_d[1] * strong_conv[3]; 559bb8a0c61SJames Wright tau_strong_conv[4] = Tau_d[2] * strong_conv[4]; 560c1a52365SJed Brown PrimitiveToConservative_fwd(s.U.density, s.Y.velocity, s.U.E_total, Rd, cv, 561c1a52365SJed Brown tau_strong_conv, tau_strong_conv_conservative); 562493642f1SJames Wright for (CeedInt j=0; j<3; j++) 563493642f1SJames Wright for (CeedInt k=0; k<5; k++) 564493642f1SJames Wright for (CeedInt l=0; l<5; l++) 565bb8a0c61SJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l]; 5663a8779fbSJames Wright 567493642f1SJames Wright for (CeedInt j=0; j<5; j++) 568493642f1SJames Wright for (CeedInt k=0; k<3; k++) 569752f40e3SJed Brown Grad_v[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 5703a8779fbSJames Wright stab[j][1] * dXdx[k][1] + 5713a8779fbSJames Wright stab[j][2] * dXdx[k][2]); 572eef2387dSJed Brown 5733a8779fbSJames Wright break; 5743a8779fbSJames Wright case STAB_SUPG: // SUPG 575c1a52365SJed Brown Tau_diagPrim(Tau_d, dXdx, s.Y.velocity, cv, context, mu, dt, s.U.density); 576bb8a0c61SJames Wright tau_strong_res[0] = Tau_d[0] * strong_res[0]; 577bb8a0c61SJames Wright tau_strong_res[1] = Tau_d[1] * strong_res[1]; 578bb8a0c61SJames Wright tau_strong_res[2] = Tau_d[1] * strong_res[2]; 579bb8a0c61SJames Wright tau_strong_res[3] = Tau_d[1] * strong_res[3]; 580bb8a0c61SJames Wright tau_strong_res[4] = Tau_d[2] * strong_res[4]; 581bb8a0c61SJames Wright // Alternate route (useful later with primitive variable code) 582bb8a0c61SJames Wright // this function was verified against PHASTA for as IC that was as close as possible 583bb8a0c61SJames Wright // computeFluxJacobian_NSp(jacob_F_conv_p, rho, u, E, Rd, cv); 584bb8a0c61SJames Wright // it has also been verified to compute a correct through the following 585bb8a0c61SJames Wright // stab[k][j] += jacob_F_conv_p[j][k][l] * tau_strong_res[l] // flux Jacobian wrt primitive 586bb8a0c61SJames Wright // applied in the triple loop below 587bb8a0c61SJames Wright // However, it is more flops than using the existing Jacobian wrt q after q_{,Y} viz 588c1a52365SJed Brown PrimitiveToConservative_fwd(s.U.density, s.Y.velocity, s.U.E_total, Rd, cv, 589c1a52365SJed Brown tau_strong_res, tau_strong_res_conservative); 590493642f1SJames Wright for (CeedInt j=0; j<3; j++) 591493642f1SJames Wright for (CeedInt k=0; k<5; k++) 592493642f1SJames Wright for (CeedInt l=0; l<5; l++) 593bb8a0c61SJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_res_conservative[l]; 5943a8779fbSJames Wright 595493642f1SJames Wright for (CeedInt j=0; j<5; j++) 596493642f1SJames Wright for (CeedInt k=0; k<3; k++) 597752f40e3SJed Brown Grad_v[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 5983a8779fbSJames Wright stab[j][1] * dXdx[k][1] + 5993a8779fbSJames Wright stab[j][2] * dXdx[k][2]); 6003a8779fbSJames Wright break; 6013a8779fbSJames Wright } 602eef2387dSJed Brown for (CeedInt j=0; j<5; j++) jac_data[j][i] = U[j]; 603eef2387dSJed Brown for (CeedInt j=0; j<6; j++) jac_data[5+j][i] = kmstress[j]; 604eef2387dSJed Brown for (CeedInt j=0; j<3; j++) jac_data[5+6+j][i] = Tau_d[j]; 6053a8779fbSJames Wright 6063a8779fbSJames Wright } // End Quadrature Point Loop 6073a8779fbSJames Wright 6083a8779fbSJames Wright // Return 6093a8779fbSJames Wright return 0; 6103a8779fbSJames Wright } 611f0b65372SJed Brown 612f0b65372SJed Brown CEED_QFUNCTION(IJacobian_Newtonian)(void *ctx, CeedInt Q, 613f0b65372SJed Brown const CeedScalar *const *in, 614f0b65372SJed Brown CeedScalar *const *out) { 615f0b65372SJed Brown // *INDENT-OFF* 616f0b65372SJed Brown // Inputs 617f0b65372SJed Brown const CeedScalar (*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 618f0b65372SJed Brown (*Grad_dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 619f0b65372SJed Brown (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 620f0b65372SJed Brown (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3], 621f0b65372SJed Brown (*jac_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 622f0b65372SJed Brown // Outputs 623f0b65372SJed Brown CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 624f0b65372SJed Brown (*Grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 625f0b65372SJed Brown // *INDENT-ON* 626f0b65372SJed Brown // Context 627f0b65372SJed Brown NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx; 628f0b65372SJed Brown const CeedScalar *g = context->g; 629f0b65372SJed Brown const CeedScalar cp = context->cp; 630f0b65372SJed Brown const CeedScalar cv = context->cv; 631f0b65372SJed Brown const CeedScalar Rd = cp - cv; 632f0b65372SJed Brown const CeedScalar gamma = cp / cv; 633f0b65372SJed Brown 634f0b65372SJed Brown CeedPragmaSIMD 635f0b65372SJed Brown // Quadrature Point Loop 636f0b65372SJed Brown for (CeedInt i=0; i<Q; i++) { 637f0b65372SJed Brown // -- Interp-to-Interp q_data 638f0b65372SJed Brown const CeedScalar wdetJ = q_data[0][i]; 639f0b65372SJed Brown // -- Interp-to-Grad q_data 640f0b65372SJed Brown // ---- Inverse of change of coordinate matrix: X_i,j 641f0b65372SJed Brown // *INDENT-OFF* 642f0b65372SJed Brown const CeedScalar dXdx[3][3] = {{q_data[1][i], 643f0b65372SJed Brown q_data[2][i], 644f0b65372SJed Brown q_data[3][i]}, 645f0b65372SJed Brown {q_data[4][i], 646f0b65372SJed Brown q_data[5][i], 647f0b65372SJed Brown q_data[6][i]}, 648f0b65372SJed Brown {q_data[7][i], 649f0b65372SJed Brown q_data[8][i], 650f0b65372SJed Brown q_data[9][i]} 651f0b65372SJed Brown }; 652f0b65372SJed Brown // *INDENT-ON* 653f0b65372SJed Brown 654f0b65372SJed Brown CeedScalar U[5], kmstress[6], Tau_d[3] __attribute((unused)); 655f0b65372SJed Brown for (int j=0; j<5; j++) U[j] = jac_data[j][i]; 656f0b65372SJed Brown for (int j=0; j<6; j++) kmstress[j] = jac_data[5+j][i]; 657f0b65372SJed Brown for (int j=0; j<3; j++) Tau_d[j] = jac_data[5+6+j][i]; 658f0b65372SJed Brown const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 659f0b65372SJed Brown State s = StateFromU(context, U, x_i); 660f0b65372SJed Brown 661f0b65372SJed Brown CeedScalar dU[5], dx0[3] = {0}; 662f0b65372SJed Brown for (int j=0; j<5; j++) dU[j] = dq[j][i]; 663f0b65372SJed Brown State ds = StateFromU_fwd(context, s, dU, x_i, dx0); 664f0b65372SJed Brown 665f0b65372SJed Brown State grad_ds[3]; 666f0b65372SJed Brown for (int j=0; j<3; j++) { 667f0b65372SJed Brown CeedScalar dUj[5]; 668f0b65372SJed Brown for (int k=0; k<5; k++) dUj[k] = Grad_dq[0][k][i] * dXdx[0][j] 669f0b65372SJed Brown + Grad_dq[1][k][i] * dXdx[1][j] 670f0b65372SJed Brown + Grad_dq[2][k][i] * dXdx[2][j]; 671f0b65372SJed Brown grad_ds[j] = StateFromU_fwd(context, s, dUj, x_i, dx0); 672f0b65372SJed Brown } 673f0b65372SJed Brown 674f0b65372SJed Brown CeedScalar dstrain_rate[6], dkmstress[6], stress[3][3], dstress[3][3], dFe[3]; 675f0b65372SJed Brown KMStrainRate(grad_ds, dstrain_rate); 676f0b65372SJed Brown NewtonianStress(context, dstrain_rate, dkmstress); 677f0b65372SJed Brown KMUnpack(dkmstress, dstress); 678f0b65372SJed Brown KMUnpack(kmstress, stress); 679f0b65372SJed Brown ViscousEnergyFlux_fwd(context, s.Y, ds.Y, grad_ds, stress, dstress, dFe); 680f0b65372SJed Brown 681f0b65372SJed Brown StateConservative dF_inviscid[3]; 682f0b65372SJed Brown FluxInviscid_fwd(context, s, ds, dF_inviscid); 683f0b65372SJed Brown 684f0b65372SJed Brown // Total flux 685f0b65372SJed Brown CeedScalar dFlux[5][3]; 686f0b65372SJed Brown for (int j=0; j<3; j++) { 687f0b65372SJed Brown dFlux[0][j] = dF_inviscid[j].density; 688f0b65372SJed Brown for (int k=0; k<3; k++) 689f0b65372SJed Brown dFlux[k+1][j] = dF_inviscid[j].momentum[k] - dstress[k][j]; 690f0b65372SJed Brown dFlux[4][j] = dF_inviscid[j].E_total + dFe[j]; 691f0b65372SJed Brown } 692f0b65372SJed Brown 693f0b65372SJed Brown for (int j=0; j<3; j++) { 694f0b65372SJed Brown for (int k=0; k<5; k++) { 695f0b65372SJed Brown Grad_v[j][k][i] = -wdetJ * (dXdx[j][0] * dFlux[k][0] + 696f0b65372SJed Brown dXdx[j][1] * dFlux[k][1] + 697f0b65372SJed Brown dXdx[j][2] * dFlux[k][2]); 698f0b65372SJed Brown } 699f0b65372SJed Brown } 700f0b65372SJed Brown 701f0b65372SJed Brown const CeedScalar dbody_force[5] = {0, ds.U.density *g[0], ds.U.density *g[1], ds.U.density *g[2], 0}; 702f0b65372SJed Brown for (int j=0; j<5; j++) 703f0b65372SJed Brown v[j][i] = wdetJ * (context->ijacobian_time_shift * dU[j] - dbody_force[j]); 704f0b65372SJed Brown 705f0b65372SJed Brown if (1) { 706f0b65372SJed Brown CeedScalar jacob_F_conv[3][5][5] = {0}; 707f0b65372SJed Brown computeFluxJacobian_NS(jacob_F_conv, s.U.density, s.Y.velocity, s.U.E_total, 708f0b65372SJed Brown gamma, g, x_i); 709f0b65372SJed Brown CeedScalar grad_dU[5][3]; 710f0b65372SJed Brown for (int j=0; j<3; j++) { 711f0b65372SJed Brown grad_dU[0][j] = grad_ds[j].U.density; 712f0b65372SJed Brown for (int k=0; k<3; k++) grad_dU[k+1][j] = grad_ds[j].U.momentum[k]; 713f0b65372SJed Brown grad_dU[4][j] = grad_ds[j].U.E_total; 714f0b65372SJed Brown } 715f0b65372SJed Brown CeedScalar dstrong_conv[5] = {0}; 716f0b65372SJed Brown for (int j=0; j<3; j++) 717f0b65372SJed Brown for (int k=0; k<5; k++) 718f0b65372SJed Brown for (int l=0; l<5; l++) 719f0b65372SJed Brown dstrong_conv[k] += jacob_F_conv[j][k][l] * grad_dU[l][j]; 720f0b65372SJed Brown CeedScalar dstrong_res[5]; 721f0b65372SJed Brown for (int j=0; j<5; j++) 722f0b65372SJed Brown dstrong_res[j] = context->ijacobian_time_shift * dU[j] + dstrong_conv[j] - 723f0b65372SJed Brown dbody_force[j]; 724f0b65372SJed Brown CeedScalar dtau_strong_res[5] = {0.}, dtau_strong_res_conservative[5] = {0}; 725f0b65372SJed Brown dtau_strong_res[0] = Tau_d[0] * dstrong_res[0]; 726f0b65372SJed Brown dtau_strong_res[1] = Tau_d[1] * dstrong_res[1]; 727f0b65372SJed Brown dtau_strong_res[2] = Tau_d[1] * dstrong_res[2]; 728f0b65372SJed Brown dtau_strong_res[3] = Tau_d[1] * dstrong_res[3]; 729f0b65372SJed Brown dtau_strong_res[4] = Tau_d[2] * dstrong_res[4]; 730f0b65372SJed Brown PrimitiveToConservative_fwd(s.U.density, s.Y.velocity, s.U.E_total, Rd, cv, 731f0b65372SJed Brown dtau_strong_res, dtau_strong_res_conservative); 732f0b65372SJed Brown CeedScalar dstab[5][3] = {0}; 733f0b65372SJed Brown for (int j=0; j<3; j++) 734f0b65372SJed Brown for (int k=0; k<5; k++) 735f0b65372SJed Brown for (int l=0; l<5; l++) 736f0b65372SJed Brown dstab[k][j] += jacob_F_conv[j][k][l] * dtau_strong_res_conservative[l]; 737f0b65372SJed Brown for (int j=0; j<5; j++) 738f0b65372SJed Brown for (int k=0; k<3; k++) 739f0b65372SJed Brown Grad_v[k][j][i] += wdetJ*(dstab[j][0] * dXdx[k][0] + 740f0b65372SJed Brown dstab[j][1] * dXdx[k][1] + 741f0b65372SJed Brown dstab[j][2] * dXdx[k][2]); 742f0b65372SJed Brown 743f0b65372SJed Brown } 744f0b65372SJed Brown } // End Quadrature Point Loop 745f0b65372SJed Brown return 0; 746f0b65372SJed Brown } 747*8085925cSJames Wright 748*8085925cSJames Wright // Compute boundary integral (ie. for strongly set inflows) 749*8085925cSJames Wright CEED_QFUNCTION(BoundaryIntegral)(void *ctx, CeedInt Q, 750*8085925cSJames Wright const CeedScalar *const *in, 751*8085925cSJames Wright CeedScalar *const *out) { 752*8085925cSJames Wright 753*8085925cSJames Wright //*INDENT-OFF* 754*8085925cSJames Wright const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0], 755*8085925cSJames Wright (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[2]; 756*8085925cSJames Wright 757*8085925cSJames Wright CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[0]; 758*8085925cSJames Wright 759*8085925cSJames Wright //*INDENT-ON* 760*8085925cSJames Wright 761*8085925cSJames Wright const NewtonianIdealGasContext newt_ctx = (NewtonianIdealGasContext) ctx; 762*8085925cSJames Wright const bool is_implicit = newt_ctx->is_implicit; 763*8085925cSJames Wright const CeedScalar cv = newt_ctx->cv; 764*8085925cSJames Wright const CeedScalar cp = newt_ctx->cp; 765*8085925cSJames Wright const CeedScalar gamma = cp/cv; 766*8085925cSJames Wright 767*8085925cSJames Wright CeedPragmaSIMD 768*8085925cSJames Wright for(CeedInt i=0; i<Q; i++) { 769*8085925cSJames Wright const CeedScalar rho = q[0][i]; 770*8085925cSJames Wright const CeedScalar u[] = {q[1][i]/rho, q[2][i]/rho, q[3][i]/rho}; 771*8085925cSJames Wright const CeedScalar E_kinetic = .5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]); 772*8085925cSJames Wright const CeedScalar E_internal = q[4][i] - E_kinetic; 773*8085925cSJames Wright const CeedScalar P = E_internal * (gamma - 1.); 774*8085925cSJames Wright 775*8085925cSJames Wright const CeedScalar wdetJb = (is_implicit ? -1. : 1.) * q_data_sur[0][i]; 776*8085925cSJames Wright // ---- Normal vect 777*8085925cSJames Wright const CeedScalar norm[3] = {q_data_sur[1][i], 778*8085925cSJames Wright q_data_sur[2][i], 779*8085925cSJames Wright q_data_sur[3][i] 780*8085925cSJames Wright }; 781*8085925cSJames Wright 782*8085925cSJames Wright const CeedScalar E = E_internal + E_kinetic; 783*8085925cSJames Wright 784*8085925cSJames Wright // Velocity normal to the boundary 785*8085925cSJames Wright const CeedScalar u_normal = norm[0]*u[0] + 786*8085925cSJames Wright norm[1]*u[1] + 787*8085925cSJames Wright norm[2]*u[2]; 788*8085925cSJames Wright // The Physics 789*8085925cSJames Wright // Zero v so all future terms can safely sum into it 790*8085925cSJames Wright for (CeedInt j=0; j<5; j++) v[j][i] = 0.; 791*8085925cSJames Wright 792*8085925cSJames Wright // The Physics 793*8085925cSJames Wright // -- Density 794*8085925cSJames Wright v[0][i] -= wdetJb * rho * u_normal; 795*8085925cSJames Wright 796*8085925cSJames Wright // -- Momentum 797*8085925cSJames Wright for (CeedInt j=0; j<3; j++) 798*8085925cSJames Wright v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + 799*8085925cSJames Wright norm[j] * P); 800*8085925cSJames Wright 801*8085925cSJames Wright // -- Total Energy Density 802*8085925cSJames Wright v[4][i] -= wdetJb * u_normal * (E + P); 803*8085925cSJames Wright } 804*8085925cSJames Wright return 0; 805*8085925cSJames Wright } 806*8085925cSJames Wright 8073a8779fbSJames Wright // ***************************************************************************** 8083a8779fbSJames Wright #endif // newtonian_h 809