13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 388b783a1SJames Wright // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 588b783a1SJames Wright // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 788b783a1SJames Wright 888b783a1SJames Wright /// @file 988b783a1SJames Wright /// Operator for Navier-Stokes example using PETSc 1088b783a1SJames Wright 1188b783a1SJames Wright 1288b783a1SJames Wright #ifndef newtonian_h 1388b783a1SJames Wright #define newtonian_h 1488b783a1SJames Wright 1588b783a1SJames Wright #include <math.h> 1688b783a1SJames Wright #include <ceed.h> 1788b783a1SJames Wright 1888b783a1SJames Wright #ifndef M_PI 1988b783a1SJames Wright #define M_PI 3.14159265358979323846 2088b783a1SJames Wright #endif 2188b783a1SJames Wright 2288b783a1SJames Wright #ifndef setup_context_struct 2388b783a1SJames Wright #define setup_context_struct 2488b783a1SJames Wright typedef struct SetupContext_ *SetupContext; 2588b783a1SJames Wright struct SetupContext_ { 2688b783a1SJames Wright CeedScalar theta0; 2788b783a1SJames Wright CeedScalar thetaC; 2888b783a1SJames Wright CeedScalar P0; 2988b783a1SJames Wright CeedScalar N; 3088b783a1SJames Wright CeedScalar cv; 3188b783a1SJames Wright CeedScalar cp; 3288626eedSJames Wright CeedScalar g[3]; 3388b783a1SJames Wright CeedScalar rc; 3488b783a1SJames Wright CeedScalar lx; 3588b783a1SJames Wright CeedScalar ly; 3688b783a1SJames Wright CeedScalar lz; 3788b783a1SJames Wright CeedScalar center[3]; 3888b783a1SJames Wright CeedScalar dc_axis[3]; 3988b783a1SJames Wright CeedScalar wind[3]; 4088b783a1SJames Wright CeedScalar time; 4188b783a1SJames Wright int wind_type; // See WindType: 0=ROTATION, 1=TRANSLATION 4288b783a1SJames Wright int bubble_type; // See BubbleType: 0=SPHERE, 1=CYLINDER 4388b783a1SJames Wright int bubble_continuity_type; // See BubbleContinuityType: 0=SMOOTH, 1=BACK_SHARP 2=THICK 4488b783a1SJames Wright }; 4588b783a1SJames Wright #endif 4688b783a1SJames Wright 4788b783a1SJames Wright #ifndef newtonian_context_struct 4888b783a1SJames Wright #define newtonian_context_struct 4988b783a1SJames Wright typedef enum { 5088b783a1SJames Wright STAB_NONE = 0, 5188b783a1SJames Wright STAB_SU = 1, // Streamline Upwind 5288b783a1SJames Wright STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin 5388b783a1SJames Wright } StabilizationType; 5488b783a1SJames Wright 5588b783a1SJames Wright typedef struct NewtonianIdealGasContext_ *NewtonianIdealGasContext; 5688b783a1SJames Wright struct NewtonianIdealGasContext_ { 5788b783a1SJames Wright CeedScalar lambda; 5888b783a1SJames Wright CeedScalar mu; 5988b783a1SJames Wright CeedScalar k; 6088b783a1SJames Wright CeedScalar cv; 6188b783a1SJames Wright CeedScalar cp; 6288626eedSJames Wright CeedScalar g[3]; 6388b783a1SJames Wright CeedScalar c_tau; 6488626eedSJames Wright CeedScalar Ctau_t; 6588626eedSJames Wright CeedScalar Ctau_v; 6688626eedSJames Wright CeedScalar Ctau_C; 6788626eedSJames Wright CeedScalar Ctau_M; 6888626eedSJames Wright CeedScalar Ctau_E; 6988626eedSJames Wright CeedScalar dt; 7088b783a1SJames Wright StabilizationType stabilization; 7188b783a1SJames Wright }; 7288b783a1SJames Wright #endif 7388b783a1SJames Wright 7488b783a1SJames Wright // ***************************************************************************** 7588b783a1SJames Wright // Helper function for computing flux Jacobian 7688b783a1SJames Wright // ***************************************************************************** 7788b783a1SJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NS(CeedScalar dF[3][5][5], 7888b783a1SJames Wright const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 7988626eedSJames Wright const CeedScalar gamma, const CeedScalar g[3], const CeedScalar x[3]) { 8088b783a1SJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 8188626eedSJames Wright CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]); 8288b783a1SJames Wright for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 8388b783a1SJames Wright for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix 8488626eedSJames Wright dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2. - e_potential)) : 0.) - 8588626eedSJames Wright u[i]*u[j]; 8688b783a1SJames Wright for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix 8788b783a1SJames Wright dF[i][0][k+1] = ((i==k) ? 1. : 0.); 8888b783a1SJames Wright dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) + 8988b783a1SJames Wright ((i==k) ? u[j] : 0.) - 9088b783a1SJames Wright ((i==j) ? u[k] : 0.) * (gamma-1.); 9188b783a1SJames Wright dF[i][4][k+1] = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) - 9288b783a1SJames Wright (gamma-1.)*u[i]*u[k]; 9388b783a1SJames Wright } 9488b783a1SJames Wright dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.); 9588b783a1SJames Wright } 9688b783a1SJames Wright dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho); 9788b783a1SJames Wright dF[i][4][4] = u[i] * gamma; 9888b783a1SJames Wright } 9988b783a1SJames Wright } 10088b783a1SJames Wright 10188b783a1SJames Wright // ***************************************************************************** 10288626eedSJames Wright // Helper function for computing flux Jacobian of Primitive variables 10388626eedSJames Wright // ***************************************************************************** 10488626eedSJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NSp(CeedScalar dF[3][5][5], 10588626eedSJames Wright const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 10688626eedSJames Wright const CeedScalar Rd, const CeedScalar cv) { 10788626eedSJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 10888626eedSJames Wright // TODO Add in gravity's contribution 10988626eedSJames Wright 11088626eedSJames Wright CeedScalar T = ( E / rho - u_sq / 2. ) / cv; 11188626eedSJames Wright CeedScalar drdT = -rho / T; 11288626eedSJames Wright CeedScalar drdP = 1. / ( Rd * T); 11388626eedSJames Wright CeedScalar etot = E / rho ; 11488626eedSJames Wright CeedScalar e2p = drdP * etot + 1. ; 11588626eedSJames Wright CeedScalar e3p = ( E + rho * Rd * T ); 11688626eedSJames Wright CeedScalar e4p = drdT * etot + rho * cv ; 11788626eedSJames Wright 11888626eedSJames Wright for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 11988626eedSJames Wright for (CeedInt j=0; j<3; j++) { // j counts F^{m_j} 12088626eedSJames Wright // [row][col] of A_i 12188626eedSJames Wright dF[i][j+1][0] = drdP * u[i] * u[j] + ((i==j) ? 1. : 0.); // F^{{m_j} wrt p 12288626eedSJames Wright for (CeedInt k=0; k<3; k++) { // k counts the wrt vel_k 123*871db79fSKenneth E. Jansen dF[i][0][k+1] = ((i==k) ? rho : 0.); // F^c wrt u_k 12488626eedSJames Wright dF[i][j+1][k+1] = (((j==k) ? u[i] : 0.) + // F^m_j wrt u_k 12588626eedSJames Wright ((i==k) ? u[j] : 0.) ) * rho; 12688626eedSJames Wright dF[i][4][k+1] = rho * u[i] * u[k] 12788626eedSJames Wright + ((i==k) ? e3p : 0.) ; // F^e wrt u_k 12888626eedSJames Wright } 12988626eedSJames Wright dF[i][j+1][4] = drdT * u[i] * u[j]; // F^{m_j} wrt T 13088626eedSJames Wright } 13188626eedSJames Wright dF[i][4][0] = u[i] * e2p; // F^e wrt p 13288626eedSJames Wright dF[i][4][4] = u[i] * e4p; // F^e wrt T 13388626eedSJames Wright dF[i][0][0] = u[i] * drdP; // F^c wrt p 13488626eedSJames Wright dF[i][0][4] = u[i] * drdT; // F^c wrt T 13588626eedSJames Wright } 13688626eedSJames Wright } 13788626eedSJames Wright 13888626eedSJames Wright CEED_QFUNCTION_HELPER void PrimitiveToConservative_fwd(const CeedScalar rho, 13988626eedSJames Wright const CeedScalar u[3], const CeedScalar E, const CeedScalar Rd, 14088626eedSJames Wright const CeedScalar cv, const CeedScalar dY[5], CeedScalar dU[5]) { 14188626eedSJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; 14288626eedSJames Wright CeedScalar T = ( E / rho - u_sq / 2. ) / cv; 14388626eedSJames Wright CeedScalar drdT = -rho / T; 14488626eedSJames Wright CeedScalar drdP = 1. / ( Rd * T); 14588626eedSJames Wright dU[0] = drdP * dY[0] + drdT * dY[4]; 14688626eedSJames Wright CeedScalar de_kinetic = 0; 14788626eedSJames Wright for (int i=0; i<3; i++) { 14888626eedSJames Wright dU[1+i] = dU[0] * u[i] + rho * dY[1+i]; 14988626eedSJames Wright de_kinetic += u[i] * dY[1+i]; 15088626eedSJames Wright } 15188626eedSJames Wright dU[4] = rho * cv * dY[4] + dU[0] * cv * T // internal energy: rho * e 15288626eedSJames Wright + rho * de_kinetic + .5 * dU[0] * u_sq; // kinetic energy: .5 * rho * |u|^2 15388626eedSJames Wright } 15488626eedSJames Wright 15588626eedSJames Wright // ***************************************************************************** 15688626eedSJames Wright // Helper function for computing Tau elements (stabilization constant) 15788626eedSJames Wright // Model from: 15888626eedSJames Wright // PHASTA 15988626eedSJames Wright // 16088626eedSJames Wright // Tau[i] = itau=0 which is diagonal-Shakib (3 values still but not spatial) 16188626eedSJames Wright // 16288626eedSJames Wright // Where NOT UPDATED YET 16388626eedSJames Wright // ***************************************************************************** 16488626eedSJames Wright CEED_QFUNCTION_HELPER void Tau_diagPrim(CeedScalar Tau_d[3], 16588626eedSJames Wright const CeedScalar dXdx[3][3], const CeedScalar u[3], 16688626eedSJames Wright const CeedScalar cv, const NewtonianIdealGasContext newt_ctx, 16788626eedSJames Wright const CeedScalar mu, const CeedScalar dt, 16888626eedSJames Wright const CeedScalar rho) { 16988626eedSJames Wright // Context 17088626eedSJames Wright const CeedScalar Ctau_t = newt_ctx->Ctau_t; 17188626eedSJames Wright const CeedScalar Ctau_v = newt_ctx->Ctau_v; 17288626eedSJames Wright const CeedScalar Ctau_C = newt_ctx->Ctau_C; 17388626eedSJames Wright const CeedScalar Ctau_M = newt_ctx->Ctau_M; 17488626eedSJames Wright const CeedScalar Ctau_E = newt_ctx->Ctau_E; 17588626eedSJames Wright CeedScalar gijd[6]; 17688626eedSJames Wright CeedScalar tau; 17788626eedSJames Wright CeedScalar dts; 17888626eedSJames Wright CeedScalar fact; 17988626eedSJames Wright 18088626eedSJames Wright //*INDENT-OFF* 18188626eedSJames Wright gijd[0] = dXdx[0][0] * dXdx[0][0] 18288626eedSJames Wright + dXdx[1][0] * dXdx[1][0] 18388626eedSJames Wright + dXdx[2][0] * dXdx[2][0]; 18488626eedSJames Wright 18588626eedSJames Wright gijd[1] = dXdx[0][0] * dXdx[0][1] 18688626eedSJames Wright + dXdx[1][0] * dXdx[1][1] 18788626eedSJames Wright + dXdx[2][0] * dXdx[2][1]; 18888626eedSJames Wright 18988626eedSJames Wright gijd[2] = dXdx[0][1] * dXdx[0][1] 19088626eedSJames Wright + dXdx[1][1] * dXdx[1][1] 19188626eedSJames Wright + dXdx[2][1] * dXdx[2][1]; 19288626eedSJames Wright 19388626eedSJames Wright gijd[3] = dXdx[0][0] * dXdx[0][2] 19488626eedSJames Wright + dXdx[1][0] * dXdx[1][2] 19588626eedSJames Wright + dXdx[2][0] * dXdx[2][2]; 19688626eedSJames Wright 19788626eedSJames Wright gijd[4] = dXdx[0][1] * dXdx[0][2] 19888626eedSJames Wright + dXdx[1][1] * dXdx[1][2] 19988626eedSJames Wright + dXdx[2][1] * dXdx[2][2]; 20088626eedSJames Wright 20188626eedSJames Wright gijd[5] = dXdx[0][2] * dXdx[0][2] 20288626eedSJames Wright + dXdx[1][2] * dXdx[1][2] 20388626eedSJames Wright + dXdx[2][2] * dXdx[2][2]; 20488626eedSJames Wright //*INDENT-ON* 20588626eedSJames Wright 20688626eedSJames Wright dts = Ctau_t / dt ; 20788626eedSJames Wright 20888626eedSJames Wright tau = rho*rho*((4. * dts * dts) 20988626eedSJames Wright + u[0] * ( u[0] * gijd[0] + 2. * ( u[1] * gijd[1] + u[2] * gijd[3])) 21088626eedSJames Wright + u[1] * ( u[1] * gijd[2] + 2. * u[2] * gijd[4]) 21188626eedSJames Wright + u[2] * u[2] * gijd[5]) 21288626eedSJames Wright + Ctau_v* mu * mu * 21388626eedSJames Wright (gijd[0]*gijd[0] + gijd[2]*gijd[2] + gijd[5]*gijd[5] + 21488626eedSJames Wright + 2. * (gijd[1]*gijd[1] + gijd[3]*gijd[3] + gijd[4]*gijd[4])); 21588626eedSJames Wright 21688626eedSJames Wright fact=sqrt(tau); 21788626eedSJames Wright 21888626eedSJames Wright Tau_d[0] = Ctau_C * fact / (rho*(gijd[0] + gijd[2] + gijd[5]))*0.125; 21988626eedSJames Wright 22088626eedSJames Wright Tau_d[1] = Ctau_M / fact; 22188626eedSJames Wright Tau_d[2] = Ctau_E / ( fact * cv ); 22288626eedSJames Wright 22388626eedSJames Wright // consider putting back the way I initially had it Ctau_E * Tau_d[1] /cv 22488626eedSJames Wright // to avoid a division if the compiler is smart enough to see that cv IS 22588626eedSJames Wright // a constant that it could invert once for all elements 22688626eedSJames Wright // but in that case energy tau is scaled by the product of Ctau_E * Ctau_M 22788626eedSJames Wright // OR we could absorb cv into Ctau_E but this puts more burden on user to 22888626eedSJames Wright // know how to change constants with a change of fluid or units. Same for 22988626eedSJames Wright // Ctau_v * mu * mu IF AND ONLY IF we don't add viscosity law =f(T) 23088626eedSJames Wright } 23188626eedSJames Wright 23288626eedSJames Wright // ***************************************************************************** 23388b783a1SJames Wright // Helper function for computing Tau elements (stabilization constant) 23488b783a1SJames Wright // Model from: 23588b783a1SJames Wright // Stabilized Methods for Compressible Flows, Hughes et al 2010 23688b783a1SJames Wright // 23788b783a1SJames Wright // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 23888b783a1SJames Wright // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 23988b783a1SJames Wright // 24088b783a1SJames Wright // Where 24188b783a1SJames Wright // c_tau = stabilization constant (0.5 is reported as "optimal") 24288b783a1SJames Wright // h[i] = 2 length(dxdX[i]) 24388b783a1SJames Wright // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 24488b783a1SJames Wright // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 24588b783a1SJames Wright // rho(A[i]) = spectral radius of the convective flux Jacobian i, 24688b783a1SJames Wright // wave speed in direction i 24788b783a1SJames Wright // ***************************************************************************** 24888b783a1SJames Wright CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], 24988b783a1SJames Wright const CeedScalar dXdx[3][3], const CeedScalar u[3], 25088626eedSJames Wright /* const CeedScalar sound_speed, const CeedScalar c_tau) { */ 25188626eedSJames Wright const CeedScalar sound_speed, const CeedScalar c_tau, 25288626eedSJames Wright const CeedScalar viscosity) { 25388626eedSJames Wright const CeedScalar mag_u_visc = sqrt(u[0]*u[0] +u[1]*u[1] +u[2]*u[2]) / 25488626eedSJames Wright (2*viscosity); 25588b783a1SJames Wright for (int i=0; i<3; i++) { 25688b783a1SJames Wright // length of element in direction i 25788b783a1SJames Wright CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] + 25888b783a1SJames Wright dXdx[2][i]*dXdx[2][i]); 25988626eedSJames Wright CeedScalar Pe = mag_u_visc*h; 26088626eedSJames Wright CeedScalar Xi = 1/tanh(Pe) - 1/Pe; 26188b783a1SJames Wright // fastest wave in direction i 26288b783a1SJames Wright CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 26388626eedSJames Wright Tau_x[i] = c_tau * h * Xi / fastest_wave; 26488b783a1SJames Wright } 26588b783a1SJames Wright } 26688b783a1SJames Wright 26788b783a1SJames Wright // ***************************************************************************** 26888b783a1SJames Wright // This QFunction sets a "still" initial condition for generic Newtonian IG problems 26988b783a1SJames Wright // ***************************************************************************** 27088b783a1SJames Wright CEED_QFUNCTION(ICsNewtonianIG)(void *ctx, CeedInt Q, 27188b783a1SJames Wright const CeedScalar *const *in, CeedScalar *const *out) { 27288b783a1SJames Wright // Inputs 27388b783a1SJames Wright const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 27488b783a1SJames Wright 27588b783a1SJames Wright // Outputs 27688b783a1SJames Wright CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 27788b783a1SJames Wright 27888626eedSJames Wright // Context 27988626eedSJames Wright const SetupContext context = (SetupContext)ctx; 28088626eedSJames Wright const CeedScalar theta0 = context->theta0; 28188626eedSJames Wright const CeedScalar P0 = context->P0; 28288626eedSJames Wright const CeedScalar cv = context->cv; 28388626eedSJames Wright const CeedScalar cp = context->cp; 28488626eedSJames Wright const CeedScalar *g = context->g; 28588626eedSJames Wright const CeedScalar Rd = cp - cv; 28688626eedSJames Wright 28788b783a1SJames Wright // Quadrature Point Loop 28888b783a1SJames Wright CeedPragmaSIMD 28988b783a1SJames Wright for (CeedInt i=0; i<Q; i++) { 29088b783a1SJames Wright CeedScalar q[5] = {0.}; 29188b783a1SJames Wright 29288b783a1SJames Wright // Setup 29388b783a1SJames Wright // -- Coordinates 29488626eedSJames Wright const CeedScalar x[3] = {X[0][i], X[1][i], X[2][i]}; 29588626eedSJames Wright const CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]); 29688b783a1SJames Wright 29788b783a1SJames Wright // -- Density 29888626eedSJames Wright const CeedScalar rho = P0 / (Rd*theta0); 29988b783a1SJames Wright 30088b783a1SJames Wright // Initial Conditions 30188b783a1SJames Wright q[0] = rho; 30288b783a1SJames Wright q[1] = 0.0; 30388b783a1SJames Wright q[2] = 0.0; 30488b783a1SJames Wright q[3] = 0.0; 30588626eedSJames Wright q[4] = rho * (cv*theta0 + e_potential); 30688b783a1SJames Wright 30788b783a1SJames Wright for (CeedInt j=0; j<5; j++) 30888b783a1SJames Wright q0[j][i] = q[j]; 30988b783a1SJames Wright } // End of Quadrature Point Loop 31088b783a1SJames Wright return 0; 31188b783a1SJames Wright } 31288b783a1SJames Wright 31388b783a1SJames Wright // ***************************************************************************** 31488b783a1SJames Wright // This QFunction implements the following formulation of Navier-Stokes with 31588b783a1SJames Wright // explicit time stepping method 31688b783a1SJames Wright // 31788b783a1SJames Wright // This is 3D compressible Navier-Stokes in conservation form with state 31888b783a1SJames Wright // variables of density, momentum density, and total energy density. 31988b783a1SJames Wright // 32088b783a1SJames Wright // State Variables: q = ( rho, U1, U2, U3, E ) 32188b783a1SJames Wright // rho - Mass Density 32288b783a1SJames Wright // Ui - Momentum Density, Ui = rho ui 32388b783a1SJames Wright // E - Total Energy Density, E = rho (cv T + (u u)/2 + g z) 32488b783a1SJames Wright // 32588b783a1SJames Wright // Navier-Stokes Equations: 32688b783a1SJames Wright // drho/dt + div( U ) = 0 32788b783a1SJames Wright // dU/dt + div( rho (u x u) + P I3 ) + rho g khat = div( Fu ) 32888b783a1SJames Wright // dE/dt + div( (E + P) u ) = div( Fe ) 32988b783a1SJames Wright // 33088b783a1SJames Wright // Viscous Stress: 33188b783a1SJames Wright // Fu = mu (grad( u ) + grad( u )^T + lambda div ( u ) I3) 33288b783a1SJames Wright // 33388b783a1SJames Wright // Thermal Stress: 33488b783a1SJames Wright // Fe = u Fu + k grad( T ) 33588626eedSJames Wright // Equation of State 33688b783a1SJames Wright // P = (gamma - 1) (E - rho (u u) / 2 - rho g z) 33788b783a1SJames Wright // 33888b783a1SJames Wright // Stabilization: 33988b783a1SJames Wright // Tau = diag(TauC, TauM, TauM, TauM, TauE) 34088b783a1SJames Wright // f1 = rho sqrt(ui uj gij) 34188b783a1SJames Wright // gij = dXi/dX * dXi/dX 34288b783a1SJames Wright // TauC = Cc f1 / (8 gii) 34388b783a1SJames Wright // TauM = min( 1 , 1 / f1 ) 34488b783a1SJames Wright // TauE = TauM / (Ce cv) 34588b783a1SJames Wright // 34688b783a1SJames Wright // SU = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) ) 34788b783a1SJames Wright // 34888b783a1SJames Wright // Constants: 34988b783a1SJames Wright // lambda = - 2 / 3, From Stokes hypothesis 35088b783a1SJames Wright // mu , Dynamic viscosity 35188b783a1SJames Wright // k , Thermal conductivity 35288b783a1SJames Wright // cv , Specific heat, constant volume 35388b783a1SJames Wright // cp , Specific heat, constant pressure 35488b783a1SJames Wright // g , Gravity 35588b783a1SJames Wright // gamma = cp / cv, Specific heat ratio 35688b783a1SJames Wright // 35788b783a1SJames Wright // We require the product of the inverse of the Jacobian (dXdx_j,k) and 35888b783a1SJames Wright // its transpose (dXdx_k,j) to properly compute integrals of the form: 35988b783a1SJames Wright // int( gradv gradu ) 36088b783a1SJames Wright // 36188b783a1SJames Wright // ***************************************************************************** 36288b783a1SJames Wright CEED_QFUNCTION(Newtonian)(void *ctx, CeedInt Q, 36388b783a1SJames Wright const CeedScalar *const *in, CeedScalar *const *out) { 36488b783a1SJames Wright // *INDENT-OFF* 36588b783a1SJames Wright // Inputs 36688b783a1SJames Wright const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 36788b783a1SJames Wright (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 36888b783a1SJames Wright (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 36988b783a1SJames Wright (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 37088b783a1SJames Wright // Outputs 37188b783a1SJames Wright CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 37288b783a1SJames Wright (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 37388b783a1SJames Wright // *INDENT-ON* 37488b783a1SJames Wright 37588b783a1SJames Wright // Context 37688b783a1SJames Wright NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx; 37788b783a1SJames Wright const CeedScalar lambda = context->lambda; 37888b783a1SJames Wright const CeedScalar mu = context->mu; 37988b783a1SJames Wright const CeedScalar k = context->k; 38088b783a1SJames Wright const CeedScalar cv = context->cv; 38188b783a1SJames Wright const CeedScalar cp = context->cp; 38288626eedSJames Wright const CeedScalar *g = context->g; 38388626eedSJames Wright const CeedScalar dt = context->dt; 38488b783a1SJames Wright const CeedScalar gamma = cp / cv; 38588626eedSJames Wright const CeedScalar Rd = cp - cv; 38688b783a1SJames Wright 38788b783a1SJames Wright CeedPragmaSIMD 38888b783a1SJames Wright // Quadrature Point Loop 38988b783a1SJames Wright for (CeedInt i=0; i<Q; i++) { 39088b783a1SJames Wright // *INDENT-OFF* 39188b783a1SJames Wright // Setup 39288b783a1SJames Wright // -- Interp in 39388b783a1SJames Wright const CeedScalar rho = q[0][i]; 39488b783a1SJames Wright const CeedScalar u[3] = {q[1][i] / rho, 39588b783a1SJames Wright q[2][i] / rho, 39688b783a1SJames Wright q[3][i] / rho 39788b783a1SJames Wright }; 39888b783a1SJames Wright const CeedScalar E = q[4][i]; 39988b783a1SJames Wright // -- Grad in 40088b783a1SJames Wright const CeedScalar drho[3] = {dq[0][0][i], 40188b783a1SJames Wright dq[1][0][i], 40288b783a1SJames Wright dq[2][0][i] 40388b783a1SJames Wright }; 40488b783a1SJames Wright const CeedScalar dU[3][3] = {{dq[0][1][i], 40588b783a1SJames Wright dq[1][1][i], 40688b783a1SJames Wright dq[2][1][i]}, 40788b783a1SJames Wright {dq[0][2][i], 40888b783a1SJames Wright dq[1][2][i], 40988b783a1SJames Wright dq[2][2][i]}, 41088b783a1SJames Wright {dq[0][3][i], 41188b783a1SJames Wright dq[1][3][i], 41288b783a1SJames Wright dq[2][3][i]} 41388b783a1SJames Wright }; 41488b783a1SJames Wright const CeedScalar dE[3] = {dq[0][4][i], 41588b783a1SJames Wright dq[1][4][i], 41688b783a1SJames Wright dq[2][4][i] 41788b783a1SJames Wright }; 41888b783a1SJames Wright // -- Interp-to-Interp q_data 41988b783a1SJames Wright const CeedScalar wdetJ = q_data[0][i]; 42088b783a1SJames Wright // -- Interp-to-Grad q_data 42188b783a1SJames Wright // ---- Inverse of change of coordinate matrix: X_i,j 42288b783a1SJames Wright // *INDENT-OFF* 42388b783a1SJames Wright const CeedScalar dXdx[3][3] = {{q_data[1][i], 42488b783a1SJames Wright q_data[2][i], 42588b783a1SJames Wright q_data[3][i]}, 42688b783a1SJames Wright {q_data[4][i], 42788b783a1SJames Wright q_data[5][i], 42888b783a1SJames Wright q_data[6][i]}, 42988b783a1SJames Wright {q_data[7][i], 43088b783a1SJames Wright q_data[8][i], 43188b783a1SJames Wright q_data[9][i]} 43288b783a1SJames Wright }; 43388626eedSJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 43488b783a1SJames Wright // *INDENT-ON* 43588b783a1SJames Wright // -- Grad-to-Grad q_data 43688b783a1SJames Wright // dU/dx 43788b783a1SJames Wright CeedScalar du[3][3] = {{0}}; 43888b783a1SJames Wright CeedScalar drhodx[3] = {0}; 43988b783a1SJames Wright CeedScalar dEdx[3] = {0}; 44088b783a1SJames Wright CeedScalar dUdx[3][3] = {{0}}; 44188b783a1SJames Wright CeedScalar dXdxdXdxT[3][3] = {{0}}; 44288b783a1SJames Wright for (int j=0; j<3; j++) { 44388b783a1SJames Wright for (int k=0; k<3; k++) { 44488b783a1SJames Wright du[j][k] = (dU[j][k] - drho[k]*u[j]) / rho; 44588b783a1SJames Wright drhodx[j] += drho[k] * dXdx[k][j]; 44688b783a1SJames Wright dEdx[j] += dE[k] * dXdx[k][j]; 44788b783a1SJames Wright for (int l=0; l<3; l++) { 44888b783a1SJames Wright dUdx[j][k] += dU[j][l] * dXdx[l][k]; 44988b783a1SJames Wright dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 45088b783a1SJames Wright } 45188b783a1SJames Wright } 45288b783a1SJames Wright } 45388b783a1SJames Wright CeedScalar dudx[3][3] = {{0}}; 45488b783a1SJames Wright for (int j=0; j<3; j++) 45588b783a1SJames Wright for (int k=0; k<3; k++) 45688b783a1SJames Wright for (int l=0; l<3; l++) 45788b783a1SJames Wright dudx[j][k] += du[j][l] * dXdx[l][k]; 45888b783a1SJames Wright // -- grad_T 45988b783a1SJames Wright const CeedScalar grad_T[3] = {(dEdx[0]/rho - E*drhodx[0]/(rho*rho) - /* *NOPAD* */ 46088626eedSJames Wright (u[0]*dudx[0][0] + u[1]*dudx[1][0] + u[2]*dudx[2][0]) + g[0])/cv, 46188b783a1SJames Wright (dEdx[1]/rho - E*drhodx[1]/(rho*rho) - /* *NOPAD* */ 46288626eedSJames Wright (u[0]*dudx[0][1] + u[1]*dudx[1][1] + u[2]*dudx[2][1]) + g[1])/cv, 46388b783a1SJames Wright (dEdx[2]/rho - E*drhodx[2]/(rho*rho) - /* *NOPAD* */ 46488626eedSJames Wright (u[0]*dudx[0][2] + u[1]*dudx[1][2] + u[2]*dudx[2][2]) + g[2])/cv 46588b783a1SJames Wright }; 46688b783a1SJames Wright 46788b783a1SJames Wright // -- Fuvisc 46888b783a1SJames Wright // ---- Symmetric 3x3 matrix 46988b783a1SJames Wright const CeedScalar Fu[6] = {mu*(dudx[0][0] * (2 + lambda) + /* *NOPAD* */ 47088b783a1SJames Wright lambda * (dudx[1][1] + dudx[2][2])), 47188b783a1SJames Wright mu*(dudx[0][1] + dudx[1][0]), /* *NOPAD* */ 47288b783a1SJames Wright mu*(dudx[0][2] + dudx[2][0]), /* *NOPAD* */ 47388b783a1SJames Wright mu*(dudx[1][1] * (2 + lambda) + /* *NOPAD* */ 47488b783a1SJames Wright lambda * (dudx[0][0] + dudx[2][2])), 47588b783a1SJames Wright mu*(dudx[1][2] + dudx[2][1]), /* *NOPAD* */ 47688b783a1SJames Wright mu*(dudx[2][2] * (2 + lambda) + /* *NOPAD* */ 47788b783a1SJames Wright lambda * (dudx[0][0] + dudx[1][1])) 47888b783a1SJames Wright }; 47988b783a1SJames Wright // -- Fevisc 48088b783a1SJames Wright const CeedScalar Fe[3] = {u[0]*Fu[0] + u[1]*Fu[1] + u[2]*Fu[2] + /* *NOPAD* */ 48188b783a1SJames Wright k*grad_T[0], /* *NOPAD* */ 48288b783a1SJames Wright u[0]*Fu[1] + u[1]*Fu[3] + u[2]*Fu[4] + /* *NOPAD* */ 48388b783a1SJames Wright k*grad_T[1], /* *NOPAD* */ 48488b783a1SJames Wright u[0]*Fu[2] + u[1]*Fu[4] + u[2]*Fu[5] + /* *NOPAD* */ 48588b783a1SJames Wright k*grad_T[2] /* *NOPAD* */ 48688b783a1SJames Wright }; 48788b783a1SJames Wright // Pressure 48888b783a1SJames Wright const CeedScalar 48988b783a1SJames Wright E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 49088626eedSJames Wright E_potential = -rho*(g[0]*x_i[0] + g[1]*x_i[1] + g[2]*x_i[2]), 49188b783a1SJames Wright E_internal = E - E_kinetic - E_potential, 49288b783a1SJames Wright P = E_internal * (gamma - 1.); // P = pressure 49388b783a1SJames Wright 49488b783a1SJames Wright // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 49588b783a1SJames Wright CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 49688626eedSJames Wright computeFluxJacobian_NS(jacob_F_conv, rho, u, E, gamma, g, x_i); 49788b783a1SJames Wright 49888b783a1SJames Wright // dqdx collects drhodx, dUdx and dEdx in one vector 49988b783a1SJames Wright CeedScalar dqdx[5][3]; 50088b783a1SJames Wright for (int j=0; j<3; j++) { 50188b783a1SJames Wright dqdx[0][j] = drhodx[j]; 50288b783a1SJames Wright dqdx[4][j] = dEdx[j]; 50388b783a1SJames Wright for (int k=0; k<3; k++) 50488b783a1SJames Wright dqdx[k+1][j] = dUdx[k][j]; 50588b783a1SJames Wright } 50688b783a1SJames Wright 50788b783a1SJames Wright // strong_conv = dF/dq * dq/dx (Strong convection) 50888b783a1SJames Wright CeedScalar strong_conv[5] = {0}; 50988b783a1SJames Wright for (int j=0; j<3; j++) 51088b783a1SJames Wright for (int k=0; k<5; k++) 51188b783a1SJames Wright for (int l=0; l<5; l++) 51288b783a1SJames Wright strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 51388b783a1SJames Wright 51488b783a1SJames Wright // Body force 51588626eedSJames Wright const CeedScalar body_force[5] = {0, rho *g[0], rho *g[1], rho *g[2], 0}; 51688b783a1SJames Wright 51788b783a1SJames Wright // The Physics 51888b783a1SJames Wright // Zero dv so all future terms can safely sum into it 51988b783a1SJames Wright for (int j=0; j<5; j++) 52088b783a1SJames Wright for (int k=0; k<3; k++) 52188b783a1SJames Wright dv[k][j][i] = 0; 52288b783a1SJames Wright 52388b783a1SJames Wright // -- Density 52488b783a1SJames Wright // ---- u rho 52588b783a1SJames Wright for (int j=0; j<3; j++) 52688b783a1SJames Wright dv[j][0][i] += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 52788b783a1SJames Wright rho*u[2]*dXdx[j][2]); 52888b783a1SJames Wright // -- Momentum 52988b783a1SJames Wright // ---- rho (u x u) + P I3 53088b783a1SJames Wright for (int j=0; j<3; j++) 53188b783a1SJames Wright for (int k=0; k<3; k++) 53288b783a1SJames Wright dv[k][j+1][i] += wdetJ*((rho*u[j]*u[0] + (j==0?P:0))*dXdx[k][0] + 53388b783a1SJames Wright (rho*u[j]*u[1] + (j==1?P:0))*dXdx[k][1] + 53488b783a1SJames Wright (rho*u[j]*u[2] + (j==2?P:0))*dXdx[k][2]); 53588b783a1SJames Wright // ---- Fuvisc 53688b783a1SJames Wright const CeedInt Fuviscidx[3][3] = {{0, 1, 2}, {1, 3, 4}, {2, 4, 5}}; // symmetric matrix indices 53788b783a1SJames Wright for (int j=0; j<3; j++) 53888b783a1SJames Wright for (int k=0; k<3; k++) 53988b783a1SJames Wright dv[k][j+1][i] -= wdetJ*(Fu[Fuviscidx[j][0]]*dXdx[k][0] + 54088b783a1SJames Wright Fu[Fuviscidx[j][1]]*dXdx[k][1] + 54188b783a1SJames Wright Fu[Fuviscidx[j][2]]*dXdx[k][2]); 54288b783a1SJames Wright // -- Total Energy Density 54388b783a1SJames Wright // ---- (E + P) u 54488b783a1SJames Wright for (int j=0; j<3; j++) 54588b783a1SJames Wright dv[j][4][i] += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 54688b783a1SJames Wright u[2]*dXdx[j][2]); 54788b783a1SJames Wright // ---- Fevisc 54888b783a1SJames Wright for (int j=0; j<3; j++) 54988b783a1SJames Wright dv[j][4][i] -= wdetJ * (Fe[0]*dXdx[j][0] + Fe[1]*dXdx[j][1] + 55088b783a1SJames Wright Fe[2]*dXdx[j][2]); 55188b783a1SJames Wright // Body Force 55288b783a1SJames Wright for (int j=0; j<5; j++) 55388b783a1SJames Wright v[j][i] = wdetJ * body_force[j]; 55488b783a1SJames Wright 55588626eedSJames Wright // Spatial Stabilization 55688626eedSJames Wright // -- Not used in favor of diagonal tau. Kept for future testing 55788626eedSJames Wright // const CeedScalar sound_speed = sqrt(gamma * P / rho); 55888626eedSJames Wright // CeedScalar Tau_x[3] = {0.}; 55988626eedSJames Wright // Tau_spatial(Tau_x, dXdx, u, sound_speed, context->c_tau, mu); 56088b783a1SJames Wright 56188626eedSJames Wright // -- Stabilization method: none, SU, or SUPG 56288626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 56388626eedSJames Wright CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0}; 56488626eedSJames Wright CeedScalar Tau_d[3] = {0.}; 56588b783a1SJames Wright switch (context->stabilization) { 56688b783a1SJames Wright case STAB_NONE: // Galerkin 56788b783a1SJames Wright break; 56888b783a1SJames Wright case STAB_SU: // SU 56988626eedSJames Wright Tau_diagPrim(Tau_d, dXdx, u, cv, context, mu, dt, rho); 57088626eedSJames Wright tau_strong_conv[0] = Tau_d[0] * strong_conv[0]; 57188626eedSJames Wright tau_strong_conv[1] = Tau_d[1] * strong_conv[1]; 57288626eedSJames Wright tau_strong_conv[2] = Tau_d[1] * strong_conv[2]; 57388626eedSJames Wright tau_strong_conv[3] = Tau_d[1] * strong_conv[3]; 57488626eedSJames Wright tau_strong_conv[4] = Tau_d[2] * strong_conv[4]; 57588626eedSJames Wright PrimitiveToConservative_fwd(rho, u, E, Rd, cv, tau_strong_conv, 57688626eedSJames Wright tau_strong_conv_conservative); 57788b783a1SJames Wright for (int j=0; j<3; j++) 57888b783a1SJames Wright for (int k=0; k<5; k++) 57988b783a1SJames Wright for (int l=0; l<5; l++) 58088626eedSJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l]; 58188b783a1SJames Wright 58288b783a1SJames Wright for (int j=0; j<5; j++) 58388b783a1SJames Wright for (int k=0; k<3; k++) 58488b783a1SJames Wright dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] + 58588b783a1SJames Wright stab[j][1] * dXdx[k][1] + 58688b783a1SJames Wright stab[j][2] * dXdx[k][2]); 58788b783a1SJames Wright break; 58888b783a1SJames Wright case STAB_SUPG: // SUPG is not implemented for explicit scheme 58988b783a1SJames Wright break; 59088b783a1SJames Wright } 59188b783a1SJames Wright 59288b783a1SJames Wright } // End Quadrature Point Loop 59388b783a1SJames Wright 59488b783a1SJames Wright // Return 59588b783a1SJames Wright return 0; 59688b783a1SJames Wright } 59788b783a1SJames Wright 59888b783a1SJames Wright // ***************************************************************************** 59988b783a1SJames Wright // This QFunction implements the Navier-Stokes equations (mentioned above) with 60088b783a1SJames Wright // implicit time stepping method 60188b783a1SJames Wright // 60288b783a1SJames Wright // SU = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) ) 60388b783a1SJames Wright // SUPG = Galerkin + grad(v) . ( Ai^T * Tau * (q_dot + Aj q,j - body force) ) 60488b783a1SJames Wright // (diffussive terms will be added later) 60588b783a1SJames Wright // 60688b783a1SJames Wright // ***************************************************************************** 60788b783a1SJames Wright CEED_QFUNCTION(IFunction_Newtonian)(void *ctx, CeedInt Q, 60888b783a1SJames Wright const CeedScalar *const *in, 60988b783a1SJames Wright CeedScalar *const *out) { 61088b783a1SJames Wright // *INDENT-OFF* 61188b783a1SJames Wright // Inputs 61288b783a1SJames Wright const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 61388b783a1SJames Wright (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 61488b783a1SJames Wright (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 61588b783a1SJames Wright (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3], 61688b783a1SJames Wright (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 61788b783a1SJames Wright // Outputs 61888b783a1SJames Wright CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 61988b783a1SJames Wright (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 62088b783a1SJames Wright // *INDENT-ON* 62188b783a1SJames Wright // Context 62288b783a1SJames Wright NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx; 62388b783a1SJames Wright const CeedScalar lambda = context->lambda; 62488b783a1SJames Wright const CeedScalar mu = context->mu; 62588b783a1SJames Wright const CeedScalar k = context->k; 62688b783a1SJames Wright const CeedScalar cv = context->cv; 62788b783a1SJames Wright const CeedScalar cp = context->cp; 62888626eedSJames Wright const CeedScalar *g = context->g; 62988626eedSJames Wright const CeedScalar dt = context->dt; 63088b783a1SJames Wright const CeedScalar gamma = cp / cv; 63188626eedSJames Wright const CeedScalar Rd = cp-cv; 63288b783a1SJames Wright 63388b783a1SJames Wright CeedPragmaSIMD 63488b783a1SJames Wright // Quadrature Point Loop 63588b783a1SJames Wright for (CeedInt i=0; i<Q; i++) { 63688b783a1SJames Wright // Setup 63788b783a1SJames Wright // -- Interp in 63888b783a1SJames Wright const CeedScalar rho = q[0][i]; 63988b783a1SJames Wright const CeedScalar u[3] = {q[1][i] / rho, 64088b783a1SJames Wright q[2][i] / rho, 64188b783a1SJames Wright q[3][i] / rho 64288b783a1SJames Wright }; 64388b783a1SJames Wright const CeedScalar E = q[4][i]; 64488b783a1SJames Wright // -- Grad in 64588b783a1SJames Wright const CeedScalar drho[3] = {dq[0][0][i], 64688b783a1SJames Wright dq[1][0][i], 64788b783a1SJames Wright dq[2][0][i] 64888b783a1SJames Wright }; 64988b783a1SJames Wright // *INDENT-OFF* 65088b783a1SJames Wright const CeedScalar dU[3][3] = {{dq[0][1][i], 65188b783a1SJames Wright dq[1][1][i], 65288b783a1SJames Wright dq[2][1][i]}, 65388b783a1SJames Wright {dq[0][2][i], 65488b783a1SJames Wright dq[1][2][i], 65588b783a1SJames Wright dq[2][2][i]}, 65688b783a1SJames Wright {dq[0][3][i], 65788b783a1SJames Wright dq[1][3][i], 65888b783a1SJames Wright dq[2][3][i]} 65988b783a1SJames Wright }; 66088b783a1SJames Wright // *INDENT-ON* 66188b783a1SJames Wright const CeedScalar dE[3] = {dq[0][4][i], 66288b783a1SJames Wright dq[1][4][i], 66388b783a1SJames Wright dq[2][4][i] 66488b783a1SJames Wright }; 66588b783a1SJames Wright // -- Interp-to-Interp q_data 66688b783a1SJames Wright const CeedScalar wdetJ = q_data[0][i]; 66788b783a1SJames Wright // -- Interp-to-Grad q_data 66888b783a1SJames Wright // ---- Inverse of change of coordinate matrix: X_i,j 66988b783a1SJames Wright // *INDENT-OFF* 67088b783a1SJames Wright const CeedScalar dXdx[3][3] = {{q_data[1][i], 67188b783a1SJames Wright q_data[2][i], 67288b783a1SJames Wright q_data[3][i]}, 67388b783a1SJames Wright {q_data[4][i], 67488b783a1SJames Wright q_data[5][i], 67588b783a1SJames Wright q_data[6][i]}, 67688b783a1SJames Wright {q_data[7][i], 67788b783a1SJames Wright q_data[8][i], 67888b783a1SJames Wright q_data[9][i]} 67988b783a1SJames Wright }; 68088626eedSJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 68188b783a1SJames Wright // *INDENT-ON* 68288b783a1SJames Wright // -- Grad-to-Grad q_data 68388b783a1SJames Wright // dU/dx 68488b783a1SJames Wright CeedScalar du[3][3] = {{0}}; 68588b783a1SJames Wright CeedScalar drhodx[3] = {0}; 68688b783a1SJames Wright CeedScalar dEdx[3] = {0}; 68788b783a1SJames Wright CeedScalar dUdx[3][3] = {{0}}; 68888b783a1SJames Wright CeedScalar dXdxdXdxT[3][3] = {{0}}; 68988b783a1SJames Wright for (int j=0; j<3; j++) { 69088b783a1SJames Wright for (int k=0; k<3; k++) { 69188b783a1SJames Wright du[j][k] = (dU[j][k] - drho[k]*u[j]) / rho; 69288b783a1SJames Wright drhodx[j] += drho[k] * dXdx[k][j]; 69388b783a1SJames Wright dEdx[j] += dE[k] * dXdx[k][j]; 69488b783a1SJames Wright for (int l=0; l<3; l++) { 69588b783a1SJames Wright dUdx[j][k] += dU[j][l] * dXdx[l][k]; 69688b783a1SJames Wright dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 69788b783a1SJames Wright } 69888b783a1SJames Wright } 69988b783a1SJames Wright } 70088b783a1SJames Wright CeedScalar dudx[3][3] = {{0}}; 70188b783a1SJames Wright for (int j=0; j<3; j++) 70288b783a1SJames Wright for (int k=0; k<3; k++) 70388b783a1SJames Wright for (int l=0; l<3; l++) 70488b783a1SJames Wright dudx[j][k] += du[j][l] * dXdx[l][k]; 70588b783a1SJames Wright // -- grad_T 70688b783a1SJames Wright const CeedScalar grad_T[3] = {(dEdx[0]/rho - E*drhodx[0]/(rho*rho) - /* *NOPAD* */ 70788626eedSJames Wright (u[0]*dudx[0][0] + u[1]*dudx[1][0] + u[2]*dudx[2][0]) + g[0])/cv, 70888b783a1SJames Wright (dEdx[1]/rho - E*drhodx[1]/(rho*rho) - /* *NOPAD* */ 70988626eedSJames Wright (u[0]*dudx[0][1] + u[1]*dudx[1][1] + u[2]*dudx[2][1]) + g[1])/cv, 71088b783a1SJames Wright (dEdx[2]/rho - E*drhodx[2]/(rho*rho) - /* *NOPAD* */ 71188626eedSJames Wright (u[0]*dudx[0][2] + u[1]*dudx[1][2] + u[2]*dudx[2][2]) + g[2])/cv 71288b783a1SJames Wright }; 71388b783a1SJames Wright // -- Fuvisc 71488b783a1SJames Wright // ---- Symmetric 3x3 matrix 71588b783a1SJames Wright const CeedScalar Fu[6] = {mu*(dudx[0][0] * (2 + lambda) + /* *NOPAD* */ 71688b783a1SJames Wright lambda * (dudx[1][1] + dudx[2][2])), 71788b783a1SJames Wright mu*(dudx[0][1] + dudx[1][0]), /* *NOPAD* */ 71888b783a1SJames Wright mu*(dudx[0][2] + dudx[2][0]), /* *NOPAD* */ 71988b783a1SJames Wright mu*(dudx[1][1] * (2 + lambda) + /* *NOPAD* */ 72088b783a1SJames Wright lambda * (dudx[0][0] + dudx[2][2])), 72188b783a1SJames Wright mu*(dudx[1][2] + dudx[2][1]), /* *NOPAD* */ 72288b783a1SJames Wright mu*(dudx[2][2] * (2 + lambda) + /* *NOPAD* */ 72388b783a1SJames Wright lambda * (dudx[0][0] + dudx[1][1])) 72488b783a1SJames Wright }; 72588b783a1SJames Wright // -- Fevisc 72688b783a1SJames Wright const CeedScalar Fe[3] = {u[0]*Fu[0] + u[1]*Fu[1] + u[2]*Fu[2] + /* *NOPAD* */ 72788b783a1SJames Wright k*grad_T[0], /* *NOPAD* */ 72888b783a1SJames Wright u[0]*Fu[1] + u[1]*Fu[3] + u[2]*Fu[4] + /* *NOPAD* */ 72988b783a1SJames Wright k*grad_T[1], /* *NOPAD* */ 73088b783a1SJames Wright u[0]*Fu[2] + u[1]*Fu[4] + u[2]*Fu[5] + /* *NOPAD* */ 73188b783a1SJames Wright k*grad_T[2] /* *NOPAD* */ 73288b783a1SJames Wright }; 73388b783a1SJames Wright // Pressure 73488b783a1SJames Wright const CeedScalar 73588b783a1SJames Wright E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 73688626eedSJames Wright E_potential = -rho*(g[0]*x_i[0] + g[1]*x_i[1] + g[2]*x_i[2]), 73788b783a1SJames Wright E_internal = E - E_kinetic - E_potential, 73888b783a1SJames Wright P = E_internal * (gamma - 1.); // P = pressure 73988b783a1SJames Wright 74088b783a1SJames Wright // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 74188b783a1SJames Wright CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 74288626eedSJames Wright computeFluxJacobian_NS(jacob_F_conv, rho, u, E, gamma, g, x_i); 74388b783a1SJames Wright 74488b783a1SJames Wright // dqdx collects drhodx, dUdx and dEdx in one vector 74588b783a1SJames Wright CeedScalar dqdx[5][3]; 74688b783a1SJames Wright for (int j=0; j<3; j++) { 74788b783a1SJames Wright dqdx[0][j] = drhodx[j]; 74888b783a1SJames Wright dqdx[4][j] = dEdx[j]; 74988b783a1SJames Wright for (int k=0; k<3; k++) 75088b783a1SJames Wright dqdx[k+1][j] = dUdx[k][j]; 75188b783a1SJames Wright } 75288b783a1SJames Wright // strong_conv = dF/dq * dq/dx (Strong convection) 75388b783a1SJames Wright CeedScalar strong_conv[5] = {0}; 75488b783a1SJames Wright for (int j=0; j<3; j++) 75588b783a1SJames Wright for (int k=0; k<5; k++) 75688b783a1SJames Wright for (int l=0; l<5; l++) 75788b783a1SJames Wright strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 75888b783a1SJames Wright 75988b783a1SJames Wright // Body force 76088626eedSJames Wright const CeedScalar body_force[5] = {0, rho *g[0], rho *g[1], rho *g[2], 0}; 76188b783a1SJames Wright 76288b783a1SJames Wright // Strong residual 76388b783a1SJames Wright CeedScalar strong_res[5]; 76488b783a1SJames Wright for (int j=0; j<5; j++) 76588b783a1SJames Wright strong_res[j] = q_dot[j][i] + strong_conv[j] - body_force[j]; 76688b783a1SJames Wright 76788b783a1SJames Wright // The Physics 76888b783a1SJames Wright //-----mass matrix 76988b783a1SJames Wright for (int j=0; j<5; j++) 77088b783a1SJames Wright v[j][i] = wdetJ*q_dot[j][i]; 77188b783a1SJames Wright 77288b783a1SJames Wright // Zero dv so all future terms can safely sum into it 77388b783a1SJames Wright for (int j=0; j<5; j++) 77488b783a1SJames Wright for (int k=0; k<3; k++) 77588b783a1SJames Wright dv[k][j][i] = 0; 77688b783a1SJames Wright 77788b783a1SJames Wright // -- Density 77888b783a1SJames Wright // ---- u rho 77988b783a1SJames Wright for (int j=0; j<3; j++) 78088b783a1SJames Wright dv[j][0][i] -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 78188b783a1SJames Wright rho*u[2]*dXdx[j][2]); 78288b783a1SJames Wright // -- Momentum 78388b783a1SJames Wright // ---- rho (u x u) + P I3 78488b783a1SJames Wright for (int j=0; j<3; j++) 78588b783a1SJames Wright for (int k=0; k<3; k++) 78688b783a1SJames Wright dv[k][j+1][i] -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0))*dXdx[k][0] + 78788b783a1SJames Wright (rho*u[j]*u[1] + (j==1?P:0))*dXdx[k][1] + 78888b783a1SJames Wright (rho*u[j]*u[2] + (j==2?P:0))*dXdx[k][2]); 78988b783a1SJames Wright // ---- Fuvisc 79088b783a1SJames Wright const CeedInt Fuviscidx[3][3] = {{0, 1, 2}, {1, 3, 4}, {2, 4, 5}}; // symmetric matrix indices 79188b783a1SJames Wright for (int j=0; j<3; j++) 79288b783a1SJames Wright for (int k=0; k<3; k++) 79388b783a1SJames Wright dv[k][j+1][i] += wdetJ*(Fu[Fuviscidx[j][0]]*dXdx[k][0] + 79488b783a1SJames Wright Fu[Fuviscidx[j][1]]*dXdx[k][1] + 79588b783a1SJames Wright Fu[Fuviscidx[j][2]]*dXdx[k][2]); 79688b783a1SJames Wright // -- Total Energy Density 79788b783a1SJames Wright // ---- (E + P) u 79888b783a1SJames Wright for (int j=0; j<3; j++) 79988b783a1SJames Wright dv[j][4][i] -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 80088b783a1SJames Wright u[2]*dXdx[j][2]); 80188b783a1SJames Wright // ---- Fevisc 80288b783a1SJames Wright for (int j=0; j<3; j++) 80388b783a1SJames Wright dv[j][4][i] += wdetJ * (Fe[0]*dXdx[j][0] + Fe[1]*dXdx[j][1] + 80488b783a1SJames Wright Fe[2]*dXdx[j][2]); 80588b783a1SJames Wright // Body Force 80688b783a1SJames Wright for (int j=0; j<5; j++) 80788b783a1SJames Wright v[j][i] -= wdetJ*body_force[j]; 80888b783a1SJames Wright 80988626eedSJames Wright // Spatial Stabilization 81088626eedSJames Wright // -- Not used in favor of diagonal tau. Kept for future testing 81188626eedSJames Wright // const CeedScalar sound_speed = sqrt(gamma * P / rho); 81288626eedSJames Wright // CeedScalar Tau_x[3] = {0.}; 81388626eedSJames Wright // Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau, mu); 81488b783a1SJames Wright 81588b783a1SJames Wright // -- Stabilization method: none, SU, or SUPG 81688626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 81788626eedSJames Wright CeedScalar tau_strong_res[5] = {0.}, tau_strong_res_conservative[5] = {0}; 81888626eedSJames Wright CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0}; 81988626eedSJames Wright CeedScalar Tau_d[3] = {0.}; 82088b783a1SJames Wright switch (context->stabilization) { 82188b783a1SJames Wright case STAB_NONE: // Galerkin 82288b783a1SJames Wright break; 82388b783a1SJames Wright case STAB_SU: // SU 82488626eedSJames Wright Tau_diagPrim(Tau_d, dXdx, u, cv, context, mu, dt, rho); 82588626eedSJames Wright tau_strong_conv[0] = Tau_d[0] * strong_conv[0]; 82688626eedSJames Wright tau_strong_conv[1] = Tau_d[1] * strong_conv[1]; 82788626eedSJames Wright tau_strong_conv[2] = Tau_d[1] * strong_conv[2]; 82888626eedSJames Wright tau_strong_conv[3] = Tau_d[1] * strong_conv[3]; 82988626eedSJames Wright tau_strong_conv[4] = Tau_d[2] * strong_conv[4]; 83088626eedSJames Wright PrimitiveToConservative_fwd(rho, u, E, Rd, cv, tau_strong_conv, 83188626eedSJames Wright tau_strong_conv_conservative); 83288b783a1SJames Wright for (int j=0; j<3; j++) 83388b783a1SJames Wright for (int k=0; k<5; k++) 83488b783a1SJames Wright for (int l=0; l<5; l++) 83588626eedSJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l]; 83688b783a1SJames Wright 83788b783a1SJames Wright for (int j=0; j<5; j++) 83888b783a1SJames Wright for (int k=0; k<3; k++) 83988b783a1SJames Wright dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 84088b783a1SJames Wright stab[j][1] * dXdx[k][1] + 84188b783a1SJames Wright stab[j][2] * dXdx[k][2]); 84288b783a1SJames Wright break; 84388b783a1SJames Wright case STAB_SUPG: // SUPG 84488626eedSJames Wright Tau_diagPrim(Tau_d, dXdx, u, cv, context, mu, dt, rho); 84588626eedSJames Wright tau_strong_res[0] = Tau_d[0] * strong_res[0]; 84688626eedSJames Wright tau_strong_res[1] = Tau_d[1] * strong_res[1]; 84788626eedSJames Wright tau_strong_res[2] = Tau_d[1] * strong_res[2]; 84888626eedSJames Wright tau_strong_res[3] = Tau_d[1] * strong_res[3]; 84988626eedSJames Wright tau_strong_res[4] = Tau_d[2] * strong_res[4]; 85088626eedSJames Wright // Alternate route (useful later with primitive variable code) 85188626eedSJames Wright // this function was verified against PHASTA for as IC that was as close as possible 85288626eedSJames Wright // computeFluxJacobian_NSp(jacob_F_conv_p, rho, u, E, Rd, cv); 85388626eedSJames Wright // it has also been verified to compute a correct through the following 85488626eedSJames Wright // stab[k][j] += jacob_F_conv_p[j][k][l] * tau_strong_res[l] // flux Jacobian wrt primitive 85588626eedSJames Wright // applied in the triple loop below 85688626eedSJames Wright // However, it is more flops than using the existing Jacobian wrt q after q_{,Y} viz 85788626eedSJames Wright PrimitiveToConservative_fwd(rho, u, E, Rd, cv, tau_strong_res, 85888626eedSJames Wright tau_strong_res_conservative); 85988b783a1SJames Wright for (int j=0; j<3; j++) 86088b783a1SJames Wright for (int k=0; k<5; k++) 86188b783a1SJames Wright for (int l=0; l<5; l++) 86288626eedSJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_res_conservative[l]; 86388b783a1SJames Wright 86488b783a1SJames Wright for (int j=0; j<5; j++) 86588b783a1SJames Wright for (int k=0; k<3; k++) 86688b783a1SJames Wright dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 86788b783a1SJames Wright stab[j][1] * dXdx[k][1] + 86888b783a1SJames Wright stab[j][2] * dXdx[k][2]); 86988b783a1SJames Wright break; 87088b783a1SJames Wright } 87188b783a1SJames Wright 87288b783a1SJames Wright } // End Quadrature Point Loop 87388b783a1SJames Wright 87488b783a1SJames Wright // Return 87588b783a1SJames Wright return 0; 87688b783a1SJames Wright } 87788b783a1SJames Wright // ***************************************************************************** 87888b783a1SJames Wright #endif // newtonian_h 879