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> 17*841e4c73SJed Brown #include "newtonian_types.h" 1888b783a1SJames Wright 1988b783a1SJames Wright #ifndef M_PI 2088b783a1SJames Wright #define M_PI 3.14159265358979323846 2188b783a1SJames Wright #endif 2288b783a1SJames Wright 2388b783a1SJames Wright #ifndef setup_context_struct 2488b783a1SJames Wright #define setup_context_struct 2588b783a1SJames Wright typedef struct SetupContext_ *SetupContext; 2688b783a1SJames Wright struct SetupContext_ { 2788b783a1SJames Wright CeedScalar theta0; 2888b783a1SJames Wright CeedScalar thetaC; 2988b783a1SJames Wright CeedScalar P0; 3088b783a1SJames Wright CeedScalar N; 3188b783a1SJames Wright CeedScalar cv; 3288b783a1SJames Wright CeedScalar cp; 3388626eedSJames Wright CeedScalar g[3]; 3488b783a1SJames Wright CeedScalar rc; 3588b783a1SJames Wright CeedScalar lx; 3688b783a1SJames Wright CeedScalar ly; 3788b783a1SJames Wright CeedScalar lz; 3888b783a1SJames Wright CeedScalar center[3]; 3988b783a1SJames Wright CeedScalar dc_axis[3]; 4088b783a1SJames Wright CeedScalar wind[3]; 4188b783a1SJames Wright CeedScalar time; 4288b783a1SJames Wright int wind_type; // See WindType: 0=ROTATION, 1=TRANSLATION 4388b783a1SJames Wright int bubble_type; // See BubbleType: 0=SPHERE, 1=CYLINDER 4488b783a1SJames Wright int bubble_continuity_type; // See BubbleContinuityType: 0=SMOOTH, 1=BACK_SHARP 2=THICK 4588b783a1SJames Wright }; 4688b783a1SJames Wright #endif 4788b783a1SJames Wright 4888b783a1SJames Wright // ***************************************************************************** 4988b783a1SJames Wright // Helper function for computing flux Jacobian 5088b783a1SJames Wright // ***************************************************************************** 5188b783a1SJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NS(CeedScalar dF[3][5][5], 5288b783a1SJames Wright const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 5388626eedSJames Wright const CeedScalar gamma, const CeedScalar g[3], const CeedScalar x[3]) { 5488b783a1SJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 5588626eedSJames Wright CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]); 5688b783a1SJames Wright for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 5788b783a1SJames Wright for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix 5888626eedSJames Wright dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2. - e_potential)) : 0.) - 5988626eedSJames Wright u[i]*u[j]; 6088b783a1SJames Wright for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix 6188b783a1SJames Wright dF[i][0][k+1] = ((i==k) ? 1. : 0.); 6288b783a1SJames Wright dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) + 6388b783a1SJames Wright ((i==k) ? u[j] : 0.) - 6488b783a1SJames Wright ((i==j) ? u[k] : 0.) * (gamma-1.); 6588b783a1SJames Wright dF[i][4][k+1] = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) - 6688b783a1SJames Wright (gamma-1.)*u[i]*u[k]; 6788b783a1SJames Wright } 6888b783a1SJames Wright dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.); 6988b783a1SJames Wright } 7088b783a1SJames Wright dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho); 7188b783a1SJames Wright dF[i][4][4] = u[i] * gamma; 7288b783a1SJames Wright } 7388b783a1SJames Wright } 7488b783a1SJames Wright 7588b783a1SJames Wright // ***************************************************************************** 7688626eedSJames Wright // Helper function for computing flux Jacobian of Primitive variables 7788626eedSJames Wright // ***************************************************************************** 7888626eedSJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NSp(CeedScalar dF[3][5][5], 7988626eedSJames Wright const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 8088626eedSJames Wright const CeedScalar Rd, const CeedScalar cv) { 8188626eedSJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 8288626eedSJames Wright // TODO Add in gravity's contribution 8388626eedSJames Wright 8488626eedSJames Wright CeedScalar T = ( E / rho - u_sq / 2. ) / cv; 8588626eedSJames Wright CeedScalar drdT = -rho / T; 8688626eedSJames Wright CeedScalar drdP = 1. / ( Rd * T); 8788626eedSJames Wright CeedScalar etot = E / rho ; 8888626eedSJames Wright CeedScalar e2p = drdP * etot + 1. ; 8988626eedSJames Wright CeedScalar e3p = ( E + rho * Rd * T ); 9088626eedSJames Wright CeedScalar e4p = drdT * etot + rho * cv ; 9188626eedSJames Wright 9288626eedSJames Wright for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 9388626eedSJames Wright for (CeedInt j=0; j<3; j++) { // j counts F^{m_j} 9488626eedSJames Wright // [row][col] of A_i 9588626eedSJames Wright dF[i][j+1][0] = drdP * u[i] * u[j] + ((i==j) ? 1. : 0.); // F^{{m_j} wrt p 9688626eedSJames Wright for (CeedInt k=0; k<3; k++) { // k counts the wrt vel_k 97871db79fSKenneth E. Jansen dF[i][0][k+1] = ((i==k) ? rho : 0.); // F^c wrt u_k 9888626eedSJames Wright dF[i][j+1][k+1] = (((j==k) ? u[i] : 0.) + // F^m_j wrt u_k 9988626eedSJames Wright ((i==k) ? u[j] : 0.) ) * rho; 10088626eedSJames Wright dF[i][4][k+1] = rho * u[i] * u[k] 10188626eedSJames Wright + ((i==k) ? e3p : 0.) ; // F^e wrt u_k 10288626eedSJames Wright } 10388626eedSJames Wright dF[i][j+1][4] = drdT * u[i] * u[j]; // F^{m_j} wrt T 10488626eedSJames Wright } 10588626eedSJames Wright dF[i][4][0] = u[i] * e2p; // F^e wrt p 10688626eedSJames Wright dF[i][4][4] = u[i] * e4p; // F^e wrt T 10788626eedSJames Wright dF[i][0][0] = u[i] * drdP; // F^c wrt p 10888626eedSJames Wright dF[i][0][4] = u[i] * drdT; // F^c wrt T 10988626eedSJames Wright } 11088626eedSJames Wright } 11188626eedSJames Wright 11288626eedSJames Wright CEED_QFUNCTION_HELPER void PrimitiveToConservative_fwd(const CeedScalar rho, 11388626eedSJames Wright const CeedScalar u[3], const CeedScalar E, const CeedScalar Rd, 11488626eedSJames Wright const CeedScalar cv, const CeedScalar dY[5], CeedScalar dU[5]) { 11588626eedSJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; 11688626eedSJames Wright CeedScalar T = ( E / rho - u_sq / 2. ) / cv; 11788626eedSJames Wright CeedScalar drdT = -rho / T; 11888626eedSJames Wright CeedScalar drdP = 1. / ( Rd * T); 11988626eedSJames Wright dU[0] = drdP * dY[0] + drdT * dY[4]; 12088626eedSJames Wright CeedScalar de_kinetic = 0; 12188626eedSJames Wright for (int i=0; i<3; i++) { 12288626eedSJames Wright dU[1+i] = dU[0] * u[i] + rho * dY[1+i]; 12388626eedSJames Wright de_kinetic += u[i] * dY[1+i]; 12488626eedSJames Wright } 12588626eedSJames Wright dU[4] = rho * cv * dY[4] + dU[0] * cv * T // internal energy: rho * e 12688626eedSJames Wright + rho * de_kinetic + .5 * dU[0] * u_sq; // kinetic energy: .5 * rho * |u|^2 12788626eedSJames Wright } 12888626eedSJames Wright 12988626eedSJames Wright // ***************************************************************************** 13088626eedSJames Wright // Helper function for computing Tau elements (stabilization constant) 13188626eedSJames Wright // Model from: 13288626eedSJames Wright // PHASTA 13388626eedSJames Wright // 13488626eedSJames Wright // Tau[i] = itau=0 which is diagonal-Shakib (3 values still but not spatial) 13588626eedSJames Wright // 13688626eedSJames Wright // Where NOT UPDATED YET 13788626eedSJames Wright // ***************************************************************************** 13888626eedSJames Wright CEED_QFUNCTION_HELPER void Tau_diagPrim(CeedScalar Tau_d[3], 13988626eedSJames Wright const CeedScalar dXdx[3][3], const CeedScalar u[3], 14088626eedSJames Wright const CeedScalar cv, const NewtonianIdealGasContext newt_ctx, 14188626eedSJames Wright const CeedScalar mu, const CeedScalar dt, 14288626eedSJames Wright const CeedScalar rho) { 14388626eedSJames Wright // Context 14488626eedSJames Wright const CeedScalar Ctau_t = newt_ctx->Ctau_t; 14588626eedSJames Wright const CeedScalar Ctau_v = newt_ctx->Ctau_v; 14688626eedSJames Wright const CeedScalar Ctau_C = newt_ctx->Ctau_C; 14788626eedSJames Wright const CeedScalar Ctau_M = newt_ctx->Ctau_M; 14888626eedSJames Wright const CeedScalar Ctau_E = newt_ctx->Ctau_E; 14988626eedSJames Wright CeedScalar gijd[6]; 15088626eedSJames Wright CeedScalar tau; 15188626eedSJames Wright CeedScalar dts; 15288626eedSJames Wright CeedScalar fact; 15388626eedSJames Wright 15488626eedSJames Wright //*INDENT-OFF* 15588626eedSJames Wright gijd[0] = dXdx[0][0] * dXdx[0][0] 15688626eedSJames Wright + dXdx[1][0] * dXdx[1][0] 15788626eedSJames Wright + dXdx[2][0] * dXdx[2][0]; 15888626eedSJames Wright 15988626eedSJames Wright gijd[1] = dXdx[0][0] * dXdx[0][1] 16088626eedSJames Wright + dXdx[1][0] * dXdx[1][1] 16188626eedSJames Wright + dXdx[2][0] * dXdx[2][1]; 16288626eedSJames Wright 16388626eedSJames Wright gijd[2] = dXdx[0][1] * dXdx[0][1] 16488626eedSJames Wright + dXdx[1][1] * dXdx[1][1] 16588626eedSJames Wright + dXdx[2][1] * dXdx[2][1]; 16688626eedSJames Wright 16788626eedSJames Wright gijd[3] = dXdx[0][0] * dXdx[0][2] 16888626eedSJames Wright + dXdx[1][0] * dXdx[1][2] 16988626eedSJames Wright + dXdx[2][0] * dXdx[2][2]; 17088626eedSJames Wright 17188626eedSJames Wright gijd[4] = dXdx[0][1] * dXdx[0][2] 17288626eedSJames Wright + dXdx[1][1] * dXdx[1][2] 17388626eedSJames Wright + dXdx[2][1] * dXdx[2][2]; 17488626eedSJames Wright 17588626eedSJames Wright gijd[5] = dXdx[0][2] * dXdx[0][2] 17688626eedSJames Wright + dXdx[1][2] * dXdx[1][2] 17788626eedSJames Wright + dXdx[2][2] * dXdx[2][2]; 17888626eedSJames Wright //*INDENT-ON* 17988626eedSJames Wright 18088626eedSJames Wright dts = Ctau_t / dt ; 18188626eedSJames Wright 18288626eedSJames Wright tau = rho*rho*((4. * dts * dts) 18388626eedSJames Wright + u[0] * ( u[0] * gijd[0] + 2. * ( u[1] * gijd[1] + u[2] * gijd[3])) 18488626eedSJames Wright + u[1] * ( u[1] * gijd[2] + 2. * u[2] * gijd[4]) 18588626eedSJames Wright + u[2] * u[2] * gijd[5]) 18688626eedSJames Wright + Ctau_v* mu * mu * 18788626eedSJames Wright (gijd[0]*gijd[0] + gijd[2]*gijd[2] + gijd[5]*gijd[5] + 18888626eedSJames Wright + 2. * (gijd[1]*gijd[1] + gijd[3]*gijd[3] + gijd[4]*gijd[4])); 18988626eedSJames Wright 19088626eedSJames Wright fact=sqrt(tau); 19188626eedSJames Wright 19288626eedSJames Wright Tau_d[0] = Ctau_C * fact / (rho*(gijd[0] + gijd[2] + gijd[5]))*0.125; 19388626eedSJames Wright 19488626eedSJames Wright Tau_d[1] = Ctau_M / fact; 19588626eedSJames Wright Tau_d[2] = Ctau_E / ( fact * cv ); 19688626eedSJames Wright 19788626eedSJames Wright // consider putting back the way I initially had it Ctau_E * Tau_d[1] /cv 19888626eedSJames Wright // to avoid a division if the compiler is smart enough to see that cv IS 19988626eedSJames Wright // a constant that it could invert once for all elements 20088626eedSJames Wright // but in that case energy tau is scaled by the product of Ctau_E * Ctau_M 20188626eedSJames Wright // OR we could absorb cv into Ctau_E but this puts more burden on user to 20288626eedSJames Wright // know how to change constants with a change of fluid or units. Same for 20388626eedSJames Wright // Ctau_v * mu * mu IF AND ONLY IF we don't add viscosity law =f(T) 20488626eedSJames Wright } 20588626eedSJames Wright 20688626eedSJames Wright // ***************************************************************************** 20788b783a1SJames Wright // Helper function for computing Tau elements (stabilization constant) 20888b783a1SJames Wright // Model from: 20988b783a1SJames Wright // Stabilized Methods for Compressible Flows, Hughes et al 2010 21088b783a1SJames Wright // 21188b783a1SJames Wright // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 21288b783a1SJames Wright // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 21388b783a1SJames Wright // 21488b783a1SJames Wright // Where 21588b783a1SJames Wright // c_tau = stabilization constant (0.5 is reported as "optimal") 21688b783a1SJames Wright // h[i] = 2 length(dxdX[i]) 21788b783a1SJames Wright // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 21888b783a1SJames Wright // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 21988b783a1SJames Wright // rho(A[i]) = spectral radius of the convective flux Jacobian i, 22088b783a1SJames Wright // wave speed in direction i 22188b783a1SJames Wright // ***************************************************************************** 22288b783a1SJames Wright CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], 22388b783a1SJames Wright const CeedScalar dXdx[3][3], const CeedScalar u[3], 22488626eedSJames Wright /* const CeedScalar sound_speed, const CeedScalar c_tau) { */ 22588626eedSJames Wright const CeedScalar sound_speed, const CeedScalar c_tau, 22688626eedSJames Wright const CeedScalar viscosity) { 22788626eedSJames Wright const CeedScalar mag_u_visc = sqrt(u[0]*u[0] +u[1]*u[1] +u[2]*u[2]) / 22888626eedSJames Wright (2*viscosity); 22988b783a1SJames Wright for (int i=0; i<3; i++) { 23088b783a1SJames Wright // length of element in direction i 23188b783a1SJames Wright CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] + 23288b783a1SJames Wright dXdx[2][i]*dXdx[2][i]); 23388626eedSJames Wright CeedScalar Pe = mag_u_visc*h; 23488626eedSJames Wright CeedScalar Xi = 1/tanh(Pe) - 1/Pe; 23588b783a1SJames Wright // fastest wave in direction i 23688b783a1SJames Wright CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 23788626eedSJames Wright Tau_x[i] = c_tau * h * Xi / fastest_wave; 23888b783a1SJames Wright } 23988b783a1SJames Wright } 24088b783a1SJames Wright 24188b783a1SJames Wright // ***************************************************************************** 24288b783a1SJames Wright // This QFunction sets a "still" initial condition for generic Newtonian IG problems 24388b783a1SJames Wright // ***************************************************************************** 24488b783a1SJames Wright CEED_QFUNCTION(ICsNewtonianIG)(void *ctx, CeedInt Q, 24588b783a1SJames Wright const CeedScalar *const *in, CeedScalar *const *out) { 24688b783a1SJames Wright // Inputs 24788b783a1SJames Wright const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 24888b783a1SJames Wright 24988b783a1SJames Wright // Outputs 25088b783a1SJames Wright CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 25188b783a1SJames Wright 25288626eedSJames Wright // Context 25388626eedSJames Wright const SetupContext context = (SetupContext)ctx; 25488626eedSJames Wright const CeedScalar theta0 = context->theta0; 25588626eedSJames Wright const CeedScalar P0 = context->P0; 25688626eedSJames Wright const CeedScalar cv = context->cv; 25788626eedSJames Wright const CeedScalar cp = context->cp; 25888626eedSJames Wright const CeedScalar *g = context->g; 25988626eedSJames Wright const CeedScalar Rd = cp - cv; 26088626eedSJames Wright 26188b783a1SJames Wright // Quadrature Point Loop 26288b783a1SJames Wright CeedPragmaSIMD 26388b783a1SJames Wright for (CeedInt i=0; i<Q; i++) { 26488b783a1SJames Wright CeedScalar q[5] = {0.}; 26588b783a1SJames Wright 26688b783a1SJames Wright // Setup 26788b783a1SJames Wright // -- Coordinates 26888626eedSJames Wright const CeedScalar x[3] = {X[0][i], X[1][i], X[2][i]}; 26988626eedSJames Wright const CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]); 27088b783a1SJames Wright 27188b783a1SJames Wright // -- Density 27288626eedSJames Wright const CeedScalar rho = P0 / (Rd*theta0); 27388b783a1SJames Wright 27488b783a1SJames Wright // Initial Conditions 27588b783a1SJames Wright q[0] = rho; 27688b783a1SJames Wright q[1] = 0.0; 27788b783a1SJames Wright q[2] = 0.0; 27888b783a1SJames Wright q[3] = 0.0; 27988626eedSJames Wright q[4] = rho * (cv*theta0 + e_potential); 28088b783a1SJames Wright 28188b783a1SJames Wright for (CeedInt j=0; j<5; j++) 28288b783a1SJames Wright q0[j][i] = q[j]; 28388b783a1SJames Wright } // End of Quadrature Point Loop 28488b783a1SJames Wright return 0; 28588b783a1SJames Wright } 28688b783a1SJames Wright 28788b783a1SJames Wright // ***************************************************************************** 28888b783a1SJames Wright // This QFunction implements the following formulation of Navier-Stokes with 28988b783a1SJames Wright // explicit time stepping method 29088b783a1SJames Wright // 29188b783a1SJames Wright // This is 3D compressible Navier-Stokes in conservation form with state 29288b783a1SJames Wright // variables of density, momentum density, and total energy density. 29388b783a1SJames Wright // 29488b783a1SJames Wright // State Variables: q = ( rho, U1, U2, U3, E ) 29588b783a1SJames Wright // rho - Mass Density 29688b783a1SJames Wright // Ui - Momentum Density, Ui = rho ui 29788b783a1SJames Wright // E - Total Energy Density, E = rho (cv T + (u u)/2 + g z) 29888b783a1SJames Wright // 29988b783a1SJames Wright // Navier-Stokes Equations: 30088b783a1SJames Wright // drho/dt + div( U ) = 0 30188b783a1SJames Wright // dU/dt + div( rho (u x u) + P I3 ) + rho g khat = div( Fu ) 30288b783a1SJames Wright // dE/dt + div( (E + P) u ) = div( Fe ) 30388b783a1SJames Wright // 30488b783a1SJames Wright // Viscous Stress: 30588b783a1SJames Wright // Fu = mu (grad( u ) + grad( u )^T + lambda div ( u ) I3) 30688b783a1SJames Wright // 30788b783a1SJames Wright // Thermal Stress: 30888b783a1SJames Wright // Fe = u Fu + k grad( T ) 30988626eedSJames Wright // Equation of State 31088b783a1SJames Wright // P = (gamma - 1) (E - rho (u u) / 2 - rho g z) 31188b783a1SJames Wright // 31288b783a1SJames Wright // Stabilization: 31388b783a1SJames Wright // Tau = diag(TauC, TauM, TauM, TauM, TauE) 31488b783a1SJames Wright // f1 = rho sqrt(ui uj gij) 31588b783a1SJames Wright // gij = dXi/dX * dXi/dX 31688b783a1SJames Wright // TauC = Cc f1 / (8 gii) 31788b783a1SJames Wright // TauM = min( 1 , 1 / f1 ) 31888b783a1SJames Wright // TauE = TauM / (Ce cv) 31988b783a1SJames Wright // 32088b783a1SJames Wright // SU = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) ) 32188b783a1SJames Wright // 32288b783a1SJames Wright // Constants: 32388b783a1SJames Wright // lambda = - 2 / 3, From Stokes hypothesis 32488b783a1SJames Wright // mu , Dynamic viscosity 32588b783a1SJames Wright // k , Thermal conductivity 32688b783a1SJames Wright // cv , Specific heat, constant volume 32788b783a1SJames Wright // cp , Specific heat, constant pressure 32888b783a1SJames Wright // g , Gravity 32988b783a1SJames Wright // gamma = cp / cv, Specific heat ratio 33088b783a1SJames Wright // 33188b783a1SJames Wright // We require the product of the inverse of the Jacobian (dXdx_j,k) and 33288b783a1SJames Wright // its transpose (dXdx_k,j) to properly compute integrals of the form: 33388b783a1SJames Wright // int( gradv gradu ) 33488b783a1SJames Wright // 33588b783a1SJames Wright // ***************************************************************************** 33688b783a1SJames Wright CEED_QFUNCTION(Newtonian)(void *ctx, CeedInt Q, 33788b783a1SJames Wright const CeedScalar *const *in, CeedScalar *const *out) { 33888b783a1SJames Wright // *INDENT-OFF* 33988b783a1SJames Wright // Inputs 34088b783a1SJames Wright const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 34188b783a1SJames Wright (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 34288b783a1SJames Wright (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 34388b783a1SJames Wright (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 34488b783a1SJames Wright // Outputs 34588b783a1SJames Wright CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 34688b783a1SJames Wright (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 34788b783a1SJames Wright // *INDENT-ON* 34888b783a1SJames Wright 34988b783a1SJames Wright // Context 35088b783a1SJames Wright NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx; 35188b783a1SJames Wright const CeedScalar lambda = context->lambda; 35288b783a1SJames Wright const CeedScalar mu = context->mu; 35388b783a1SJames Wright const CeedScalar k = context->k; 35488b783a1SJames Wright const CeedScalar cv = context->cv; 35588b783a1SJames Wright const CeedScalar cp = context->cp; 35688626eedSJames Wright const CeedScalar *g = context->g; 35788626eedSJames Wright const CeedScalar dt = context->dt; 35888b783a1SJames Wright const CeedScalar gamma = cp / cv; 35988626eedSJames Wright const CeedScalar Rd = cp - cv; 36088b783a1SJames Wright 36188b783a1SJames Wright CeedPragmaSIMD 36288b783a1SJames Wright // Quadrature Point Loop 36388b783a1SJames Wright for (CeedInt i=0; i<Q; i++) { 36488b783a1SJames Wright // *INDENT-OFF* 36588b783a1SJames Wright // Setup 36688b783a1SJames Wright // -- Interp in 36788b783a1SJames Wright const CeedScalar rho = q[0][i]; 36888b783a1SJames Wright const CeedScalar u[3] = {q[1][i] / rho, 36988b783a1SJames Wright q[2][i] / rho, 37088b783a1SJames Wright q[3][i] / rho 37188b783a1SJames Wright }; 37288b783a1SJames Wright const CeedScalar E = q[4][i]; 37388b783a1SJames Wright // -- Grad in 37488b783a1SJames Wright const CeedScalar drho[3] = {dq[0][0][i], 37588b783a1SJames Wright dq[1][0][i], 37688b783a1SJames Wright dq[2][0][i] 37788b783a1SJames Wright }; 37888b783a1SJames Wright const CeedScalar dU[3][3] = {{dq[0][1][i], 37988b783a1SJames Wright dq[1][1][i], 38088b783a1SJames Wright dq[2][1][i]}, 38188b783a1SJames Wright {dq[0][2][i], 38288b783a1SJames Wright dq[1][2][i], 38388b783a1SJames Wright dq[2][2][i]}, 38488b783a1SJames Wright {dq[0][3][i], 38588b783a1SJames Wright dq[1][3][i], 38688b783a1SJames Wright dq[2][3][i]} 38788b783a1SJames Wright }; 38888b783a1SJames Wright const CeedScalar dE[3] = {dq[0][4][i], 38988b783a1SJames Wright dq[1][4][i], 39088b783a1SJames Wright dq[2][4][i] 39188b783a1SJames Wright }; 39288b783a1SJames Wright // -- Interp-to-Interp q_data 39388b783a1SJames Wright const CeedScalar wdetJ = q_data[0][i]; 39488b783a1SJames Wright // -- Interp-to-Grad q_data 39588b783a1SJames Wright // ---- Inverse of change of coordinate matrix: X_i,j 39688b783a1SJames Wright // *INDENT-OFF* 39788b783a1SJames Wright const CeedScalar dXdx[3][3] = {{q_data[1][i], 39888b783a1SJames Wright q_data[2][i], 39988b783a1SJames Wright q_data[3][i]}, 40088b783a1SJames Wright {q_data[4][i], 40188b783a1SJames Wright q_data[5][i], 40288b783a1SJames Wright q_data[6][i]}, 40388b783a1SJames Wright {q_data[7][i], 40488b783a1SJames Wright q_data[8][i], 40588b783a1SJames Wright q_data[9][i]} 40688b783a1SJames Wright }; 40788626eedSJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 40888b783a1SJames Wright // *INDENT-ON* 40988b783a1SJames Wright // -- Grad-to-Grad q_data 41088b783a1SJames Wright // dU/dx 41188b783a1SJames Wright CeedScalar du[3][3] = {{0}}; 41288b783a1SJames Wright CeedScalar drhodx[3] = {0}; 41388b783a1SJames Wright CeedScalar dEdx[3] = {0}; 41488b783a1SJames Wright CeedScalar dUdx[3][3] = {{0}}; 41588b783a1SJames Wright CeedScalar dXdxdXdxT[3][3] = {{0}}; 41688b783a1SJames Wright for (int j=0; j<3; j++) { 41788b783a1SJames Wright for (int k=0; k<3; k++) { 41888b783a1SJames Wright du[j][k] = (dU[j][k] - drho[k]*u[j]) / rho; 41988b783a1SJames Wright drhodx[j] += drho[k] * dXdx[k][j]; 42088b783a1SJames Wright dEdx[j] += dE[k] * dXdx[k][j]; 42188b783a1SJames Wright for (int l=0; l<3; l++) { 42288b783a1SJames Wright dUdx[j][k] += dU[j][l] * dXdx[l][k]; 42388b783a1SJames Wright dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 42488b783a1SJames Wright } 42588b783a1SJames Wright } 42688b783a1SJames Wright } 42788b783a1SJames Wright CeedScalar dudx[3][3] = {{0}}; 42888b783a1SJames Wright for (int j=0; j<3; j++) 42988b783a1SJames Wright for (int k=0; k<3; k++) 43088b783a1SJames Wright for (int l=0; l<3; l++) 43188b783a1SJames Wright dudx[j][k] += du[j][l] * dXdx[l][k]; 43288b783a1SJames Wright // -- grad_T 43388b783a1SJames Wright const CeedScalar grad_T[3] = {(dEdx[0]/rho - E*drhodx[0]/(rho*rho) - /* *NOPAD* */ 43488626eedSJames Wright (u[0]*dudx[0][0] + u[1]*dudx[1][0] + u[2]*dudx[2][0]) + g[0])/cv, 43588b783a1SJames Wright (dEdx[1]/rho - E*drhodx[1]/(rho*rho) - /* *NOPAD* */ 43688626eedSJames Wright (u[0]*dudx[0][1] + u[1]*dudx[1][1] + u[2]*dudx[2][1]) + g[1])/cv, 43788b783a1SJames Wright (dEdx[2]/rho - E*drhodx[2]/(rho*rho) - /* *NOPAD* */ 43888626eedSJames Wright (u[0]*dudx[0][2] + u[1]*dudx[1][2] + u[2]*dudx[2][2]) + g[2])/cv 43988b783a1SJames Wright }; 44088b783a1SJames Wright 44188b783a1SJames Wright // -- Fuvisc 44288b783a1SJames Wright // ---- Symmetric 3x3 matrix 44388b783a1SJames Wright const CeedScalar Fu[6] = {mu*(dudx[0][0] * (2 + lambda) + /* *NOPAD* */ 44488b783a1SJames Wright lambda * (dudx[1][1] + dudx[2][2])), 44588b783a1SJames Wright mu*(dudx[0][1] + dudx[1][0]), /* *NOPAD* */ 44688b783a1SJames Wright mu*(dudx[0][2] + dudx[2][0]), /* *NOPAD* */ 44788b783a1SJames Wright mu*(dudx[1][1] * (2 + lambda) + /* *NOPAD* */ 44888b783a1SJames Wright lambda * (dudx[0][0] + dudx[2][2])), 44988b783a1SJames Wright mu*(dudx[1][2] + dudx[2][1]), /* *NOPAD* */ 45088b783a1SJames Wright mu*(dudx[2][2] * (2 + lambda) + /* *NOPAD* */ 45188b783a1SJames Wright lambda * (dudx[0][0] + dudx[1][1])) 45288b783a1SJames Wright }; 45388b783a1SJames Wright // -- Fevisc 45488b783a1SJames Wright const CeedScalar Fe[3] = {u[0]*Fu[0] + u[1]*Fu[1] + u[2]*Fu[2] + /* *NOPAD* */ 45588b783a1SJames Wright k*grad_T[0], /* *NOPAD* */ 45688b783a1SJames Wright u[0]*Fu[1] + u[1]*Fu[3] + u[2]*Fu[4] + /* *NOPAD* */ 45788b783a1SJames Wright k*grad_T[1], /* *NOPAD* */ 45888b783a1SJames Wright u[0]*Fu[2] + u[1]*Fu[4] + u[2]*Fu[5] + /* *NOPAD* */ 45988b783a1SJames Wright k*grad_T[2] /* *NOPAD* */ 46088b783a1SJames Wright }; 46188b783a1SJames Wright // Pressure 46288b783a1SJames Wright const CeedScalar 46388b783a1SJames Wright E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 46488626eedSJames Wright E_potential = -rho*(g[0]*x_i[0] + g[1]*x_i[1] + g[2]*x_i[2]), 46588b783a1SJames Wright E_internal = E - E_kinetic - E_potential, 46688b783a1SJames Wright P = E_internal * (gamma - 1.); // P = pressure 46788b783a1SJames Wright 46888b783a1SJames Wright // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 46988b783a1SJames Wright CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 47088626eedSJames Wright computeFluxJacobian_NS(jacob_F_conv, rho, u, E, gamma, g, x_i); 47188b783a1SJames Wright 47288b783a1SJames Wright // dqdx collects drhodx, dUdx and dEdx in one vector 47388b783a1SJames Wright CeedScalar dqdx[5][3]; 47488b783a1SJames Wright for (int j=0; j<3; j++) { 47588b783a1SJames Wright dqdx[0][j] = drhodx[j]; 47688b783a1SJames Wright dqdx[4][j] = dEdx[j]; 47788b783a1SJames Wright for (int k=0; k<3; k++) 47888b783a1SJames Wright dqdx[k+1][j] = dUdx[k][j]; 47988b783a1SJames Wright } 48088b783a1SJames Wright 48188b783a1SJames Wright // strong_conv = dF/dq * dq/dx (Strong convection) 48288b783a1SJames Wright CeedScalar strong_conv[5] = {0}; 48388b783a1SJames Wright for (int j=0; j<3; j++) 48488b783a1SJames Wright for (int k=0; k<5; k++) 48588b783a1SJames Wright for (int l=0; l<5; l++) 48688b783a1SJames Wright strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 48788b783a1SJames Wright 48888b783a1SJames Wright // Body force 48988626eedSJames Wright const CeedScalar body_force[5] = {0, rho *g[0], rho *g[1], rho *g[2], 0}; 49088b783a1SJames Wright 49188b783a1SJames Wright // The Physics 49288b783a1SJames Wright // Zero dv so all future terms can safely sum into it 49388b783a1SJames Wright for (int j=0; j<5; j++) 49488b783a1SJames Wright for (int k=0; k<3; k++) 49588b783a1SJames Wright dv[k][j][i] = 0; 49688b783a1SJames Wright 49788b783a1SJames Wright // -- Density 49888b783a1SJames Wright // ---- u rho 49988b783a1SJames Wright for (int j=0; j<3; j++) 50088b783a1SJames Wright dv[j][0][i] += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 50188b783a1SJames Wright rho*u[2]*dXdx[j][2]); 50288b783a1SJames Wright // -- Momentum 50388b783a1SJames Wright // ---- rho (u x u) + P I3 50488b783a1SJames Wright for (int j=0; j<3; j++) 50588b783a1SJames Wright for (int k=0; k<3; k++) 50688b783a1SJames Wright dv[k][j+1][i] += wdetJ*((rho*u[j]*u[0] + (j==0?P:0))*dXdx[k][0] + 50788b783a1SJames Wright (rho*u[j]*u[1] + (j==1?P:0))*dXdx[k][1] + 50888b783a1SJames Wright (rho*u[j]*u[2] + (j==2?P:0))*dXdx[k][2]); 50988b783a1SJames Wright // ---- Fuvisc 51088b783a1SJames Wright const CeedInt Fuviscidx[3][3] = {{0, 1, 2}, {1, 3, 4}, {2, 4, 5}}; // symmetric matrix indices 51188b783a1SJames Wright for (int j=0; j<3; j++) 51288b783a1SJames Wright for (int k=0; k<3; k++) 51388b783a1SJames Wright dv[k][j+1][i] -= wdetJ*(Fu[Fuviscidx[j][0]]*dXdx[k][0] + 51488b783a1SJames Wright Fu[Fuviscidx[j][1]]*dXdx[k][1] + 51588b783a1SJames Wright Fu[Fuviscidx[j][2]]*dXdx[k][2]); 51688b783a1SJames Wright // -- Total Energy Density 51788b783a1SJames Wright // ---- (E + P) u 51888b783a1SJames Wright for (int j=0; j<3; j++) 51988b783a1SJames Wright dv[j][4][i] += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 52088b783a1SJames Wright u[2]*dXdx[j][2]); 52188b783a1SJames Wright // ---- Fevisc 52288b783a1SJames Wright for (int j=0; j<3; j++) 52388b783a1SJames Wright dv[j][4][i] -= wdetJ * (Fe[0]*dXdx[j][0] + Fe[1]*dXdx[j][1] + 52488b783a1SJames Wright Fe[2]*dXdx[j][2]); 52588b783a1SJames Wright // Body Force 52688b783a1SJames Wright for (int j=0; j<5; j++) 52788b783a1SJames Wright v[j][i] = wdetJ * body_force[j]; 52888b783a1SJames Wright 52988626eedSJames Wright // Spatial Stabilization 53088626eedSJames Wright // -- Not used in favor of diagonal tau. Kept for future testing 53188626eedSJames Wright // const CeedScalar sound_speed = sqrt(gamma * P / rho); 53288626eedSJames Wright // CeedScalar Tau_x[3] = {0.}; 53388626eedSJames Wright // Tau_spatial(Tau_x, dXdx, u, sound_speed, context->c_tau, mu); 53488b783a1SJames Wright 53588626eedSJames Wright // -- Stabilization method: none, SU, or SUPG 53688626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 53788626eedSJames Wright CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0}; 53888626eedSJames Wright CeedScalar Tau_d[3] = {0.}; 53988b783a1SJames Wright switch (context->stabilization) { 54088b783a1SJames Wright case STAB_NONE: // Galerkin 54188b783a1SJames Wright break; 54288b783a1SJames Wright case STAB_SU: // SU 54388626eedSJames Wright Tau_diagPrim(Tau_d, dXdx, u, cv, context, mu, dt, rho); 54488626eedSJames Wright tau_strong_conv[0] = Tau_d[0] * strong_conv[0]; 54588626eedSJames Wright tau_strong_conv[1] = Tau_d[1] * strong_conv[1]; 54688626eedSJames Wright tau_strong_conv[2] = Tau_d[1] * strong_conv[2]; 54788626eedSJames Wright tau_strong_conv[3] = Tau_d[1] * strong_conv[3]; 54888626eedSJames Wright tau_strong_conv[4] = Tau_d[2] * strong_conv[4]; 54988626eedSJames Wright PrimitiveToConservative_fwd(rho, u, E, Rd, cv, tau_strong_conv, 55088626eedSJames Wright tau_strong_conv_conservative); 55188b783a1SJames Wright for (int j=0; j<3; j++) 55288b783a1SJames Wright for (int k=0; k<5; k++) 55388b783a1SJames Wright for (int l=0; l<5; l++) 55488626eedSJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l]; 55588b783a1SJames Wright 55688b783a1SJames Wright for (int j=0; j<5; j++) 55788b783a1SJames Wright for (int k=0; k<3; k++) 55888b783a1SJames Wright dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] + 55988b783a1SJames Wright stab[j][1] * dXdx[k][1] + 56088b783a1SJames Wright stab[j][2] * dXdx[k][2]); 56188b783a1SJames Wright break; 56288b783a1SJames Wright case STAB_SUPG: // SUPG is not implemented for explicit scheme 56388b783a1SJames Wright break; 56488b783a1SJames Wright } 56588b783a1SJames Wright 56688b783a1SJames Wright } // End Quadrature Point Loop 56788b783a1SJames Wright 56888b783a1SJames Wright // Return 56988b783a1SJames Wright return 0; 57088b783a1SJames Wright } 57188b783a1SJames Wright 57288b783a1SJames Wright // ***************************************************************************** 57388b783a1SJames Wright // This QFunction implements the Navier-Stokes equations (mentioned above) with 57488b783a1SJames Wright // implicit time stepping method 57588b783a1SJames Wright // 57688b783a1SJames Wright // SU = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) ) 57788b783a1SJames Wright // SUPG = Galerkin + grad(v) . ( Ai^T * Tau * (q_dot + Aj q,j - body force) ) 57888b783a1SJames Wright // (diffussive terms will be added later) 57988b783a1SJames Wright // 58088b783a1SJames Wright // ***************************************************************************** 58188b783a1SJames Wright CEED_QFUNCTION(IFunction_Newtonian)(void *ctx, CeedInt Q, 58288b783a1SJames Wright const CeedScalar *const *in, 58388b783a1SJames Wright CeedScalar *const *out) { 58488b783a1SJames Wright // *INDENT-OFF* 58588b783a1SJames Wright // Inputs 58688b783a1SJames Wright const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 58788b783a1SJames Wright (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 58888b783a1SJames Wright (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 58988b783a1SJames Wright (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3], 59088b783a1SJames Wright (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 59188b783a1SJames Wright // Outputs 59288b783a1SJames Wright CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 59388b783a1SJames Wright (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 59488b783a1SJames Wright // *INDENT-ON* 59588b783a1SJames Wright // Context 59688b783a1SJames Wright NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx; 59788b783a1SJames Wright const CeedScalar lambda = context->lambda; 59888b783a1SJames Wright const CeedScalar mu = context->mu; 59988b783a1SJames Wright const CeedScalar k = context->k; 60088b783a1SJames Wright const CeedScalar cv = context->cv; 60188b783a1SJames Wright const CeedScalar cp = context->cp; 60288626eedSJames Wright const CeedScalar *g = context->g; 60388626eedSJames Wright const CeedScalar dt = context->dt; 60488b783a1SJames Wright const CeedScalar gamma = cp / cv; 60588626eedSJames Wright const CeedScalar Rd = cp-cv; 60688b783a1SJames Wright 60788b783a1SJames Wright CeedPragmaSIMD 60888b783a1SJames Wright // Quadrature Point Loop 60988b783a1SJames Wright for (CeedInt i=0; i<Q; i++) { 61088b783a1SJames Wright // Setup 61188b783a1SJames Wright // -- Interp in 61288b783a1SJames Wright const CeedScalar rho = q[0][i]; 61388b783a1SJames Wright const CeedScalar u[3] = {q[1][i] / rho, 61488b783a1SJames Wright q[2][i] / rho, 61588b783a1SJames Wright q[3][i] / rho 61688b783a1SJames Wright }; 61788b783a1SJames Wright const CeedScalar E = q[4][i]; 61888b783a1SJames Wright // -- Grad in 61988b783a1SJames Wright const CeedScalar drho[3] = {dq[0][0][i], 62088b783a1SJames Wright dq[1][0][i], 62188b783a1SJames Wright dq[2][0][i] 62288b783a1SJames Wright }; 62388b783a1SJames Wright // *INDENT-OFF* 62488b783a1SJames Wright const CeedScalar dU[3][3] = {{dq[0][1][i], 62588b783a1SJames Wright dq[1][1][i], 62688b783a1SJames Wright dq[2][1][i]}, 62788b783a1SJames Wright {dq[0][2][i], 62888b783a1SJames Wright dq[1][2][i], 62988b783a1SJames Wright dq[2][2][i]}, 63088b783a1SJames Wright {dq[0][3][i], 63188b783a1SJames Wright dq[1][3][i], 63288b783a1SJames Wright dq[2][3][i]} 63388b783a1SJames Wright }; 63488b783a1SJames Wright // *INDENT-ON* 63588b783a1SJames Wright const CeedScalar dE[3] = {dq[0][4][i], 63688b783a1SJames Wright dq[1][4][i], 63788b783a1SJames Wright dq[2][4][i] 63888b783a1SJames Wright }; 63988b783a1SJames Wright // -- Interp-to-Interp q_data 64088b783a1SJames Wright const CeedScalar wdetJ = q_data[0][i]; 64188b783a1SJames Wright // -- Interp-to-Grad q_data 64288b783a1SJames Wright // ---- Inverse of change of coordinate matrix: X_i,j 64388b783a1SJames Wright // *INDENT-OFF* 64488b783a1SJames Wright const CeedScalar dXdx[3][3] = {{q_data[1][i], 64588b783a1SJames Wright q_data[2][i], 64688b783a1SJames Wright q_data[3][i]}, 64788b783a1SJames Wright {q_data[4][i], 64888b783a1SJames Wright q_data[5][i], 64988b783a1SJames Wright q_data[6][i]}, 65088b783a1SJames Wright {q_data[7][i], 65188b783a1SJames Wright q_data[8][i], 65288b783a1SJames Wright q_data[9][i]} 65388b783a1SJames Wright }; 65488626eedSJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 65588b783a1SJames Wright // *INDENT-ON* 65688b783a1SJames Wright // -- Grad-to-Grad q_data 65788b783a1SJames Wright // dU/dx 65888b783a1SJames Wright CeedScalar du[3][3] = {{0}}; 65988b783a1SJames Wright CeedScalar drhodx[3] = {0}; 66088b783a1SJames Wright CeedScalar dEdx[3] = {0}; 66188b783a1SJames Wright CeedScalar dUdx[3][3] = {{0}}; 66288b783a1SJames Wright CeedScalar dXdxdXdxT[3][3] = {{0}}; 66388b783a1SJames Wright for (int j=0; j<3; j++) { 66488b783a1SJames Wright for (int k=0; k<3; k++) { 66588b783a1SJames Wright du[j][k] = (dU[j][k] - drho[k]*u[j]) / rho; 66688b783a1SJames Wright drhodx[j] += drho[k] * dXdx[k][j]; 66788b783a1SJames Wright dEdx[j] += dE[k] * dXdx[k][j]; 66888b783a1SJames Wright for (int l=0; l<3; l++) { 66988b783a1SJames Wright dUdx[j][k] += dU[j][l] * dXdx[l][k]; 67088b783a1SJames Wright dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 67188b783a1SJames Wright } 67288b783a1SJames Wright } 67388b783a1SJames Wright } 67488b783a1SJames Wright CeedScalar dudx[3][3] = {{0}}; 67588b783a1SJames Wright for (int j=0; j<3; j++) 67688b783a1SJames Wright for (int k=0; k<3; k++) 67788b783a1SJames Wright for (int l=0; l<3; l++) 67888b783a1SJames Wright dudx[j][k] += du[j][l] * dXdx[l][k]; 67988b783a1SJames Wright // -- grad_T 68088b783a1SJames Wright const CeedScalar grad_T[3] = {(dEdx[0]/rho - E*drhodx[0]/(rho*rho) - /* *NOPAD* */ 68188626eedSJames Wright (u[0]*dudx[0][0] + u[1]*dudx[1][0] + u[2]*dudx[2][0]) + g[0])/cv, 68288b783a1SJames Wright (dEdx[1]/rho - E*drhodx[1]/(rho*rho) - /* *NOPAD* */ 68388626eedSJames Wright (u[0]*dudx[0][1] + u[1]*dudx[1][1] + u[2]*dudx[2][1]) + g[1])/cv, 68488b783a1SJames Wright (dEdx[2]/rho - E*drhodx[2]/(rho*rho) - /* *NOPAD* */ 68588626eedSJames Wright (u[0]*dudx[0][2] + u[1]*dudx[1][2] + u[2]*dudx[2][2]) + g[2])/cv 68688b783a1SJames Wright }; 68788b783a1SJames Wright // -- Fuvisc 68888b783a1SJames Wright // ---- Symmetric 3x3 matrix 68988b783a1SJames Wright const CeedScalar Fu[6] = {mu*(dudx[0][0] * (2 + lambda) + /* *NOPAD* */ 69088b783a1SJames Wright lambda * (dudx[1][1] + dudx[2][2])), 69188b783a1SJames Wright mu*(dudx[0][1] + dudx[1][0]), /* *NOPAD* */ 69288b783a1SJames Wright mu*(dudx[0][2] + dudx[2][0]), /* *NOPAD* */ 69388b783a1SJames Wright mu*(dudx[1][1] * (2 + lambda) + /* *NOPAD* */ 69488b783a1SJames Wright lambda * (dudx[0][0] + dudx[2][2])), 69588b783a1SJames Wright mu*(dudx[1][2] + dudx[2][1]), /* *NOPAD* */ 69688b783a1SJames Wright mu*(dudx[2][2] * (2 + lambda) + /* *NOPAD* */ 69788b783a1SJames Wright lambda * (dudx[0][0] + dudx[1][1])) 69888b783a1SJames Wright }; 69988b783a1SJames Wright // -- Fevisc 70088b783a1SJames Wright const CeedScalar Fe[3] = {u[0]*Fu[0] + u[1]*Fu[1] + u[2]*Fu[2] + /* *NOPAD* */ 70188b783a1SJames Wright k*grad_T[0], /* *NOPAD* */ 70288b783a1SJames Wright u[0]*Fu[1] + u[1]*Fu[3] + u[2]*Fu[4] + /* *NOPAD* */ 70388b783a1SJames Wright k*grad_T[1], /* *NOPAD* */ 70488b783a1SJames Wright u[0]*Fu[2] + u[1]*Fu[4] + u[2]*Fu[5] + /* *NOPAD* */ 70588b783a1SJames Wright k*grad_T[2] /* *NOPAD* */ 70688b783a1SJames Wright }; 70788b783a1SJames Wright // Pressure 70888b783a1SJames Wright const CeedScalar 70988b783a1SJames Wright E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 71088626eedSJames Wright E_potential = -rho*(g[0]*x_i[0] + g[1]*x_i[1] + g[2]*x_i[2]), 71188b783a1SJames Wright E_internal = E - E_kinetic - E_potential, 71288b783a1SJames Wright P = E_internal * (gamma - 1.); // P = pressure 71388b783a1SJames Wright 71488b783a1SJames Wright // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 71588b783a1SJames Wright CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 71688626eedSJames Wright computeFluxJacobian_NS(jacob_F_conv, rho, u, E, gamma, g, x_i); 71788b783a1SJames Wright 71888b783a1SJames Wright // dqdx collects drhodx, dUdx and dEdx in one vector 71988b783a1SJames Wright CeedScalar dqdx[5][3]; 72088b783a1SJames Wright for (int j=0; j<3; j++) { 72188b783a1SJames Wright dqdx[0][j] = drhodx[j]; 72288b783a1SJames Wright dqdx[4][j] = dEdx[j]; 72388b783a1SJames Wright for (int k=0; k<3; k++) 72488b783a1SJames Wright dqdx[k+1][j] = dUdx[k][j]; 72588b783a1SJames Wright } 72688b783a1SJames Wright // strong_conv = dF/dq * dq/dx (Strong convection) 72788b783a1SJames Wright CeedScalar strong_conv[5] = {0}; 72888b783a1SJames Wright for (int j=0; j<3; j++) 72988b783a1SJames Wright for (int k=0; k<5; k++) 73088b783a1SJames Wright for (int l=0; l<5; l++) 73188b783a1SJames Wright strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 73288b783a1SJames Wright 73388b783a1SJames Wright // Body force 73488626eedSJames Wright const CeedScalar body_force[5] = {0, rho *g[0], rho *g[1], rho *g[2], 0}; 73588b783a1SJames Wright 73688b783a1SJames Wright // Strong residual 73788b783a1SJames Wright CeedScalar strong_res[5]; 73888b783a1SJames Wright for (int j=0; j<5; j++) 73988b783a1SJames Wright strong_res[j] = q_dot[j][i] + strong_conv[j] - body_force[j]; 74088b783a1SJames Wright 74188b783a1SJames Wright // The Physics 74288b783a1SJames Wright //-----mass matrix 74388b783a1SJames Wright for (int j=0; j<5; j++) 74488b783a1SJames Wright v[j][i] = wdetJ*q_dot[j][i]; 74588b783a1SJames Wright 74688b783a1SJames Wright // Zero dv so all future terms can safely sum into it 74788b783a1SJames Wright for (int j=0; j<5; j++) 74888b783a1SJames Wright for (int k=0; k<3; k++) 74988b783a1SJames Wright dv[k][j][i] = 0; 75088b783a1SJames Wright 75188b783a1SJames Wright // -- Density 75288b783a1SJames Wright // ---- u rho 75388b783a1SJames Wright for (int j=0; j<3; j++) 75488b783a1SJames Wright dv[j][0][i] -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 75588b783a1SJames Wright rho*u[2]*dXdx[j][2]); 75688b783a1SJames Wright // -- Momentum 75788b783a1SJames Wright // ---- rho (u x u) + P I3 75888b783a1SJames Wright for (int j=0; j<3; j++) 75988b783a1SJames Wright for (int k=0; k<3; k++) 76088b783a1SJames Wright dv[k][j+1][i] -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0))*dXdx[k][0] + 76188b783a1SJames Wright (rho*u[j]*u[1] + (j==1?P:0))*dXdx[k][1] + 76288b783a1SJames Wright (rho*u[j]*u[2] + (j==2?P:0))*dXdx[k][2]); 76388b783a1SJames Wright // ---- Fuvisc 76488b783a1SJames Wright const CeedInt Fuviscidx[3][3] = {{0, 1, 2}, {1, 3, 4}, {2, 4, 5}}; // symmetric matrix indices 76588b783a1SJames Wright for (int j=0; j<3; j++) 76688b783a1SJames Wright for (int k=0; k<3; k++) 76788b783a1SJames Wright dv[k][j+1][i] += wdetJ*(Fu[Fuviscidx[j][0]]*dXdx[k][0] + 76888b783a1SJames Wright Fu[Fuviscidx[j][1]]*dXdx[k][1] + 76988b783a1SJames Wright Fu[Fuviscidx[j][2]]*dXdx[k][2]); 77088b783a1SJames Wright // -- Total Energy Density 77188b783a1SJames Wright // ---- (E + P) u 77288b783a1SJames Wright for (int j=0; j<3; j++) 77388b783a1SJames Wright dv[j][4][i] -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 77488b783a1SJames Wright u[2]*dXdx[j][2]); 77588b783a1SJames Wright // ---- Fevisc 77688b783a1SJames Wright for (int j=0; j<3; j++) 77788b783a1SJames Wright dv[j][4][i] += wdetJ * (Fe[0]*dXdx[j][0] + Fe[1]*dXdx[j][1] + 77888b783a1SJames Wright Fe[2]*dXdx[j][2]); 77988b783a1SJames Wright // Body Force 78088b783a1SJames Wright for (int j=0; j<5; j++) 78188b783a1SJames Wright v[j][i] -= wdetJ*body_force[j]; 78288b783a1SJames Wright 78388626eedSJames Wright // Spatial Stabilization 78488626eedSJames Wright // -- Not used in favor of diagonal tau. Kept for future testing 78588626eedSJames Wright // const CeedScalar sound_speed = sqrt(gamma * P / rho); 78688626eedSJames Wright // CeedScalar Tau_x[3] = {0.}; 78788626eedSJames Wright // Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau, mu); 78888b783a1SJames Wright 78988b783a1SJames Wright // -- Stabilization method: none, SU, or SUPG 79088626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 79188626eedSJames Wright CeedScalar tau_strong_res[5] = {0.}, tau_strong_res_conservative[5] = {0}; 79288626eedSJames Wright CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0}; 79388626eedSJames Wright CeedScalar Tau_d[3] = {0.}; 79488b783a1SJames Wright switch (context->stabilization) { 79588b783a1SJames Wright case STAB_NONE: // Galerkin 79688b783a1SJames Wright break; 79788b783a1SJames Wright case STAB_SU: // SU 79888626eedSJames Wright Tau_diagPrim(Tau_d, dXdx, u, cv, context, mu, dt, rho); 79988626eedSJames Wright tau_strong_conv[0] = Tau_d[0] * strong_conv[0]; 80088626eedSJames Wright tau_strong_conv[1] = Tau_d[1] * strong_conv[1]; 80188626eedSJames Wright tau_strong_conv[2] = Tau_d[1] * strong_conv[2]; 80288626eedSJames Wright tau_strong_conv[3] = Tau_d[1] * strong_conv[3]; 80388626eedSJames Wright tau_strong_conv[4] = Tau_d[2] * strong_conv[4]; 80488626eedSJames Wright PrimitiveToConservative_fwd(rho, u, E, Rd, cv, tau_strong_conv, 80588626eedSJames Wright tau_strong_conv_conservative); 80688b783a1SJames Wright for (int j=0; j<3; j++) 80788b783a1SJames Wright for (int k=0; k<5; k++) 80888b783a1SJames Wright for (int l=0; l<5; l++) 80988626eedSJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l]; 81088b783a1SJames Wright 81188b783a1SJames Wright for (int j=0; j<5; j++) 81288b783a1SJames Wright for (int k=0; k<3; k++) 81388b783a1SJames Wright dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 81488b783a1SJames Wright stab[j][1] * dXdx[k][1] + 81588b783a1SJames Wright stab[j][2] * dXdx[k][2]); 81688b783a1SJames Wright break; 81788b783a1SJames Wright case STAB_SUPG: // SUPG 81888626eedSJames Wright Tau_diagPrim(Tau_d, dXdx, u, cv, context, mu, dt, rho); 81988626eedSJames Wright tau_strong_res[0] = Tau_d[0] * strong_res[0]; 82088626eedSJames Wright tau_strong_res[1] = Tau_d[1] * strong_res[1]; 82188626eedSJames Wright tau_strong_res[2] = Tau_d[1] * strong_res[2]; 82288626eedSJames Wright tau_strong_res[3] = Tau_d[1] * strong_res[3]; 82388626eedSJames Wright tau_strong_res[4] = Tau_d[2] * strong_res[4]; 82488626eedSJames Wright // Alternate route (useful later with primitive variable code) 82588626eedSJames Wright // this function was verified against PHASTA for as IC that was as close as possible 82688626eedSJames Wright // computeFluxJacobian_NSp(jacob_F_conv_p, rho, u, E, Rd, cv); 82788626eedSJames Wright // it has also been verified to compute a correct through the following 82888626eedSJames Wright // stab[k][j] += jacob_F_conv_p[j][k][l] * tau_strong_res[l] // flux Jacobian wrt primitive 82988626eedSJames Wright // applied in the triple loop below 83088626eedSJames Wright // However, it is more flops than using the existing Jacobian wrt q after q_{,Y} viz 83188626eedSJames Wright PrimitiveToConservative_fwd(rho, u, E, Rd, cv, tau_strong_res, 83288626eedSJames Wright tau_strong_res_conservative); 83388b783a1SJames Wright for (int j=0; j<3; j++) 83488b783a1SJames Wright for (int k=0; k<5; k++) 83588b783a1SJames Wright for (int l=0; l<5; l++) 83688626eedSJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_res_conservative[l]; 83788b783a1SJames Wright 83888b783a1SJames Wright for (int j=0; j<5; j++) 83988b783a1SJames Wright for (int k=0; k<3; k++) 84088b783a1SJames Wright dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 84188b783a1SJames Wright stab[j][1] * dXdx[k][1] + 84288b783a1SJames Wright stab[j][2] * dXdx[k][2]); 84388b783a1SJames Wright break; 84488b783a1SJames Wright } 84588b783a1SJames Wright 84688b783a1SJames Wright } // End Quadrature Point Loop 84788b783a1SJames Wright 84888b783a1SJames Wright // Return 84988b783a1SJames Wright return 0; 85088b783a1SJames Wright } 85188b783a1SJames Wright // ***************************************************************************** 85288b783a1SJames Wright #endif // newtonian_h 853