1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 33a8779fbSJames Wright // 4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 53a8779fbSJames Wright // 6727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 73a8779fbSJames Wright 83a8779fbSJames Wright /// @file 93a8779fbSJames Wright /// Operator for Navier-Stokes example using PETSc 103a8779fbSJames Wright 113a8779fbSJames Wright 123a8779fbSJames Wright #ifndef newtonian_h 133a8779fbSJames Wright #define newtonian_h 143a8779fbSJames Wright 153a8779fbSJames Wright #include <math.h> 163a8779fbSJames Wright #include <ceed.h> 173a8779fbSJames Wright 183a8779fbSJames Wright #ifndef M_PI 193a8779fbSJames Wright #define M_PI 3.14159265358979323846 203a8779fbSJames Wright #endif 213a8779fbSJames Wright 223a8779fbSJames Wright #ifndef setup_context_struct 233a8779fbSJames Wright #define setup_context_struct 243a8779fbSJames Wright typedef struct SetupContext_ *SetupContext; 253a8779fbSJames Wright struct SetupContext_ { 263a8779fbSJames Wright CeedScalar theta0; 273a8779fbSJames Wright CeedScalar thetaC; 283a8779fbSJames Wright CeedScalar P0; 293a8779fbSJames Wright CeedScalar N; 303a8779fbSJames Wright CeedScalar cv; 313a8779fbSJames Wright CeedScalar cp; 32*bb8a0c61SJames Wright CeedScalar g[3]; 333a8779fbSJames Wright CeedScalar rc; 343a8779fbSJames Wright CeedScalar lx; 353a8779fbSJames Wright CeedScalar ly; 363a8779fbSJames Wright CeedScalar lz; 373a8779fbSJames Wright CeedScalar center[3]; 383a8779fbSJames Wright CeedScalar dc_axis[3]; 393a8779fbSJames Wright CeedScalar wind[3]; 403a8779fbSJames Wright CeedScalar time; 413a8779fbSJames Wright int wind_type; // See WindType: 0=ROTATION, 1=TRANSLATION 423a8779fbSJames Wright int bubble_type; // See BubbleType: 0=SPHERE, 1=CYLINDER 433a8779fbSJames Wright int bubble_continuity_type; // See BubbleContinuityType: 0=SMOOTH, 1=BACK_SHARP 2=THICK 443a8779fbSJames Wright }; 453a8779fbSJames Wright #endif 463a8779fbSJames Wright 473a8779fbSJames Wright #ifndef newtonian_context_struct 483a8779fbSJames Wright #define newtonian_context_struct 493a8779fbSJames Wright typedef enum { 503a8779fbSJames Wright STAB_NONE = 0, 513a8779fbSJames Wright STAB_SU = 1, // Streamline Upwind 523a8779fbSJames Wright STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin 533a8779fbSJames Wright } StabilizationType; 543a8779fbSJames Wright 553a8779fbSJames Wright typedef struct NewtonianIdealGasContext_ *NewtonianIdealGasContext; 563a8779fbSJames Wright struct NewtonianIdealGasContext_ { 573a8779fbSJames Wright CeedScalar lambda; 583a8779fbSJames Wright CeedScalar mu; 593a8779fbSJames Wright CeedScalar k; 603a8779fbSJames Wright CeedScalar cv; 613a8779fbSJames Wright CeedScalar cp; 62*bb8a0c61SJames Wright CeedScalar g[3]; 633a8779fbSJames Wright CeedScalar c_tau; 64*bb8a0c61SJames Wright CeedScalar Ctau_t; 65*bb8a0c61SJames Wright CeedScalar Ctau_v; 66*bb8a0c61SJames Wright CeedScalar Ctau_C; 67*bb8a0c61SJames Wright CeedScalar Ctau_M; 68*bb8a0c61SJames Wright CeedScalar Ctau_E; 69*bb8a0c61SJames Wright CeedScalar dt; 703a8779fbSJames Wright StabilizationType stabilization; 713a8779fbSJames Wright }; 723a8779fbSJames Wright #endif 733a8779fbSJames Wright 743a8779fbSJames Wright // ***************************************************************************** 753a8779fbSJames Wright // Helper function for computing flux Jacobian 763a8779fbSJames Wright // ***************************************************************************** 773a8779fbSJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NS(CeedScalar dF[3][5][5], 783a8779fbSJames Wright const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 79*bb8a0c61SJames Wright const CeedScalar gamma, const CeedScalar g[3], const CeedScalar x[3]) { 803a8779fbSJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 81*bb8a0c61SJames Wright CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]); 823a8779fbSJames Wright for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 833a8779fbSJames Wright for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix 84*bb8a0c61SJames Wright dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2. - e_potential)) : 0.) - 85*bb8a0c61SJames Wright u[i]*u[j]; 863a8779fbSJames Wright for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix 873a8779fbSJames Wright dF[i][0][k+1] = ((i==k) ? 1. : 0.); 883a8779fbSJames Wright dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) + 893a8779fbSJames Wright ((i==k) ? u[j] : 0.) - 903a8779fbSJames Wright ((i==j) ? u[k] : 0.) * (gamma-1.); 913a8779fbSJames Wright dF[i][4][k+1] = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) - 923a8779fbSJames Wright (gamma-1.)*u[i]*u[k]; 933a8779fbSJames Wright } 943a8779fbSJames Wright dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.); 953a8779fbSJames Wright } 963a8779fbSJames Wright dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho); 973a8779fbSJames Wright dF[i][4][4] = u[i] * gamma; 983a8779fbSJames Wright } 993a8779fbSJames Wright } 1003a8779fbSJames Wright 1013a8779fbSJames Wright // ***************************************************************************** 102*bb8a0c61SJames Wright // Helper function for computing flux Jacobian of Primitive variables 103*bb8a0c61SJames Wright // ***************************************************************************** 104*bb8a0c61SJames Wright CEED_QFUNCTION_HELPER void computeFluxJacobian_NSp(CeedScalar dF[3][5][5], 105*bb8a0c61SJames Wright const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 106*bb8a0c61SJames Wright const CeedScalar Rd, const CeedScalar cv) { 107*bb8a0c61SJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 108*bb8a0c61SJames Wright // TODO Add in gravity's contribution 109*bb8a0c61SJames Wright 110*bb8a0c61SJames Wright CeedScalar T = ( E / rho - u_sq / 2. ) / cv; 111*bb8a0c61SJames Wright CeedScalar drdT = -rho / T; 112*bb8a0c61SJames Wright CeedScalar drdP = 1. / ( Rd * T); 113*bb8a0c61SJames Wright CeedScalar etot = E / rho ; 114*bb8a0c61SJames Wright CeedScalar e2p = drdP * etot + 1. ; 115*bb8a0c61SJames Wright CeedScalar e3p = ( E + rho * Rd * T ); 116*bb8a0c61SJames Wright CeedScalar e4p = drdT * etot + rho * cv ; 117*bb8a0c61SJames Wright 118*bb8a0c61SJames Wright for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 119*bb8a0c61SJames Wright for (CeedInt j=0; j<3; j++) { // j counts F^{m_j} 120*bb8a0c61SJames Wright // [row][col] of A_i 121*bb8a0c61SJames Wright dF[i][j+1][0] = drdP * u[i] * u[j] + ((i==j) ? 1. : 0.); // F^{{m_j} wrt p 122*bb8a0c61SJames Wright for (CeedInt k=0; k<3; k++) { // k counts the wrt vel_k 123*bb8a0c61SJames Wright // this loop handles middle columns for all 5 rows 124*bb8a0c61SJames Wright dF[i][0][k+1] = ((i==k) ? rho : 0.); // F^c wrt vel_k 125*bb8a0c61SJames Wright dF[i][j+1][k+1] = (((j==k) ? u[i] : 0.) + // F^m_j wrt u_k 126*bb8a0c61SJames Wright ((i==k) ? u[j] : 0.) ) * rho; 127*bb8a0c61SJames Wright dF[i][4][k+1] = rho * u[i] * u[k] 128*bb8a0c61SJames Wright + ((i==k) ? e3p : 0.) ; // F^e wrt u_k 129*bb8a0c61SJames Wright } 130*bb8a0c61SJames Wright dF[i][j+1][4] = drdT * u[i] * u[j]; // F^{m_j} wrt T 131*bb8a0c61SJames Wright } 132*bb8a0c61SJames Wright dF[i][4][0] = u[i] * e2p; // F^e wrt p 133*bb8a0c61SJames Wright dF[i][4][4] = u[i] * e4p; // F^e wrt T 134*bb8a0c61SJames Wright dF[i][0][0] = u[i] * drdP; // F^c wrt p 135*bb8a0c61SJames Wright dF[i][0][4] = u[i] * drdT; // F^c wrt T 136*bb8a0c61SJames Wright } 137*bb8a0c61SJames Wright } 138*bb8a0c61SJames Wright 139*bb8a0c61SJames Wright CEED_QFUNCTION_HELPER void PrimitiveToConservative_fwd(const CeedScalar rho, 140*bb8a0c61SJames Wright const CeedScalar u[3], const CeedScalar E, const CeedScalar Rd, 141*bb8a0c61SJames Wright const CeedScalar cv, const CeedScalar dY[5], CeedScalar dU[5]) { 142*bb8a0c61SJames Wright CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; 143*bb8a0c61SJames Wright CeedScalar T = ( E / rho - u_sq / 2. ) / cv; 144*bb8a0c61SJames Wright CeedScalar drdT = -rho / T; 145*bb8a0c61SJames Wright CeedScalar drdP = 1. / ( Rd * T); 146*bb8a0c61SJames Wright dU[0] = drdP * dY[0] + drdT * dY[4]; 147*bb8a0c61SJames Wright CeedScalar de_kinetic = 0; 148*bb8a0c61SJames Wright for (int i=0; i<3; i++) { 149*bb8a0c61SJames Wright dU[1+i] = dU[0] * u[i] + rho * dY[1+i]; 150*bb8a0c61SJames Wright de_kinetic += u[i] * dY[1+i]; 151*bb8a0c61SJames Wright } 152*bb8a0c61SJames Wright dU[4] = rho * cv * dY[4] + dU[0] * cv * T // internal energy: rho * e 153*bb8a0c61SJames Wright + rho * de_kinetic + .5 * dU[0] * u_sq; // kinetic energy: .5 * rho * |u|^2 154*bb8a0c61SJames Wright } 155*bb8a0c61SJames Wright 156*bb8a0c61SJames Wright // ***************************************************************************** 157*bb8a0c61SJames Wright // Helper function for computing Tau elements (stabilization constant) 158*bb8a0c61SJames Wright // Model from: 159*bb8a0c61SJames Wright // PHASTA 160*bb8a0c61SJames Wright // 161*bb8a0c61SJames Wright // Tau[i] = itau=0 which is diagonal-Shakib (3 values still but not spatial) 162*bb8a0c61SJames Wright // 163*bb8a0c61SJames Wright // Where NOT UPDATED YET 164*bb8a0c61SJames Wright // ***************************************************************************** 165*bb8a0c61SJames Wright CEED_QFUNCTION_HELPER void Tau_diagPrim(CeedScalar Tau_d[3], 166*bb8a0c61SJames Wright const CeedScalar dXdx[3][3], const CeedScalar u[3], 167*bb8a0c61SJames Wright const CeedScalar cv, const NewtonianIdealGasContext newt_ctx, 168*bb8a0c61SJames Wright const CeedScalar mu, const CeedScalar dt, 169*bb8a0c61SJames Wright const CeedScalar rho) { 170*bb8a0c61SJames Wright // Context 171*bb8a0c61SJames Wright const CeedScalar Ctau_t = newt_ctx->Ctau_t; 172*bb8a0c61SJames Wright const CeedScalar Ctau_v = newt_ctx->Ctau_v; 173*bb8a0c61SJames Wright const CeedScalar Ctau_C = newt_ctx->Ctau_C; 174*bb8a0c61SJames Wright const CeedScalar Ctau_M = newt_ctx->Ctau_M; 175*bb8a0c61SJames Wright const CeedScalar Ctau_E = newt_ctx->Ctau_E; 176*bb8a0c61SJames Wright CeedScalar gijd[6]; 177*bb8a0c61SJames Wright CeedScalar tau; 178*bb8a0c61SJames Wright CeedScalar dts; 179*bb8a0c61SJames Wright CeedScalar fact; 180*bb8a0c61SJames Wright 181*bb8a0c61SJames Wright //*INDENT-OFF* 182*bb8a0c61SJames Wright gijd[0] = dXdx[0][0] * dXdx[0][0] 183*bb8a0c61SJames Wright + dXdx[1][0] * dXdx[1][0] 184*bb8a0c61SJames Wright + dXdx[2][0] * dXdx[2][0]; 185*bb8a0c61SJames Wright 186*bb8a0c61SJames Wright gijd[1] = dXdx[0][0] * dXdx[0][1] 187*bb8a0c61SJames Wright + dXdx[1][0] * dXdx[1][1] 188*bb8a0c61SJames Wright + dXdx[2][0] * dXdx[2][1]; 189*bb8a0c61SJames Wright 190*bb8a0c61SJames Wright gijd[2] = dXdx[0][1] * dXdx[0][1] 191*bb8a0c61SJames Wright + dXdx[1][1] * dXdx[1][1] 192*bb8a0c61SJames Wright + dXdx[2][1] * dXdx[2][1]; 193*bb8a0c61SJames Wright 194*bb8a0c61SJames Wright gijd[3] = dXdx[0][0] * dXdx[0][2] 195*bb8a0c61SJames Wright + dXdx[1][0] * dXdx[1][2] 196*bb8a0c61SJames Wright + dXdx[2][0] * dXdx[2][2]; 197*bb8a0c61SJames Wright 198*bb8a0c61SJames Wright gijd[4] = dXdx[0][1] * dXdx[0][2] 199*bb8a0c61SJames Wright + dXdx[1][1] * dXdx[1][2] 200*bb8a0c61SJames Wright + dXdx[2][1] * dXdx[2][2]; 201*bb8a0c61SJames Wright 202*bb8a0c61SJames Wright gijd[5] = dXdx[0][2] * dXdx[0][2] 203*bb8a0c61SJames Wright + dXdx[1][2] * dXdx[1][2] 204*bb8a0c61SJames Wright + dXdx[2][2] * dXdx[2][2]; 205*bb8a0c61SJames Wright //*INDENT-ON* 206*bb8a0c61SJames Wright 207*bb8a0c61SJames Wright dts = Ctau_t / dt ; 208*bb8a0c61SJames Wright 209*bb8a0c61SJames Wright tau = rho*rho*((4. * dts * dts) 210*bb8a0c61SJames Wright + u[0] * ( u[0] * gijd[0] + 2. * ( u[1] * gijd[1] + u[2] * gijd[3])) 211*bb8a0c61SJames Wright + u[1] * ( u[1] * gijd[2] + 2. * u[2] * gijd[4]) 212*bb8a0c61SJames Wright + u[2] * u[2] * gijd[5]) 213*bb8a0c61SJames Wright + Ctau_v* mu * mu * 214*bb8a0c61SJames Wright (gijd[0]*gijd[0] + gijd[2]*gijd[2] + gijd[5]*gijd[5] + 215*bb8a0c61SJames Wright + 2. * (gijd[1]*gijd[1] + gijd[3]*gijd[3] + gijd[4]*gijd[4])); 216*bb8a0c61SJames Wright 217*bb8a0c61SJames Wright fact=sqrt(tau); 218*bb8a0c61SJames Wright 219*bb8a0c61SJames Wright Tau_d[0] = Ctau_C * fact / (rho*(gijd[0] + gijd[2] + gijd[5]))*0.125; 220*bb8a0c61SJames Wright 221*bb8a0c61SJames Wright Tau_d[1] = Ctau_M / fact; 222*bb8a0c61SJames Wright Tau_d[2] = Ctau_E / ( fact * cv ); 223*bb8a0c61SJames Wright 224*bb8a0c61SJames Wright // consider putting back the way I initially had it Ctau_E * Tau_d[1] /cv 225*bb8a0c61SJames Wright // to avoid a division if the compiler is smart enough to see that cv IS 226*bb8a0c61SJames Wright // a constant that it could invert once for all elements 227*bb8a0c61SJames Wright // but in that case energy tau is scaled by the product of Ctau_E * Ctau_M 228*bb8a0c61SJames Wright // OR we could absorb cv into Ctau_E but this puts more burden on user to 229*bb8a0c61SJames Wright // know how to change constants with a change of fluid or units. Same for 230*bb8a0c61SJames Wright // Ctau_v * mu * mu IF AND ONLY IF we don't add viscosity law =f(T) 231*bb8a0c61SJames Wright } 232*bb8a0c61SJames Wright 233*bb8a0c61SJames Wright // ***************************************************************************** 2343a8779fbSJames Wright // Helper function for computing Tau elements (stabilization constant) 2353a8779fbSJames Wright // Model from: 2363a8779fbSJames Wright // Stabilized Methods for Compressible Flows, Hughes et al 2010 2373a8779fbSJames Wright // 2383a8779fbSJames Wright // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 2393a8779fbSJames Wright // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 2403a8779fbSJames Wright // 2413a8779fbSJames Wright // Where 2423a8779fbSJames Wright // c_tau = stabilization constant (0.5 is reported as "optimal") 2433a8779fbSJames Wright // h[i] = 2 length(dxdX[i]) 2443a8779fbSJames Wright // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 2453a8779fbSJames Wright // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 2463a8779fbSJames Wright // rho(A[i]) = spectral radius of the convective flux Jacobian i, 2473a8779fbSJames Wright // wave speed in direction i 2483a8779fbSJames Wright // ***************************************************************************** 2493a8779fbSJames Wright CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], 2503a8779fbSJames Wright const CeedScalar dXdx[3][3], const CeedScalar u[3], 251*bb8a0c61SJames Wright /* const CeedScalar sound_speed, const CeedScalar c_tau) { */ 252*bb8a0c61SJames Wright const CeedScalar sound_speed, const CeedScalar c_tau, 253*bb8a0c61SJames Wright const CeedScalar viscosity) { 254*bb8a0c61SJames Wright const CeedScalar mag_u_visc = sqrt(u[0]*u[0] +u[1]*u[1] +u[2]*u[2]) / 255*bb8a0c61SJames Wright (2*viscosity); 2563a8779fbSJames Wright for (int i=0; i<3; i++) { 2573a8779fbSJames Wright // length of element in direction i 2583a8779fbSJames Wright CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] + 2593a8779fbSJames Wright dXdx[2][i]*dXdx[2][i]); 260*bb8a0c61SJames Wright CeedScalar Pe = mag_u_visc*h; 261*bb8a0c61SJames Wright CeedScalar Xi = 1/tanh(Pe) - 1/Pe; 2623a8779fbSJames Wright // fastest wave in direction i 2633a8779fbSJames Wright CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 264*bb8a0c61SJames Wright Tau_x[i] = c_tau * h * Xi / fastest_wave; 2653a8779fbSJames Wright } 2663a8779fbSJames Wright } 2673a8779fbSJames Wright 2683a8779fbSJames Wright // ***************************************************************************** 2693a8779fbSJames Wright // This QFunction sets a "still" initial condition for generic Newtonian IG problems 2703a8779fbSJames Wright // ***************************************************************************** 2713a8779fbSJames Wright CEED_QFUNCTION(ICsNewtonianIG)(void *ctx, CeedInt Q, 2723a8779fbSJames Wright const CeedScalar *const *in, CeedScalar *const *out) { 2733a8779fbSJames Wright // Inputs 2743a8779fbSJames Wright const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 2753a8779fbSJames Wright 2763a8779fbSJames Wright // Outputs 2773a8779fbSJames Wright CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 2783a8779fbSJames Wright 279*bb8a0c61SJames Wright // Context 280*bb8a0c61SJames Wright const SetupContext context = (SetupContext)ctx; 281*bb8a0c61SJames Wright const CeedScalar theta0 = context->theta0; 282*bb8a0c61SJames Wright const CeedScalar P0 = context->P0; 283*bb8a0c61SJames Wright const CeedScalar cv = context->cv; 284*bb8a0c61SJames Wright const CeedScalar cp = context->cp; 285*bb8a0c61SJames Wright const CeedScalar *g = context->g; 286*bb8a0c61SJames Wright const CeedScalar Rd = cp - cv; 287*bb8a0c61SJames Wright 2883a8779fbSJames Wright // Quadrature Point Loop 2893a8779fbSJames Wright CeedPragmaSIMD 2903a8779fbSJames Wright for (CeedInt i=0; i<Q; i++) { 2913a8779fbSJames Wright CeedScalar q[5] = {0.}; 2923a8779fbSJames Wright 2933a8779fbSJames Wright // Setup 2943a8779fbSJames Wright // -- Coordinates 295*bb8a0c61SJames Wright const CeedScalar x[3] = {X[0][i], X[1][i], X[2][i]}; 296*bb8a0c61SJames Wright const CeedScalar e_potential = -(g[0]*x[0] + g[1]*x[1] + g[2]*x[2]); 2973a8779fbSJames Wright 2983a8779fbSJames Wright // -- Density 299*bb8a0c61SJames Wright const CeedScalar rho = P0 / (Rd*theta0); 3003a8779fbSJames Wright 3013a8779fbSJames Wright // Initial Conditions 3023a8779fbSJames Wright q[0] = rho; 3033a8779fbSJames Wright q[1] = 0.0; 3043a8779fbSJames Wright q[2] = 0.0; 3053a8779fbSJames Wright q[3] = 0.0; 306*bb8a0c61SJames Wright q[4] = rho * (cv*theta0 + e_potential); 3073a8779fbSJames Wright 3083a8779fbSJames Wright for (CeedInt j=0; j<5; j++) 3093a8779fbSJames Wright q0[j][i] = q[j]; 3103a8779fbSJames Wright } // End of Quadrature Point Loop 3113a8779fbSJames Wright return 0; 3123a8779fbSJames Wright } 3133a8779fbSJames Wright 3143a8779fbSJames Wright // ***************************************************************************** 3153a8779fbSJames Wright // This QFunction implements the following formulation of Navier-Stokes with 3163a8779fbSJames Wright // explicit time stepping method 3173a8779fbSJames Wright // 3183a8779fbSJames Wright // This is 3D compressible Navier-Stokes in conservation form with state 3193a8779fbSJames Wright // variables of density, momentum density, and total energy density. 3203a8779fbSJames Wright // 3213a8779fbSJames Wright // State Variables: q = ( rho, U1, U2, U3, E ) 3223a8779fbSJames Wright // rho - Mass Density 3233a8779fbSJames Wright // Ui - Momentum Density, Ui = rho ui 3243a8779fbSJames Wright // E - Total Energy Density, E = rho (cv T + (u u)/2 + g z) 3253a8779fbSJames Wright // 3263a8779fbSJames Wright // Navier-Stokes Equations: 3273a8779fbSJames Wright // drho/dt + div( U ) = 0 3283a8779fbSJames Wright // dU/dt + div( rho (u x u) + P I3 ) + rho g khat = div( Fu ) 3293a8779fbSJames Wright // dE/dt + div( (E + P) u ) = div( Fe ) 3303a8779fbSJames Wright // 3313a8779fbSJames Wright // Viscous Stress: 3323a8779fbSJames Wright // Fu = mu (grad( u ) + grad( u )^T + lambda div ( u ) I3) 3333a8779fbSJames Wright // 3343a8779fbSJames Wright // Thermal Stress: 3353a8779fbSJames Wright // Fe = u Fu + k grad( T ) 336*bb8a0c61SJames Wright // Equation of State 3373a8779fbSJames Wright // P = (gamma - 1) (E - rho (u u) / 2 - rho g z) 3383a8779fbSJames Wright // 3393a8779fbSJames Wright // Stabilization: 3403a8779fbSJames Wright // Tau = diag(TauC, TauM, TauM, TauM, TauE) 3413a8779fbSJames Wright // f1 = rho sqrt(ui uj gij) 3423a8779fbSJames Wright // gij = dXi/dX * dXi/dX 3433a8779fbSJames Wright // TauC = Cc f1 / (8 gii) 3443a8779fbSJames Wright // TauM = min( 1 , 1 / f1 ) 3453a8779fbSJames Wright // TauE = TauM / (Ce cv) 3463a8779fbSJames Wright // 3473a8779fbSJames Wright // SU = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) ) 3483a8779fbSJames Wright // 3493a8779fbSJames Wright // Constants: 3503a8779fbSJames Wright // lambda = - 2 / 3, From Stokes hypothesis 3513a8779fbSJames Wright // mu , Dynamic viscosity 3523a8779fbSJames Wright // k , Thermal conductivity 3533a8779fbSJames Wright // cv , Specific heat, constant volume 3543a8779fbSJames Wright // cp , Specific heat, constant pressure 3553a8779fbSJames Wright // g , Gravity 3563a8779fbSJames Wright // gamma = cp / cv, Specific heat ratio 3573a8779fbSJames Wright // 3583a8779fbSJames Wright // We require the product of the inverse of the Jacobian (dXdx_j,k) and 3593a8779fbSJames Wright // its transpose (dXdx_k,j) to properly compute integrals of the form: 3603a8779fbSJames Wright // int( gradv gradu ) 3613a8779fbSJames Wright // 3623a8779fbSJames Wright // ***************************************************************************** 3633a8779fbSJames Wright CEED_QFUNCTION(Newtonian)(void *ctx, CeedInt Q, 3643a8779fbSJames Wright const CeedScalar *const *in, CeedScalar *const *out) { 3653a8779fbSJames Wright // *INDENT-OFF* 3663a8779fbSJames Wright // Inputs 3673a8779fbSJames Wright const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 3683a8779fbSJames Wright (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 3693a8779fbSJames Wright (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 3703a8779fbSJames Wright (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 3713a8779fbSJames Wright // Outputs 3723a8779fbSJames Wright CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 3733a8779fbSJames Wright (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 3743a8779fbSJames Wright // *INDENT-ON* 3753a8779fbSJames Wright 3763a8779fbSJames Wright // Context 3773a8779fbSJames Wright NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx; 3783a8779fbSJames Wright const CeedScalar lambda = context->lambda; 3793a8779fbSJames Wright const CeedScalar mu = context->mu; 3803a8779fbSJames Wright const CeedScalar k = context->k; 3813a8779fbSJames Wright const CeedScalar cv = context->cv; 3823a8779fbSJames Wright const CeedScalar cp = context->cp; 383*bb8a0c61SJames Wright const CeedScalar *g = context->g; 384*bb8a0c61SJames Wright const CeedScalar dt = context->dt; 3853a8779fbSJames Wright const CeedScalar gamma = cp / cv; 386*bb8a0c61SJames Wright const CeedScalar Rd = cp - cv; 3873a8779fbSJames Wright 3883a8779fbSJames Wright CeedPragmaSIMD 3893a8779fbSJames Wright // Quadrature Point Loop 3903a8779fbSJames Wright for (CeedInt i=0; i<Q; i++) { 3913a8779fbSJames Wright // *INDENT-OFF* 3923a8779fbSJames Wright // Setup 3933a8779fbSJames Wright // -- Interp in 3943a8779fbSJames Wright const CeedScalar rho = q[0][i]; 3953a8779fbSJames Wright const CeedScalar u[3] = {q[1][i] / rho, 3963a8779fbSJames Wright q[2][i] / rho, 3973a8779fbSJames Wright q[3][i] / rho 3983a8779fbSJames Wright }; 3993a8779fbSJames Wright const CeedScalar E = q[4][i]; 4003a8779fbSJames Wright // -- Grad in 4013a8779fbSJames Wright const CeedScalar drho[3] = {dq[0][0][i], 4023a8779fbSJames Wright dq[1][0][i], 4033a8779fbSJames Wright dq[2][0][i] 4043a8779fbSJames Wright }; 4053a8779fbSJames Wright const CeedScalar dU[3][3] = {{dq[0][1][i], 4063a8779fbSJames Wright dq[1][1][i], 4073a8779fbSJames Wright dq[2][1][i]}, 4083a8779fbSJames Wright {dq[0][2][i], 4093a8779fbSJames Wright dq[1][2][i], 4103a8779fbSJames Wright dq[2][2][i]}, 4113a8779fbSJames Wright {dq[0][3][i], 4123a8779fbSJames Wright dq[1][3][i], 4133a8779fbSJames Wright dq[2][3][i]} 4143a8779fbSJames Wright }; 4153a8779fbSJames Wright const CeedScalar dE[3] = {dq[0][4][i], 4163a8779fbSJames Wright dq[1][4][i], 4173a8779fbSJames Wright dq[2][4][i] 4183a8779fbSJames Wright }; 4193a8779fbSJames Wright // -- Interp-to-Interp q_data 4203a8779fbSJames Wright const CeedScalar wdetJ = q_data[0][i]; 4213a8779fbSJames Wright // -- Interp-to-Grad q_data 4223a8779fbSJames Wright // ---- Inverse of change of coordinate matrix: X_i,j 4233a8779fbSJames Wright // *INDENT-OFF* 4243a8779fbSJames Wright const CeedScalar dXdx[3][3] = {{q_data[1][i], 4253a8779fbSJames Wright q_data[2][i], 4263a8779fbSJames Wright q_data[3][i]}, 4273a8779fbSJames Wright {q_data[4][i], 4283a8779fbSJames Wright q_data[5][i], 4293a8779fbSJames Wright q_data[6][i]}, 4303a8779fbSJames Wright {q_data[7][i], 4313a8779fbSJames Wright q_data[8][i], 4323a8779fbSJames Wright q_data[9][i]} 4333a8779fbSJames Wright }; 434*bb8a0c61SJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 4353a8779fbSJames Wright // *INDENT-ON* 4363a8779fbSJames Wright // -- Grad-to-Grad q_data 4373a8779fbSJames Wright // dU/dx 4383a8779fbSJames Wright CeedScalar du[3][3] = {{0}}; 4393a8779fbSJames Wright CeedScalar drhodx[3] = {0}; 4403a8779fbSJames Wright CeedScalar dEdx[3] = {0}; 4413a8779fbSJames Wright CeedScalar dUdx[3][3] = {{0}}; 4423a8779fbSJames Wright CeedScalar dXdxdXdxT[3][3] = {{0}}; 4433a8779fbSJames Wright for (int j=0; j<3; j++) { 4443a8779fbSJames Wright for (int k=0; k<3; k++) { 4453a8779fbSJames Wright du[j][k] = (dU[j][k] - drho[k]*u[j]) / rho; 4463a8779fbSJames Wright drhodx[j] += drho[k] * dXdx[k][j]; 4473a8779fbSJames Wright dEdx[j] += dE[k] * dXdx[k][j]; 4483a8779fbSJames Wright for (int l=0; l<3; l++) { 4493a8779fbSJames Wright dUdx[j][k] += dU[j][l] * dXdx[l][k]; 4503a8779fbSJames Wright dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 4513a8779fbSJames Wright } 4523a8779fbSJames Wright } 4533a8779fbSJames Wright } 4543a8779fbSJames Wright CeedScalar dudx[3][3] = {{0}}; 4553a8779fbSJames Wright for (int j=0; j<3; j++) 4563a8779fbSJames Wright for (int k=0; k<3; k++) 4573a8779fbSJames Wright for (int l=0; l<3; l++) 4583a8779fbSJames Wright dudx[j][k] += du[j][l] * dXdx[l][k]; 4593a8779fbSJames Wright // -- grad_T 4603a8779fbSJames Wright const CeedScalar grad_T[3] = {(dEdx[0]/rho - E*drhodx[0]/(rho*rho) - /* *NOPAD* */ 461*bb8a0c61SJames Wright (u[0]*dudx[0][0] + u[1]*dudx[1][0] + u[2]*dudx[2][0]) + g[0])/cv, 4623a8779fbSJames Wright (dEdx[1]/rho - E*drhodx[1]/(rho*rho) - /* *NOPAD* */ 463*bb8a0c61SJames Wright (u[0]*dudx[0][1] + u[1]*dudx[1][1] + u[2]*dudx[2][1]) + g[1])/cv, 4643a8779fbSJames Wright (dEdx[2]/rho - E*drhodx[2]/(rho*rho) - /* *NOPAD* */ 465*bb8a0c61SJames Wright (u[0]*dudx[0][2] + u[1]*dudx[1][2] + u[2]*dudx[2][2]) + g[2])/cv 4663a8779fbSJames Wright }; 4673a8779fbSJames Wright 4683a8779fbSJames Wright // -- Fuvisc 4693a8779fbSJames Wright // ---- Symmetric 3x3 matrix 4703a8779fbSJames Wright const CeedScalar Fu[6] = {mu*(dudx[0][0] * (2 + lambda) + /* *NOPAD* */ 4713a8779fbSJames Wright lambda * (dudx[1][1] + dudx[2][2])), 4723a8779fbSJames Wright mu*(dudx[0][1] + dudx[1][0]), /* *NOPAD* */ 4733a8779fbSJames Wright mu*(dudx[0][2] + dudx[2][0]), /* *NOPAD* */ 4743a8779fbSJames Wright mu*(dudx[1][1] * (2 + lambda) + /* *NOPAD* */ 4753a8779fbSJames Wright lambda * (dudx[0][0] + dudx[2][2])), 4763a8779fbSJames Wright mu*(dudx[1][2] + dudx[2][1]), /* *NOPAD* */ 4773a8779fbSJames Wright mu*(dudx[2][2] * (2 + lambda) + /* *NOPAD* */ 4783a8779fbSJames Wright lambda * (dudx[0][0] + dudx[1][1])) 4793a8779fbSJames Wright }; 4803a8779fbSJames Wright // -- Fevisc 4813a8779fbSJames Wright const CeedScalar Fe[3] = {u[0]*Fu[0] + u[1]*Fu[1] + u[2]*Fu[2] + /* *NOPAD* */ 4823a8779fbSJames Wright k*grad_T[0], /* *NOPAD* */ 4833a8779fbSJames Wright u[0]*Fu[1] + u[1]*Fu[3] + u[2]*Fu[4] + /* *NOPAD* */ 4843a8779fbSJames Wright k*grad_T[1], /* *NOPAD* */ 4853a8779fbSJames Wright u[0]*Fu[2] + u[1]*Fu[4] + u[2]*Fu[5] + /* *NOPAD* */ 4863a8779fbSJames Wright k*grad_T[2] /* *NOPAD* */ 4873a8779fbSJames Wright }; 4883a8779fbSJames Wright // Pressure 4893a8779fbSJames Wright const CeedScalar 4903a8779fbSJames Wright E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 491*bb8a0c61SJames Wright E_potential = -rho*(g[0]*x_i[0] + g[1]*x_i[1] + g[2]*x_i[2]), 4923a8779fbSJames Wright E_internal = E - E_kinetic - E_potential, 4933a8779fbSJames Wright P = E_internal * (gamma - 1.); // P = pressure 4943a8779fbSJames Wright 4953a8779fbSJames Wright // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 4963a8779fbSJames Wright CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 497*bb8a0c61SJames Wright computeFluxJacobian_NS(jacob_F_conv, rho, u, E, gamma, g, x_i); 4983a8779fbSJames Wright 4993a8779fbSJames Wright // dqdx collects drhodx, dUdx and dEdx in one vector 5003a8779fbSJames Wright CeedScalar dqdx[5][3]; 5013a8779fbSJames Wright for (int j=0; j<3; j++) { 5023a8779fbSJames Wright dqdx[0][j] = drhodx[j]; 5033a8779fbSJames Wright dqdx[4][j] = dEdx[j]; 5043a8779fbSJames Wright for (int k=0; k<3; k++) 5053a8779fbSJames Wright dqdx[k+1][j] = dUdx[k][j]; 5063a8779fbSJames Wright } 5073a8779fbSJames Wright 5083a8779fbSJames Wright // strong_conv = dF/dq * dq/dx (Strong convection) 5093a8779fbSJames Wright CeedScalar strong_conv[5] = {0}; 5103a8779fbSJames Wright for (int j=0; j<3; j++) 5113a8779fbSJames Wright for (int k=0; k<5; k++) 5123a8779fbSJames Wright for (int l=0; l<5; l++) 5133a8779fbSJames Wright strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 5143a8779fbSJames Wright 5153a8779fbSJames Wright // Body force 516*bb8a0c61SJames Wright const CeedScalar body_force[5] = {0, rho *g[0], rho *g[1], rho *g[2], 0}; 5173a8779fbSJames Wright 5183a8779fbSJames Wright // The Physics 5193a8779fbSJames Wright // Zero dv so all future terms can safely sum into it 5203a8779fbSJames Wright for (int j=0; j<5; j++) 5213a8779fbSJames Wright for (int k=0; k<3; k++) 5223a8779fbSJames Wright dv[k][j][i] = 0; 5233a8779fbSJames Wright 5243a8779fbSJames Wright // -- Density 5253a8779fbSJames Wright // ---- u rho 5263a8779fbSJames Wright for (int j=0; j<3; j++) 5273a8779fbSJames Wright dv[j][0][i] += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 5283a8779fbSJames Wright rho*u[2]*dXdx[j][2]); 5293a8779fbSJames Wright // -- Momentum 5303a8779fbSJames Wright // ---- rho (u x u) + P I3 5313a8779fbSJames Wright for (int j=0; j<3; j++) 5323a8779fbSJames Wright for (int k=0; k<3; k++) 5333a8779fbSJames Wright dv[k][j+1][i] += wdetJ*((rho*u[j]*u[0] + (j==0?P:0))*dXdx[k][0] + 5343a8779fbSJames Wright (rho*u[j]*u[1] + (j==1?P:0))*dXdx[k][1] + 5353a8779fbSJames Wright (rho*u[j]*u[2] + (j==2?P:0))*dXdx[k][2]); 5363a8779fbSJames Wright // ---- Fuvisc 5373a8779fbSJames Wright const CeedInt Fuviscidx[3][3] = {{0, 1, 2}, {1, 3, 4}, {2, 4, 5}}; // symmetric matrix indices 5383a8779fbSJames Wright for (int j=0; j<3; j++) 5393a8779fbSJames Wright for (int k=0; k<3; k++) 5403a8779fbSJames Wright dv[k][j+1][i] -= wdetJ*(Fu[Fuviscidx[j][0]]*dXdx[k][0] + 5413a8779fbSJames Wright Fu[Fuviscidx[j][1]]*dXdx[k][1] + 5423a8779fbSJames Wright Fu[Fuviscidx[j][2]]*dXdx[k][2]); 5433a8779fbSJames Wright // -- Total Energy Density 5443a8779fbSJames Wright // ---- (E + P) u 5453a8779fbSJames Wright for (int j=0; j<3; j++) 5463a8779fbSJames Wright dv[j][4][i] += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 5473a8779fbSJames Wright u[2]*dXdx[j][2]); 5483a8779fbSJames Wright // ---- Fevisc 5493a8779fbSJames Wright for (int j=0; j<3; j++) 5503a8779fbSJames Wright dv[j][4][i] -= wdetJ * (Fe[0]*dXdx[j][0] + Fe[1]*dXdx[j][1] + 5513a8779fbSJames Wright Fe[2]*dXdx[j][2]); 5523a8779fbSJames Wright // Body Force 5533a8779fbSJames Wright for (int j=0; j<5; j++) 5543a8779fbSJames Wright v[j][i] = wdetJ * body_force[j]; 5553a8779fbSJames Wright 556*bb8a0c61SJames Wright // Spatial Stabilization 557*bb8a0c61SJames Wright // -- Not used in favor of diagonal tau. Kept for future testing 558*bb8a0c61SJames Wright // const CeedScalar sound_speed = sqrt(gamma * P / rho); 559*bb8a0c61SJames Wright // CeedScalar Tau_x[3] = {0.}; 560*bb8a0c61SJames Wright // Tau_spatial(Tau_x, dXdx, u, sound_speed, context->c_tau, mu); 5613a8779fbSJames Wright 562*bb8a0c61SJames Wright // -- Stabilization method: none, SU, or SUPG 563*bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 564*bb8a0c61SJames Wright CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0}; 565*bb8a0c61SJames Wright CeedScalar Tau_d[3] = {0.}; 5663a8779fbSJames Wright switch (context->stabilization) { 5673a8779fbSJames Wright case STAB_NONE: // Galerkin 5683a8779fbSJames Wright break; 5693a8779fbSJames Wright case STAB_SU: // SU 570*bb8a0c61SJames Wright Tau_diagPrim(Tau_d, dXdx, u, cv, context, mu, dt, rho); 571*bb8a0c61SJames Wright tau_strong_conv[0] = Tau_d[0] * strong_conv[0]; 572*bb8a0c61SJames Wright tau_strong_conv[1] = Tau_d[1] * strong_conv[1]; 573*bb8a0c61SJames Wright tau_strong_conv[2] = Tau_d[1] * strong_conv[2]; 574*bb8a0c61SJames Wright tau_strong_conv[3] = Tau_d[1] * strong_conv[3]; 575*bb8a0c61SJames Wright tau_strong_conv[4] = Tau_d[2] * strong_conv[4]; 576*bb8a0c61SJames Wright PrimitiveToConservative_fwd(rho, u, E, Rd, cv, tau_strong_conv, 577*bb8a0c61SJames Wright tau_strong_conv_conservative); 5783a8779fbSJames Wright for (int j=0; j<3; j++) 5793a8779fbSJames Wright for (int k=0; k<5; k++) 5803a8779fbSJames Wright for (int l=0; l<5; l++) 581*bb8a0c61SJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l]; 5823a8779fbSJames Wright 5833a8779fbSJames Wright for (int j=0; j<5; j++) 5843a8779fbSJames Wright for (int k=0; k<3; k++) 5853a8779fbSJames Wright dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] + 5863a8779fbSJames Wright stab[j][1] * dXdx[k][1] + 5873a8779fbSJames Wright stab[j][2] * dXdx[k][2]); 5883a8779fbSJames Wright break; 5893a8779fbSJames Wright case STAB_SUPG: // SUPG is not implemented for explicit scheme 5903a8779fbSJames Wright break; 5913a8779fbSJames Wright } 5923a8779fbSJames Wright 5933a8779fbSJames Wright } // End Quadrature Point Loop 5943a8779fbSJames Wright 5953a8779fbSJames Wright // Return 5963a8779fbSJames Wright return 0; 5973a8779fbSJames Wright } 5983a8779fbSJames Wright 5993a8779fbSJames Wright // ***************************************************************************** 6003a8779fbSJames Wright // This QFunction implements the Navier-Stokes equations (mentioned above) with 6013a8779fbSJames Wright // implicit time stepping method 6023a8779fbSJames Wright // 6033a8779fbSJames Wright // SU = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) ) 6043a8779fbSJames Wright // SUPG = Galerkin + grad(v) . ( Ai^T * Tau * (q_dot + Aj q,j - body force) ) 6053a8779fbSJames Wright // (diffussive terms will be added later) 6063a8779fbSJames Wright // 6073a8779fbSJames Wright // ***************************************************************************** 6083a8779fbSJames Wright CEED_QFUNCTION(IFunction_Newtonian)(void *ctx, CeedInt Q, 6093a8779fbSJames Wright const CeedScalar *const *in, 6103a8779fbSJames Wright CeedScalar *const *out) { 6113a8779fbSJames Wright // *INDENT-OFF* 6123a8779fbSJames Wright // Inputs 6133a8779fbSJames Wright const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 6143a8779fbSJames Wright (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 6153a8779fbSJames Wright (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 6163a8779fbSJames Wright (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3], 6173a8779fbSJames Wright (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 6183a8779fbSJames Wright // Outputs 6193a8779fbSJames Wright CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 6203a8779fbSJames Wright (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 6213a8779fbSJames Wright // *INDENT-ON* 6223a8779fbSJames Wright // Context 6233a8779fbSJames Wright NewtonianIdealGasContext context = (NewtonianIdealGasContext)ctx; 6243a8779fbSJames Wright const CeedScalar lambda = context->lambda; 6253a8779fbSJames Wright const CeedScalar mu = context->mu; 6263a8779fbSJames Wright const CeedScalar k = context->k; 6273a8779fbSJames Wright const CeedScalar cv = context->cv; 6283a8779fbSJames Wright const CeedScalar cp = context->cp; 629*bb8a0c61SJames Wright const CeedScalar *g = context->g; 630*bb8a0c61SJames Wright const CeedScalar dt = context->dt; 6313a8779fbSJames Wright const CeedScalar gamma = cp / cv; 632*bb8a0c61SJames Wright const CeedScalar Rd = cp-cv; 6333a8779fbSJames Wright 6343a8779fbSJames Wright CeedPragmaSIMD 6353a8779fbSJames Wright // Quadrature Point Loop 6363a8779fbSJames Wright for (CeedInt i=0; i<Q; i++) { 6373a8779fbSJames Wright // Setup 6383a8779fbSJames Wright // -- Interp in 6393a8779fbSJames Wright const CeedScalar rho = q[0][i]; 6403a8779fbSJames Wright const CeedScalar u[3] = {q[1][i] / rho, 6413a8779fbSJames Wright q[2][i] / rho, 6423a8779fbSJames Wright q[3][i] / rho 6433a8779fbSJames Wright }; 6443a8779fbSJames Wright const CeedScalar E = q[4][i]; 6453a8779fbSJames Wright // -- Grad in 6463a8779fbSJames Wright const CeedScalar drho[3] = {dq[0][0][i], 6473a8779fbSJames Wright dq[1][0][i], 6483a8779fbSJames Wright dq[2][0][i] 6493a8779fbSJames Wright }; 6503a8779fbSJames Wright // *INDENT-OFF* 6513a8779fbSJames Wright const CeedScalar dU[3][3] = {{dq[0][1][i], 6523a8779fbSJames Wright dq[1][1][i], 6533a8779fbSJames Wright dq[2][1][i]}, 6543a8779fbSJames Wright {dq[0][2][i], 6553a8779fbSJames Wright dq[1][2][i], 6563a8779fbSJames Wright dq[2][2][i]}, 6573a8779fbSJames Wright {dq[0][3][i], 6583a8779fbSJames Wright dq[1][3][i], 6593a8779fbSJames Wright dq[2][3][i]} 6603a8779fbSJames Wright }; 6613a8779fbSJames Wright // *INDENT-ON* 6623a8779fbSJames Wright const CeedScalar dE[3] = {dq[0][4][i], 6633a8779fbSJames Wright dq[1][4][i], 6643a8779fbSJames Wright dq[2][4][i] 6653a8779fbSJames Wright }; 6663a8779fbSJames Wright // -- Interp-to-Interp q_data 6673a8779fbSJames Wright const CeedScalar wdetJ = q_data[0][i]; 6683a8779fbSJames Wright // -- Interp-to-Grad q_data 6693a8779fbSJames Wright // ---- Inverse of change of coordinate matrix: X_i,j 6703a8779fbSJames Wright // *INDENT-OFF* 6713a8779fbSJames Wright const CeedScalar dXdx[3][3] = {{q_data[1][i], 6723a8779fbSJames Wright q_data[2][i], 6733a8779fbSJames Wright q_data[3][i]}, 6743a8779fbSJames Wright {q_data[4][i], 6753a8779fbSJames Wright q_data[5][i], 6763a8779fbSJames Wright q_data[6][i]}, 6773a8779fbSJames Wright {q_data[7][i], 6783a8779fbSJames Wright q_data[8][i], 6793a8779fbSJames Wright q_data[9][i]} 6803a8779fbSJames Wright }; 681*bb8a0c61SJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 6823a8779fbSJames Wright // *INDENT-ON* 6833a8779fbSJames Wright // -- Grad-to-Grad q_data 6843a8779fbSJames Wright // dU/dx 6853a8779fbSJames Wright CeedScalar du[3][3] = {{0}}; 6863a8779fbSJames Wright CeedScalar drhodx[3] = {0}; 6873a8779fbSJames Wright CeedScalar dEdx[3] = {0}; 6883a8779fbSJames Wright CeedScalar dUdx[3][3] = {{0}}; 6893a8779fbSJames Wright CeedScalar dXdxdXdxT[3][3] = {{0}}; 6903a8779fbSJames Wright for (int j=0; j<3; j++) { 6913a8779fbSJames Wright for (int k=0; k<3; k++) { 6923a8779fbSJames Wright du[j][k] = (dU[j][k] - drho[k]*u[j]) / rho; 6933a8779fbSJames Wright drhodx[j] += drho[k] * dXdx[k][j]; 6943a8779fbSJames Wright dEdx[j] += dE[k] * dXdx[k][j]; 6953a8779fbSJames Wright for (int l=0; l<3; l++) { 6963a8779fbSJames Wright dUdx[j][k] += dU[j][l] * dXdx[l][k]; 6973a8779fbSJames Wright dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 6983a8779fbSJames Wright } 6993a8779fbSJames Wright } 7003a8779fbSJames Wright } 7013a8779fbSJames Wright CeedScalar dudx[3][3] = {{0}}; 7023a8779fbSJames Wright for (int j=0; j<3; j++) 7033a8779fbSJames Wright for (int k=0; k<3; k++) 7043a8779fbSJames Wright for (int l=0; l<3; l++) 7053a8779fbSJames Wright dudx[j][k] += du[j][l] * dXdx[l][k]; 7063a8779fbSJames Wright // -- grad_T 7073a8779fbSJames Wright const CeedScalar grad_T[3] = {(dEdx[0]/rho - E*drhodx[0]/(rho*rho) - /* *NOPAD* */ 708*bb8a0c61SJames Wright (u[0]*dudx[0][0] + u[1]*dudx[1][0] + u[2]*dudx[2][0]) + g[0])/cv, 7093a8779fbSJames Wright (dEdx[1]/rho - E*drhodx[1]/(rho*rho) - /* *NOPAD* */ 710*bb8a0c61SJames Wright (u[0]*dudx[0][1] + u[1]*dudx[1][1] + u[2]*dudx[2][1]) + g[1])/cv, 7113a8779fbSJames Wright (dEdx[2]/rho - E*drhodx[2]/(rho*rho) - /* *NOPAD* */ 712*bb8a0c61SJames Wright (u[0]*dudx[0][2] + u[1]*dudx[1][2] + u[2]*dudx[2][2]) + g[2])/cv 7133a8779fbSJames Wright }; 7143a8779fbSJames Wright // -- Fuvisc 7153a8779fbSJames Wright // ---- Symmetric 3x3 matrix 7163a8779fbSJames Wright const CeedScalar Fu[6] = {mu*(dudx[0][0] * (2 + lambda) + /* *NOPAD* */ 7173a8779fbSJames Wright lambda * (dudx[1][1] + dudx[2][2])), 7183a8779fbSJames Wright mu*(dudx[0][1] + dudx[1][0]), /* *NOPAD* */ 7193a8779fbSJames Wright mu*(dudx[0][2] + dudx[2][0]), /* *NOPAD* */ 7203a8779fbSJames Wright mu*(dudx[1][1] * (2 + lambda) + /* *NOPAD* */ 7213a8779fbSJames Wright lambda * (dudx[0][0] + dudx[2][2])), 7223a8779fbSJames Wright mu*(dudx[1][2] + dudx[2][1]), /* *NOPAD* */ 7233a8779fbSJames Wright mu*(dudx[2][2] * (2 + lambda) + /* *NOPAD* */ 7243a8779fbSJames Wright lambda * (dudx[0][0] + dudx[1][1])) 7253a8779fbSJames Wright }; 7263a8779fbSJames Wright // -- Fevisc 7273a8779fbSJames Wright const CeedScalar Fe[3] = {u[0]*Fu[0] + u[1]*Fu[1] + u[2]*Fu[2] + /* *NOPAD* */ 7283a8779fbSJames Wright k*grad_T[0], /* *NOPAD* */ 7293a8779fbSJames Wright u[0]*Fu[1] + u[1]*Fu[3] + u[2]*Fu[4] + /* *NOPAD* */ 7303a8779fbSJames Wright k*grad_T[1], /* *NOPAD* */ 7313a8779fbSJames Wright u[0]*Fu[2] + u[1]*Fu[4] + u[2]*Fu[5] + /* *NOPAD* */ 7323a8779fbSJames Wright k*grad_T[2] /* *NOPAD* */ 7333a8779fbSJames Wright }; 7343a8779fbSJames Wright // Pressure 7353a8779fbSJames Wright const CeedScalar 7363a8779fbSJames Wright E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 737*bb8a0c61SJames Wright E_potential = -rho*(g[0]*x_i[0] + g[1]*x_i[1] + g[2]*x_i[2]), 7383a8779fbSJames Wright E_internal = E - E_kinetic - E_potential, 7393a8779fbSJames Wright P = E_internal * (gamma - 1.); // P = pressure 7403a8779fbSJames Wright 7413a8779fbSJames Wright // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 7423a8779fbSJames Wright CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 743*bb8a0c61SJames Wright computeFluxJacobian_NS(jacob_F_conv, rho, u, E, gamma, g, x_i); 7443a8779fbSJames Wright 7453a8779fbSJames Wright // dqdx collects drhodx, dUdx and dEdx in one vector 7463a8779fbSJames Wright CeedScalar dqdx[5][3]; 7473a8779fbSJames Wright for (int j=0; j<3; j++) { 7483a8779fbSJames Wright dqdx[0][j] = drhodx[j]; 7493a8779fbSJames Wright dqdx[4][j] = dEdx[j]; 7503a8779fbSJames Wright for (int k=0; k<3; k++) 7513a8779fbSJames Wright dqdx[k+1][j] = dUdx[k][j]; 7523a8779fbSJames Wright } 7533a8779fbSJames Wright // strong_conv = dF/dq * dq/dx (Strong convection) 7543a8779fbSJames Wright CeedScalar strong_conv[5] = {0}; 7553a8779fbSJames Wright for (int j=0; j<3; j++) 7563a8779fbSJames Wright for (int k=0; k<5; k++) 7573a8779fbSJames Wright for (int l=0; l<5; l++) 7583a8779fbSJames Wright strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 7593a8779fbSJames Wright 7603a8779fbSJames Wright // Body force 761*bb8a0c61SJames Wright const CeedScalar body_force[5] = {0, rho *g[0], rho *g[1], rho *g[2], 0}; 7623a8779fbSJames Wright 7633a8779fbSJames Wright // Strong residual 7643a8779fbSJames Wright CeedScalar strong_res[5]; 7653a8779fbSJames Wright for (int j=0; j<5; j++) 7663a8779fbSJames Wright strong_res[j] = q_dot[j][i] + strong_conv[j] - body_force[j]; 7673a8779fbSJames Wright 7683a8779fbSJames Wright // The Physics 7693a8779fbSJames Wright //-----mass matrix 7703a8779fbSJames Wright for (int j=0; j<5; j++) 7713a8779fbSJames Wright v[j][i] = wdetJ*q_dot[j][i]; 7723a8779fbSJames Wright 7733a8779fbSJames Wright // Zero dv so all future terms can safely sum into it 7743a8779fbSJames Wright for (int j=0; j<5; j++) 7753a8779fbSJames Wright for (int k=0; k<3; k++) 7763a8779fbSJames Wright dv[k][j][i] = 0; 7773a8779fbSJames Wright 7783a8779fbSJames Wright // -- Density 7793a8779fbSJames Wright // ---- u rho 7803a8779fbSJames Wright for (int j=0; j<3; j++) 7813a8779fbSJames Wright dv[j][0][i] -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 7823a8779fbSJames Wright rho*u[2]*dXdx[j][2]); 7833a8779fbSJames Wright // -- Momentum 7843a8779fbSJames Wright // ---- rho (u x u) + P I3 7853a8779fbSJames Wright for (int j=0; j<3; j++) 7863a8779fbSJames Wright for (int k=0; k<3; k++) 7873a8779fbSJames Wright dv[k][j+1][i] -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0))*dXdx[k][0] + 7883a8779fbSJames Wright (rho*u[j]*u[1] + (j==1?P:0))*dXdx[k][1] + 7893a8779fbSJames Wright (rho*u[j]*u[2] + (j==2?P:0))*dXdx[k][2]); 7903a8779fbSJames Wright // ---- Fuvisc 7913a8779fbSJames Wright const CeedInt Fuviscidx[3][3] = {{0, 1, 2}, {1, 3, 4}, {2, 4, 5}}; // symmetric matrix indices 7923a8779fbSJames Wright for (int j=0; j<3; j++) 7933a8779fbSJames Wright for (int k=0; k<3; k++) 7943a8779fbSJames Wright dv[k][j+1][i] += wdetJ*(Fu[Fuviscidx[j][0]]*dXdx[k][0] + 7953a8779fbSJames Wright Fu[Fuviscidx[j][1]]*dXdx[k][1] + 7963a8779fbSJames Wright Fu[Fuviscidx[j][2]]*dXdx[k][2]); 7973a8779fbSJames Wright // -- Total Energy Density 7983a8779fbSJames Wright // ---- (E + P) u 7993a8779fbSJames Wright for (int j=0; j<3; j++) 8003a8779fbSJames Wright dv[j][4][i] -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 8013a8779fbSJames Wright u[2]*dXdx[j][2]); 8023a8779fbSJames Wright // ---- Fevisc 8033a8779fbSJames Wright for (int j=0; j<3; j++) 8043a8779fbSJames Wright dv[j][4][i] += wdetJ * (Fe[0]*dXdx[j][0] + Fe[1]*dXdx[j][1] + 8053a8779fbSJames Wright Fe[2]*dXdx[j][2]); 8063a8779fbSJames Wright // Body Force 8073a8779fbSJames Wright for (int j=0; j<5; j++) 8083a8779fbSJames Wright v[j][i] -= wdetJ*body_force[j]; 8093a8779fbSJames Wright 810*bb8a0c61SJames Wright // Spatial Stabilization 811*bb8a0c61SJames Wright // -- Not used in favor of diagonal tau. Kept for future testing 812*bb8a0c61SJames Wright // const CeedScalar sound_speed = sqrt(gamma * P / rho); 813*bb8a0c61SJames Wright // CeedScalar Tau_x[3] = {0.}; 814*bb8a0c61SJames Wright // Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau, mu); 8153a8779fbSJames Wright 8163a8779fbSJames Wright // -- Stabilization method: none, SU, or SUPG 817*bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 818*bb8a0c61SJames Wright CeedScalar tau_strong_res[5] = {0.}, tau_strong_res_conservative[5] = {0}; 819*bb8a0c61SJames Wright CeedScalar tau_strong_conv[5] = {0.}, tau_strong_conv_conservative[5] = {0}; 820*bb8a0c61SJames Wright CeedScalar Tau_d[3] = {0.}; 8213a8779fbSJames Wright switch (context->stabilization) { 8223a8779fbSJames Wright case STAB_NONE: // Galerkin 8233a8779fbSJames Wright break; 8243a8779fbSJames Wright case STAB_SU: // SU 825*bb8a0c61SJames Wright Tau_diagPrim(Tau_d, dXdx, u, cv, context, mu, dt, rho); 826*bb8a0c61SJames Wright tau_strong_conv[0] = Tau_d[0] * strong_conv[0]; 827*bb8a0c61SJames Wright tau_strong_conv[1] = Tau_d[1] * strong_conv[1]; 828*bb8a0c61SJames Wright tau_strong_conv[2] = Tau_d[1] * strong_conv[2]; 829*bb8a0c61SJames Wright tau_strong_conv[3] = Tau_d[1] * strong_conv[3]; 830*bb8a0c61SJames Wright tau_strong_conv[4] = Tau_d[2] * strong_conv[4]; 831*bb8a0c61SJames Wright PrimitiveToConservative_fwd(rho, u, E, Rd, cv, tau_strong_conv, 832*bb8a0c61SJames Wright tau_strong_conv_conservative); 8333a8779fbSJames Wright for (int j=0; j<3; j++) 8343a8779fbSJames Wright for (int k=0; k<5; k++) 8353a8779fbSJames Wright for (int l=0; l<5; l++) 836*bb8a0c61SJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_conv_conservative[l]; 8373a8779fbSJames Wright 8383a8779fbSJames Wright for (int j=0; j<5; j++) 8393a8779fbSJames Wright for (int k=0; k<3; k++) 8403a8779fbSJames Wright dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 8413a8779fbSJames Wright stab[j][1] * dXdx[k][1] + 8423a8779fbSJames Wright stab[j][2] * dXdx[k][2]); 8433a8779fbSJames Wright break; 8443a8779fbSJames Wright case STAB_SUPG: // SUPG 845*bb8a0c61SJames Wright Tau_diagPrim(Tau_d, dXdx, u, cv, context, mu, dt, rho); 846*bb8a0c61SJames Wright tau_strong_res[0] = Tau_d[0] * strong_res[0]; 847*bb8a0c61SJames Wright tau_strong_res[1] = Tau_d[1] * strong_res[1]; 848*bb8a0c61SJames Wright tau_strong_res[2] = Tau_d[1] * strong_res[2]; 849*bb8a0c61SJames Wright tau_strong_res[3] = Tau_d[1] * strong_res[3]; 850*bb8a0c61SJames Wright tau_strong_res[4] = Tau_d[2] * strong_res[4]; 851*bb8a0c61SJames Wright // Alternate route (useful later with primitive variable code) 852*bb8a0c61SJames Wright // this function was verified against PHASTA for as IC that was as close as possible 853*bb8a0c61SJames Wright // computeFluxJacobian_NSp(jacob_F_conv_p, rho, u, E, Rd, cv); 854*bb8a0c61SJames Wright // it has also been verified to compute a correct through the following 855*bb8a0c61SJames Wright // stab[k][j] += jacob_F_conv_p[j][k][l] * tau_strong_res[l] // flux Jacobian wrt primitive 856*bb8a0c61SJames Wright // applied in the triple loop below 857*bb8a0c61SJames Wright // However, it is more flops than using the existing Jacobian wrt q after q_{,Y} viz 858*bb8a0c61SJames Wright PrimitiveToConservative_fwd(rho, u, E, Rd, cv, tau_strong_res, 859*bb8a0c61SJames Wright tau_strong_res_conservative); 8603a8779fbSJames Wright for (int j=0; j<3; j++) 8613a8779fbSJames Wright for (int k=0; k<5; k++) 8623a8779fbSJames Wright for (int l=0; l<5; l++) 863*bb8a0c61SJames Wright stab[k][j] += jacob_F_conv[j][k][l] * tau_strong_res_conservative[l]; 8643a8779fbSJames Wright 8653a8779fbSJames Wright for (int j=0; j<5; j++) 8663a8779fbSJames Wright for (int k=0; k<3; k++) 8673a8779fbSJames Wright dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 8683a8779fbSJames Wright stab[j][1] * dXdx[k][1] + 8693a8779fbSJames Wright stab[j][2] * dXdx[k][2]); 8703a8779fbSJames Wright break; 8713a8779fbSJames Wright } 8723a8779fbSJames Wright 8733a8779fbSJames Wright } // End Quadrature Point Loop 8743a8779fbSJames Wright 8753a8779fbSJames Wright // Return 8763a8779fbSJames Wright return 0; 8773a8779fbSJames Wright } 8783a8779fbSJames Wright // ***************************************************************************** 8793a8779fbSJames Wright #endif // newtonian_h 880