177841947SLeila Ghaffari // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 277841947SLeila Ghaffari // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 377841947SLeila Ghaffari // reserved. See files LICENSE and NOTICE for details. 477841947SLeila Ghaffari // 577841947SLeila Ghaffari // This file is part of CEED, a collection of benchmarks, miniapps, software 677841947SLeila Ghaffari // libraries and APIs for efficient high-order finite element and spectral 777841947SLeila Ghaffari // element discretizations for exascale applications. For more information and 877841947SLeila Ghaffari // source code availability see http://github.com/ceed. 977841947SLeila Ghaffari // 1077841947SLeila Ghaffari // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 1177841947SLeila Ghaffari // a collaborative effort of two U.S. Department of Energy organizations (Office 1277841947SLeila Ghaffari // of Science and the National Nuclear Security Administration) responsible for 1377841947SLeila Ghaffari // the planning and preparation of a capable exascale ecosystem, including 1477841947SLeila Ghaffari // software, applications, hardware, advanced system engineering and early 1577841947SLeila Ghaffari // testbed platforms, in support of the nation's exascale computing imperative. 1677841947SLeila Ghaffari 1777841947SLeila Ghaffari /// @file 1877841947SLeila Ghaffari /// Density current initial condition and operator for Navier-Stokes example using PETSc 1977841947SLeila Ghaffari 2077841947SLeila Ghaffari // Model from: 2177841947SLeila Ghaffari // Semi-Implicit Formulations of the Navier-Stokes Equations: Application to 2277841947SLeila Ghaffari // Nonhydrostatic Atmospheric Modeling, Giraldo, Restelli, and Lauter (2010). 2377841947SLeila Ghaffari 2477841947SLeila Ghaffari #ifndef densitycurrent_h 2577841947SLeila Ghaffari #define densitycurrent_h 2677841947SLeila Ghaffari 2777841947SLeila Ghaffari #include <math.h> 2877841947SLeila Ghaffari 2977841947SLeila Ghaffari #ifndef M_PI 3077841947SLeila Ghaffari #define M_PI 3.14159265358979323846 3177841947SLeila Ghaffari #endif 3277841947SLeila Ghaffari 3377841947SLeila Ghaffari #ifndef setup_context_struct 3477841947SLeila Ghaffari #define setup_context_struct 3577841947SLeila Ghaffari typedef struct SetupContext_ *SetupContext; 3677841947SLeila Ghaffari struct SetupContext_ { 3777841947SLeila Ghaffari CeedScalar theta0; 3877841947SLeila Ghaffari CeedScalar thetaC; 3977841947SLeila Ghaffari CeedScalar P0; 4077841947SLeila Ghaffari CeedScalar N; 4177841947SLeila Ghaffari CeedScalar cv; 4277841947SLeila Ghaffari CeedScalar cp; 4377841947SLeila Ghaffari CeedScalar g; 4477841947SLeila Ghaffari CeedScalar rc; 4577841947SLeila Ghaffari CeedScalar lx; 4677841947SLeila Ghaffari CeedScalar ly; 4777841947SLeila Ghaffari CeedScalar lz; 4877841947SLeila Ghaffari CeedScalar center[3]; 4977841947SLeila Ghaffari CeedScalar dc_axis[3]; 5077841947SLeila Ghaffari CeedScalar wind[3]; 5177841947SLeila Ghaffari CeedScalar time; 5277841947SLeila Ghaffari int wind_type; // See WindType: 0=ROTATION, 1=TRANSLATION 5377841947SLeila Ghaffari int bubble_type; // See BubbleType: 0=SPHERE, 1=CYLINDER 5477841947SLeila Ghaffari int bubble_continuity_type; // See BubbleContinuityType: 0=SMOOTH, 1=BACK_SHARP 2=THICK 5577841947SLeila Ghaffari }; 5677841947SLeila Ghaffari #endif 5777841947SLeila Ghaffari 5877841947SLeila Ghaffari #ifndef dc_context_struct 5977841947SLeila Ghaffari #define dc_context_struct 6077841947SLeila Ghaffari typedef struct DCContext_ *DCContext; 6177841947SLeila Ghaffari struct DCContext_ { 6277841947SLeila Ghaffari CeedScalar lambda; 6377841947SLeila Ghaffari CeedScalar mu; 6477841947SLeila Ghaffari CeedScalar k; 6577841947SLeila Ghaffari CeedScalar cv; 6677841947SLeila Ghaffari CeedScalar cp; 6777841947SLeila Ghaffari CeedScalar g; 6877841947SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 6977841947SLeila Ghaffari }; 7077841947SLeila Ghaffari #endif 7177841947SLeila Ghaffari 7277841947SLeila Ghaffari // ***************************************************************************** 7377841947SLeila Ghaffari // This function sets the initial conditions and the boundary conditions 7477841947SLeila Ghaffari // 7577841947SLeila Ghaffari // These initial conditions are given in terms of potential temperature and 7677841947SLeila Ghaffari // Exner pressure and then converted to density and total energy. 7777841947SLeila Ghaffari // Initial momentum density is zero. 7877841947SLeila Ghaffari // 7977841947SLeila Ghaffari // Initial Conditions: 8077841947SLeila Ghaffari // Potential Temperature: 8177841947SLeila Ghaffari // theta = thetabar + delta_theta 8277841947SLeila Ghaffari // thetabar = theta0 exp( N**2 z / g ) 8377841947SLeila Ghaffari // delta_theta = r <= rc : thetaC(1 + cos(pi r/rc)) / 2 8477841947SLeila Ghaffari // r > rc : 0 8577841947SLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 + (z - zc)**2 ) 8677841947SLeila Ghaffari // with (xc,yc,zc) center of domain, rc characteristic radius of thermal bubble 8777841947SLeila Ghaffari // Exner Pressure: 8877841947SLeila Ghaffari // Pi = Pibar + deltaPi 8977841947SLeila Ghaffari // Pibar = 1. + g**2 (exp( - N**2 z / g ) - 1) / (cp theta0 N**2) 9077841947SLeila Ghaffari // deltaPi = 0 (hydrostatic balance) 9177841947SLeila Ghaffari // Velocity/Momentum Density: 9277841947SLeila Ghaffari // Ui = ui = 0 9377841947SLeila Ghaffari // 9477841947SLeila Ghaffari // Conversion to Conserved Variables: 9577841947SLeila Ghaffari // rho = P0 Pi**(cv/Rd) / (Rd theta) 9677841947SLeila Ghaffari // E = rho (cv T + (u u)/2 + g z) 9777841947SLeila Ghaffari // 9877841947SLeila Ghaffari // Boundary Conditions: 9977841947SLeila Ghaffari // Mass Density: 10077841947SLeila Ghaffari // 0.0 flux 10177841947SLeila Ghaffari // Momentum Density: 10277841947SLeila Ghaffari // 0.0 10377841947SLeila Ghaffari // Energy Density: 10477841947SLeila Ghaffari // 0.0 flux 10577841947SLeila Ghaffari // 10677841947SLeila Ghaffari // Constants: 10777841947SLeila Ghaffari // theta0 , Potential temperature constant 10877841947SLeila Ghaffari // thetaC , Potential temperature perturbation 10977841947SLeila Ghaffari // P0 , Pressure at the surface 11077841947SLeila Ghaffari // N , Brunt-Vaisala frequency 11177841947SLeila Ghaffari // cv , Specific heat, constant volume 11277841947SLeila Ghaffari // cp , Specific heat, constant pressure 11377841947SLeila Ghaffari // Rd = cp - cv, Specific heat difference 11477841947SLeila Ghaffari // g , Gravity 11577841947SLeila Ghaffari // rc , Characteristic radius of thermal bubble 11677841947SLeila Ghaffari // center , Location of bubble center 11777841947SLeila Ghaffari // dc_axis , Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric 11877841947SLeila Ghaffari // ***************************************************************************** 11977841947SLeila Ghaffari 12077841947SLeila Ghaffari // ***************************************************************************** 12177841947SLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 12277841947SLeila Ghaffari // (currently not implemented) and IC formulation for density current 12377841947SLeila Ghaffari // ***************************************************************************** 12477841947SLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_DC(CeedInt dim, CeedScalar time, 12577841947SLeila Ghaffari const CeedScalar X[], CeedInt Nf, CeedScalar q[], 12677841947SLeila Ghaffari void *ctx) { 12777841947SLeila Ghaffari // Context 12877841947SLeila Ghaffari const SetupContext context = (SetupContext)ctx; 12977841947SLeila Ghaffari const CeedScalar theta0 = context->theta0; 13077841947SLeila Ghaffari const CeedScalar thetaC = context->thetaC; 13177841947SLeila Ghaffari const CeedScalar P0 = context->P0; 13277841947SLeila Ghaffari const CeedScalar N = context->N; 13377841947SLeila Ghaffari const CeedScalar cv = context->cv; 13477841947SLeila Ghaffari const CeedScalar cp = context->cp; 13577841947SLeila Ghaffari const CeedScalar g = context->g; 13677841947SLeila Ghaffari const CeedScalar rc = context->rc; 13777841947SLeila Ghaffari const CeedScalar *center = context->center; 13877841947SLeila Ghaffari const CeedScalar *dc_axis = context->dc_axis; 139*e6225c47SLeila Ghaffari const CeedScalar Rd = cp - cv; 14077841947SLeila Ghaffari 14177841947SLeila Ghaffari // Setup 14277841947SLeila Ghaffari // -- Coordinates 14377841947SLeila Ghaffari const CeedScalar x = X[0]; 14477841947SLeila Ghaffari const CeedScalar y = X[1]; 14577841947SLeila Ghaffari const CeedScalar z = X[2]; 14677841947SLeila Ghaffari 14777841947SLeila Ghaffari // -- Potential temperature, density current 14877841947SLeila Ghaffari CeedScalar rr[3] = {x - center[0], y - center[1], z - center[2]}; 14977841947SLeila Ghaffari // (I - q q^T) r: distance from dc_axis (or from center if dc_axis is the zero vector) 15077841947SLeila Ghaffari for (CeedInt i=0; i<3; i++) 15177841947SLeila Ghaffari rr[i] -= dc_axis[i] * 15277841947SLeila Ghaffari (dc_axis[0]*rr[0] + dc_axis[1]*rr[1] + dc_axis[2]*rr[2]); 15377841947SLeila Ghaffari const CeedScalar r = sqrt(rr[0]*rr[0] + rr[1]*rr[1] + rr[2]*rr[2]); 15477841947SLeila Ghaffari const CeedScalar delta_theta = r <= rc ? thetaC*(1. + cos(M_PI*r/rc))/2. : 0.; 15577841947SLeila Ghaffari const CeedScalar theta = theta0*exp(N*N*z/g) + delta_theta; 15677841947SLeila Ghaffari 15777841947SLeila Ghaffari // -- Exner pressure, hydrostatic balance 15877841947SLeila Ghaffari const CeedScalar Pi = 1. + g*g*(exp(-N*N*z/g) - 1.) / (cp*theta0*N*N); 15977841947SLeila Ghaffari // -- Density 16077841947SLeila Ghaffari 16177841947SLeila Ghaffari const CeedScalar rho = P0 * pow(Pi, cv/Rd) / (Rd*theta); 16277841947SLeila Ghaffari 16377841947SLeila Ghaffari // Initial Conditions 16477841947SLeila Ghaffari q[0] = rho; 16577841947SLeila Ghaffari q[1] = 0.0; 16677841947SLeila Ghaffari q[2] = 0.0; 16777841947SLeila Ghaffari q[3] = 0.0; 16877841947SLeila Ghaffari q[4] = rho * (cv*theta*Pi + g*z); 16977841947SLeila Ghaffari 17077841947SLeila Ghaffari return 0; 17177841947SLeila Ghaffari } 17277841947SLeila Ghaffari 17377841947SLeila Ghaffari // ***************************************************************************** 174*e6225c47SLeila Ghaffari // Helper function for computing flux Jacobian 175*e6225c47SLeila Ghaffari // ***************************************************************************** 176*e6225c47SLeila Ghaffari CEED_QFUNCTION_HELPER void computeFluxJacobian_NS(CeedScalar dF[3][5][5], 177*e6225c47SLeila Ghaffari const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 178*e6225c47SLeila Ghaffari const CeedScalar gamma, const CeedScalar g, CeedScalar z) { 179*e6225c47SLeila Ghaffari CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 180*e6225c47SLeila Ghaffari for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 181*e6225c47SLeila Ghaffari for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix 182*e6225c47SLeila Ghaffari dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2. - g*z)) : 0.) - u[i]*u[j]; 183*e6225c47SLeila Ghaffari for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix 184*e6225c47SLeila Ghaffari dF[i][0][k+1] = ((i==k) ? 1. : 0.); 185*e6225c47SLeila Ghaffari dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) + 186*e6225c47SLeila Ghaffari ((i==k) ? u[j] : 0.) - 187*e6225c47SLeila Ghaffari ((i==j) ? u[k] : 0.) * (gamma-1.); 188*e6225c47SLeila Ghaffari dF[i][4][k+1] = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) - 189*e6225c47SLeila Ghaffari (gamma-1.)*u[i]*u[k]; 190*e6225c47SLeila Ghaffari } 191*e6225c47SLeila Ghaffari dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.); 192*e6225c47SLeila Ghaffari } 193*e6225c47SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho); 194*e6225c47SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 195*e6225c47SLeila Ghaffari } 196*e6225c47SLeila Ghaffari } 197*e6225c47SLeila Ghaffari 198*e6225c47SLeila Ghaffari // ***************************************************************************** 19977841947SLeila Ghaffari // This QFunction sets the initial conditions for density current 20077841947SLeila Ghaffari // ***************************************************************************** 20177841947SLeila Ghaffari CEED_QFUNCTION(ICsDC)(void *ctx, CeedInt Q, 20277841947SLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 20377841947SLeila Ghaffari // Inputs 20477841947SLeila Ghaffari const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 20577841947SLeila Ghaffari 20677841947SLeila Ghaffari // Outputs 20777841947SLeila Ghaffari CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 20877841947SLeila Ghaffari 20977841947SLeila Ghaffari CeedPragmaSIMD 21077841947SLeila Ghaffari // Quadrature Point Loop 21177841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 21277841947SLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 213*e6225c47SLeila Ghaffari CeedScalar q[5] = {0.}; 21477841947SLeila Ghaffari 21577841947SLeila Ghaffari Exact_DC(3, 0., x, 5, q, ctx); 21677841947SLeila Ghaffari 21777841947SLeila Ghaffari for (CeedInt j=0; j<5; j++) 21877841947SLeila Ghaffari q0[j][i] = q[j]; 21977841947SLeila Ghaffari } // End of Quadrature Point Loop 22077841947SLeila Ghaffari 22177841947SLeila Ghaffari // Return 22277841947SLeila Ghaffari return 0; 22377841947SLeila Ghaffari } 22477841947SLeila Ghaffari 22577841947SLeila Ghaffari // ***************************************************************************** 22677841947SLeila Ghaffari // This QFunction implements the following formulation of Navier-Stokes with 22777841947SLeila Ghaffari // explicit time stepping method 22877841947SLeila Ghaffari // 22977841947SLeila Ghaffari // This is 3D compressible Navier-Stokes in conservation form with state 23077841947SLeila Ghaffari // variables of density, momentum density, and total energy density. 23177841947SLeila Ghaffari // 23277841947SLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 23377841947SLeila Ghaffari // rho - Mass Density 23477841947SLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 23577841947SLeila Ghaffari // E - Total Energy Density, E = rho (cv T + (u u)/2 + g z) 23677841947SLeila Ghaffari // 23777841947SLeila Ghaffari // Navier-Stokes Equations: 23877841947SLeila Ghaffari // drho/dt + div( U ) = 0 23977841947SLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) + rho g khat = div( Fu ) 24077841947SLeila Ghaffari // dE/dt + div( (E + P) u ) = div( Fe ) 24177841947SLeila Ghaffari // 24277841947SLeila Ghaffari // Viscous Stress: 24377841947SLeila Ghaffari // Fu = mu (grad( u ) + grad( u )^T + lambda div ( u ) I3) 24477841947SLeila Ghaffari // 24577841947SLeila Ghaffari // Thermal Stress: 24677841947SLeila Ghaffari // Fe = u Fu + k grad( T ) 24777841947SLeila Ghaffari // 24877841947SLeila Ghaffari // Equation of State: 24977841947SLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2 - rho g z) 25077841947SLeila Ghaffari // 25177841947SLeila Ghaffari // Stabilization: 25277841947SLeila Ghaffari // Tau = diag(TauC, TauM, TauM, TauM, TauE) 25377841947SLeila Ghaffari // f1 = rho sqrt(ui uj gij) 25477841947SLeila Ghaffari // gij = dXi/dX * dXi/dX 25577841947SLeila Ghaffari // TauC = Cc f1 / (8 gii) 25677841947SLeila Ghaffari // TauM = min( 1 , 1 / f1 ) 25777841947SLeila Ghaffari // TauE = TauM / (Ce cv) 25877841947SLeila Ghaffari // 25977841947SLeila Ghaffari // SU = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) ) 26077841947SLeila Ghaffari // 26177841947SLeila Ghaffari // Constants: 26277841947SLeila Ghaffari // lambda = - 2 / 3, From Stokes hypothesis 26377841947SLeila Ghaffari // mu , Dynamic viscosity 26477841947SLeila Ghaffari // k , Thermal conductivity 26577841947SLeila Ghaffari // cv , Specific heat, constant volume 26677841947SLeila Ghaffari // cp , Specific heat, constant pressure 26777841947SLeila Ghaffari // g , Gravity 26877841947SLeila Ghaffari // gamma = cp / cv, Specific heat ratio 26977841947SLeila Ghaffari // 27077841947SLeila Ghaffari // We require the product of the inverse of the Jacobian (dXdx_j,k) and 27177841947SLeila Ghaffari // its transpose (dXdx_k,j) to properly compute integrals of the form: 27277841947SLeila Ghaffari // int( gradv gradu ) 27377841947SLeila Ghaffari // 27477841947SLeila Ghaffari // ***************************************************************************** 27577841947SLeila Ghaffari CEED_QFUNCTION(DC)(void *ctx, CeedInt Q, 27677841947SLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 27777841947SLeila Ghaffari // *INDENT-OFF* 27877841947SLeila Ghaffari // Inputs 27977841947SLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 28077841947SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 28177841947SLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 28277841947SLeila Ghaffari (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 28377841947SLeila Ghaffari // Outputs 28477841947SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 28577841947SLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 28677841947SLeila Ghaffari // *INDENT-ON* 28777841947SLeila Ghaffari 28877841947SLeila Ghaffari // Context 28977841947SLeila Ghaffari DCContext context = (DCContext)ctx; 29077841947SLeila Ghaffari const CeedScalar lambda = context->lambda; 29177841947SLeila Ghaffari const CeedScalar mu = context->mu; 29277841947SLeila Ghaffari const CeedScalar k = context->k; 29377841947SLeila Ghaffari const CeedScalar cv = context->cv; 29477841947SLeila Ghaffari const CeedScalar cp = context->cp; 29577841947SLeila Ghaffari const CeedScalar g = context->g; 29677841947SLeila Ghaffari const CeedScalar gamma = cp / cv; 29777841947SLeila Ghaffari 29877841947SLeila Ghaffari CeedPragmaSIMD 29977841947SLeila Ghaffari // Quadrature Point Loop 30077841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 30177841947SLeila Ghaffari // *INDENT-OFF* 30277841947SLeila Ghaffari // Setup 30377841947SLeila Ghaffari // -- Interp in 30477841947SLeila Ghaffari const CeedScalar rho = q[0][i]; 30577841947SLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 30677841947SLeila Ghaffari q[2][i] / rho, 30777841947SLeila Ghaffari q[3][i] / rho 30877841947SLeila Ghaffari }; 30977841947SLeila Ghaffari const CeedScalar E = q[4][i]; 31077841947SLeila Ghaffari // -- Grad in 31177841947SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 31277841947SLeila Ghaffari dq[1][0][i], 31377841947SLeila Ghaffari dq[2][0][i] 31477841947SLeila Ghaffari }; 31577841947SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 31677841947SLeila Ghaffari dq[1][1][i], 31777841947SLeila Ghaffari dq[2][1][i]}, 31877841947SLeila Ghaffari {dq[0][2][i], 31977841947SLeila Ghaffari dq[1][2][i], 32077841947SLeila Ghaffari dq[2][2][i]}, 32177841947SLeila Ghaffari {dq[0][3][i], 32277841947SLeila Ghaffari dq[1][3][i], 32377841947SLeila Ghaffari dq[2][3][i]} 32477841947SLeila Ghaffari }; 32577841947SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 32677841947SLeila Ghaffari dq[1][4][i], 32777841947SLeila Ghaffari dq[2][4][i] 32877841947SLeila Ghaffari }; 32977841947SLeila Ghaffari // -- Interp-to-Interp q_data 33077841947SLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 33177841947SLeila Ghaffari // -- Interp-to-Grad q_data 33277841947SLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 33377841947SLeila Ghaffari // *INDENT-OFF* 33477841947SLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 33577841947SLeila Ghaffari q_data[2][i], 33677841947SLeila Ghaffari q_data[3][i]}, 33777841947SLeila Ghaffari {q_data[4][i], 33877841947SLeila Ghaffari q_data[5][i], 33977841947SLeila Ghaffari q_data[6][i]}, 34077841947SLeila Ghaffari {q_data[7][i], 34177841947SLeila Ghaffari q_data[8][i], 34277841947SLeila Ghaffari q_data[9][i]} 34377841947SLeila Ghaffari }; 34477841947SLeila Ghaffari // *INDENT-ON* 34577841947SLeila Ghaffari // -- Grad-to-Grad q_data 34677841947SLeila Ghaffari // dU/dx 34777841947SLeila Ghaffari CeedScalar du[3][3] = {{0}}; 34877841947SLeila Ghaffari CeedScalar drhodx[3] = {0}; 34977841947SLeila Ghaffari CeedScalar dEdx[3] = {0}; 35077841947SLeila Ghaffari CeedScalar dUdx[3][3] = {{0}}; 35177841947SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0}}; 35277841947SLeila Ghaffari for (int j=0; j<3; j++) { 35377841947SLeila Ghaffari for (int k=0; k<3; k++) { 35477841947SLeila Ghaffari du[j][k] = (dU[j][k] - drho[k]*u[j]) / rho; 35577841947SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 35677841947SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 35777841947SLeila Ghaffari for (int l=0; l<3; l++) { 35877841947SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 35977841947SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 36077841947SLeila Ghaffari } 36177841947SLeila Ghaffari } 36277841947SLeila Ghaffari } 36377841947SLeila Ghaffari CeedScalar dudx[3][3] = {{0}}; 36477841947SLeila Ghaffari for (int j=0; j<3; j++) 36577841947SLeila Ghaffari for (int k=0; k<3; k++) 36677841947SLeila Ghaffari for (int l=0; l<3; l++) 36777841947SLeila Ghaffari dudx[j][k] += du[j][l] * dXdx[l][k]; 36877841947SLeila Ghaffari // -- grad_T 36977841947SLeila Ghaffari const CeedScalar grad_T[3] = {(dEdx[0]/rho - E*drhodx[0]/(rho*rho) - /* *NOPAD* */ 37077841947SLeila Ghaffari (u[0]*dudx[0][0] + u[1]*dudx[1][0] + u[2]*dudx[2][0]))/cv, 37177841947SLeila Ghaffari (dEdx[1]/rho - E*drhodx[1]/(rho*rho) - /* *NOPAD* */ 37277841947SLeila Ghaffari (u[0]*dudx[0][1] + u[1]*dudx[1][1] + u[2]*dudx[2][1]))/cv, 37377841947SLeila Ghaffari (dEdx[2]/rho - E*drhodx[2]/(rho*rho) - /* *NOPAD* */ 37477841947SLeila Ghaffari (u[0]*dudx[0][2] + u[1]*dudx[1][2] + u[2]*dudx[2][2]) - g)/cv 37577841947SLeila Ghaffari }; 37677841947SLeila Ghaffari 37777841947SLeila Ghaffari // -- Fuvisc 37877841947SLeila Ghaffari // ---- Symmetric 3x3 matrix 37977841947SLeila Ghaffari const CeedScalar Fu[6] = {mu*(dudx[0][0] * (2 + lambda) + /* *NOPAD* */ 38077841947SLeila Ghaffari lambda * (dudx[1][1] + dudx[2][2])), 38177841947SLeila Ghaffari mu*(dudx[0][1] + dudx[1][0]), /* *NOPAD* */ 38277841947SLeila Ghaffari mu*(dudx[0][2] + dudx[2][0]), /* *NOPAD* */ 38377841947SLeila Ghaffari mu*(dudx[1][1] * (2 + lambda) + /* *NOPAD* */ 38477841947SLeila Ghaffari lambda * (dudx[0][0] + dudx[2][2])), 38577841947SLeila Ghaffari mu*(dudx[1][2] + dudx[2][1]), /* *NOPAD* */ 38677841947SLeila Ghaffari mu*(dudx[2][2] * (2 + lambda) + /* *NOPAD* */ 38777841947SLeila Ghaffari lambda * (dudx[0][0] + dudx[1][1])) 38877841947SLeila Ghaffari }; 38977841947SLeila Ghaffari // -- Fevisc 39077841947SLeila Ghaffari const CeedScalar Fe[3] = {u[0]*Fu[0] + u[1]*Fu[1] + u[2]*Fu[2] + /* *NOPAD* */ 39177841947SLeila Ghaffari k*grad_T[0], /* *NOPAD* */ 39277841947SLeila Ghaffari u[0]*Fu[1] + u[1]*Fu[3] + u[2]*Fu[4] + /* *NOPAD* */ 39377841947SLeila Ghaffari k*grad_T[1], /* *NOPAD* */ 39477841947SLeila Ghaffari u[0]*Fu[2] + u[1]*Fu[4] + u[2]*Fu[5] + /* *NOPAD* */ 39577841947SLeila Ghaffari k*grad_T[2] /* *NOPAD* */ 39677841947SLeila Ghaffari }; 397*e6225c47SLeila Ghaffari // Pressure 398*e6225c47SLeila Ghaffari const CeedScalar 399*e6225c47SLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 400*e6225c47SLeila Ghaffari E_potential = rho*g*x[2][i], 401*e6225c47SLeila Ghaffari E_internal = E - E_kinetic - E_potential, 402*e6225c47SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 403*e6225c47SLeila Ghaffari 40477841947SLeila Ghaffari // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 405*e6225c47SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 406*e6225c47SLeila Ghaffari computeFluxJacobian_NS(jacob_F_conv, rho, u, E, gamma, g, x[2][i]); 40777841947SLeila Ghaffari 40877841947SLeila Ghaffari // jacob_F_conv_T = jacob_F_conv^T 40977841947SLeila Ghaffari CeedScalar jacob_F_conv_T[3][5][5]; 41077841947SLeila Ghaffari for (int j=0; j<3; j++) 41177841947SLeila Ghaffari for (int k=0; k<5; k++) 41277841947SLeila Ghaffari for (int l=0; l<5; l++) 41377841947SLeila Ghaffari jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k]; 41477841947SLeila Ghaffari 41577841947SLeila Ghaffari // dqdx collects drhodx, dUdx and dEdx in one vector 41677841947SLeila Ghaffari CeedScalar dqdx[5][3]; 41777841947SLeila Ghaffari for (int j=0; j<3; j++) { 41877841947SLeila Ghaffari dqdx[0][j] = drhodx[j]; 41977841947SLeila Ghaffari dqdx[4][j] = dEdx[j]; 42077841947SLeila Ghaffari for (int k=0; k<3; k++) 42177841947SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 42277841947SLeila Ghaffari } 42377841947SLeila Ghaffari 42477841947SLeila Ghaffari // strong_conv = dF/dq * dq/dx (Strong convection) 42577841947SLeila Ghaffari CeedScalar strong_conv[5] = {0}; 42677841947SLeila Ghaffari for (int j=0; j<3; j++) 42777841947SLeila Ghaffari for (int k=0; k<5; k++) 42877841947SLeila Ghaffari for (int l=0; l<5; l++) 42977841947SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 43077841947SLeila Ghaffari 43177841947SLeila Ghaffari // Body force 43277841947SLeila Ghaffari const CeedScalar body_force[5] = {0, 0, 0, -rho*g, 0}; 43377841947SLeila Ghaffari 43477841947SLeila Ghaffari // The Physics 43577841947SLeila Ghaffari // Zero dv so all future terms can safely sum into it 43677841947SLeila Ghaffari for (int j=0; j<5; j++) 43777841947SLeila Ghaffari for (int k=0; k<3; k++) 43877841947SLeila Ghaffari dv[k][j][i] = 0; 43977841947SLeila Ghaffari 44077841947SLeila Ghaffari // -- Density 44177841947SLeila Ghaffari // ---- u rho 44277841947SLeila Ghaffari for (int j=0; j<3; j++) 44377841947SLeila Ghaffari dv[j][0][i] += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 44477841947SLeila Ghaffari rho*u[2]*dXdx[j][2]); 44577841947SLeila Ghaffari // -- Momentum 44677841947SLeila Ghaffari // ---- rho (u x u) + P I3 44777841947SLeila Ghaffari for (int j=0; j<3; j++) 44877841947SLeila Ghaffari for (int k=0; k<3; k++) 44977841947SLeila Ghaffari dv[k][j+1][i] += wdetJ*((rho*u[j]*u[0] + (j==0?P:0))*dXdx[k][0] + 45077841947SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0))*dXdx[k][1] + 45177841947SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0))*dXdx[k][2]); 45277841947SLeila Ghaffari // ---- Fuvisc 45377841947SLeila Ghaffari const CeedInt Fuviscidx[3][3] = {{0, 1, 2}, {1, 3, 4}, {2, 4, 5}}; // symmetric matrix indices 45477841947SLeila Ghaffari for (int j=0; j<3; j++) 45577841947SLeila Ghaffari for (int k=0; k<3; k++) 45677841947SLeila Ghaffari dv[k][j+1][i] -= wdetJ*(Fu[Fuviscidx[j][0]]*dXdx[k][0] + 45777841947SLeila Ghaffari Fu[Fuviscidx[j][1]]*dXdx[k][1] + 45877841947SLeila Ghaffari Fu[Fuviscidx[j][2]]*dXdx[k][2]); 45977841947SLeila Ghaffari // -- Total Energy Density 46077841947SLeila Ghaffari // ---- (E + P) u 46177841947SLeila Ghaffari for (int j=0; j<3; j++) 46277841947SLeila Ghaffari dv[j][4][i] += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 46377841947SLeila Ghaffari u[2]*dXdx[j][2]); 46477841947SLeila Ghaffari // ---- Fevisc 46577841947SLeila Ghaffari for (int j=0; j<3; j++) 46677841947SLeila Ghaffari dv[j][4][i] -= wdetJ * (Fe[0]*dXdx[j][0] + Fe[1]*dXdx[j][1] + 46777841947SLeila Ghaffari Fe[2]*dXdx[j][2]); 46877841947SLeila Ghaffari // Body Force 46977841947SLeila Ghaffari for (int j=0; j<5; j++) 47077841947SLeila Ghaffari v[j][i] = wdetJ * body_force[j]; 47177841947SLeila Ghaffari 47277841947SLeila Ghaffari //Stabilization 47377841947SLeila Ghaffari CeedScalar uX[3]; 474*e6225c47SLeila Ghaffari for (int j=0; j<3; j++) 475*e6225c47SLeila Ghaffari uX[j] = dXdx[j][0]*u[0] + dXdx[j][1]*u[1] + dXdx[j][2]*u[2]; 47677841947SLeila Ghaffari const CeedScalar uiujgij = uX[0]*uX[0] + uX[1]*uX[1] + uX[2]*uX[2]; 47777841947SLeila Ghaffari const CeedScalar Cc = 1.; 47877841947SLeila Ghaffari const CeedScalar Ce = 1.; 47977841947SLeila Ghaffari const CeedScalar f1 = rho * sqrt(uiujgij); 48077841947SLeila Ghaffari const CeedScalar TauC = (Cc * f1) / 48177841947SLeila Ghaffari (8 * (dXdxdXdxT[0][0] + dXdxdXdxT[1][1] + dXdxdXdxT[2][2])); 48277841947SLeila Ghaffari const CeedScalar TauM = 1. / (f1>1. ? f1 : 1.); 48377841947SLeila Ghaffari const CeedScalar TauE = TauM / (Ce * cv); 48477841947SLeila Ghaffari const CeedScalar Tau[5] = {TauC, TauM, TauM, TauM, TauE}; 48577841947SLeila Ghaffari CeedScalar stab[5][3]; 48677841947SLeila Ghaffari switch (context->stabilization) { 48777841947SLeila Ghaffari case 0: // Galerkin 48877841947SLeila Ghaffari break; 48977841947SLeila Ghaffari case 1: // SU 49077841947SLeila Ghaffari for (int j=0; j<3; j++) 49177841947SLeila Ghaffari for (int k=0; k<5; k++) 49277841947SLeila Ghaffari for (int l=0; l<5; l++) 49377841947SLeila Ghaffari stab[k][j] = jacob_F_conv_T[j][k][l] * Tau[l] * strong_conv[l]; 49477841947SLeila Ghaffari 49577841947SLeila Ghaffari for (int j=0; j<5; j++) 49677841947SLeila Ghaffari for (int k=0; k<3; k++) 49777841947SLeila Ghaffari dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] + 49877841947SLeila Ghaffari stab[j][1] * dXdx[k][1] + 49977841947SLeila Ghaffari stab[j][2] * dXdx[k][2]); 50077841947SLeila Ghaffari break; 50177841947SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 50277841947SLeila Ghaffari break; 50377841947SLeila Ghaffari } 50477841947SLeila Ghaffari 50577841947SLeila Ghaffari } // End Quadrature Point Loop 50677841947SLeila Ghaffari 50777841947SLeila Ghaffari // Return 50877841947SLeila Ghaffari return 0; 50977841947SLeila Ghaffari } 51077841947SLeila Ghaffari 51177841947SLeila Ghaffari // ***************************************************************************** 51277841947SLeila Ghaffari // This QFunction implements the Navier-Stokes equations (mentioned above) with 51377841947SLeila Ghaffari // implicit time stepping method 51477841947SLeila Ghaffari // 51577841947SLeila Ghaffari // SU = Galerkin + grad(v) . ( Ai^T * Tau * (Aj q,j) ) 51677841947SLeila Ghaffari // SUPG = Galerkin + grad(v) . ( Ai^T * Tau * (q_dot + Aj q,j - body force) ) 51777841947SLeila Ghaffari // (diffussive terms will be added later) 51877841947SLeila Ghaffari // 51977841947SLeila Ghaffari // ***************************************************************************** 52077841947SLeila Ghaffari CEED_QFUNCTION(IFunction_DC)(void *ctx, CeedInt Q, 52177841947SLeila Ghaffari const CeedScalar *const *in, 52277841947SLeila Ghaffari CeedScalar *const *out) { 52377841947SLeila Ghaffari // *INDENT-OFF* 52477841947SLeila Ghaffari // Inputs 52577841947SLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 52677841947SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 52777841947SLeila Ghaffari (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 52877841947SLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3], 52977841947SLeila Ghaffari (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 53077841947SLeila Ghaffari // Outputs 53177841947SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 53277841947SLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 53377841947SLeila Ghaffari // *INDENT-ON* 53477841947SLeila Ghaffari // Context 53577841947SLeila Ghaffari DCContext context = (DCContext)ctx; 53677841947SLeila Ghaffari const CeedScalar lambda = context->lambda; 53777841947SLeila Ghaffari const CeedScalar mu = context->mu; 53877841947SLeila Ghaffari const CeedScalar k = context->k; 53977841947SLeila Ghaffari const CeedScalar cv = context->cv; 54077841947SLeila Ghaffari const CeedScalar cp = context->cp; 54177841947SLeila Ghaffari const CeedScalar g = context->g; 54277841947SLeila Ghaffari const CeedScalar gamma = cp / cv; 54377841947SLeila Ghaffari 54477841947SLeila Ghaffari CeedPragmaSIMD 54577841947SLeila Ghaffari // Quadrature Point Loop 54677841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 54777841947SLeila Ghaffari // Setup 54877841947SLeila Ghaffari // -- Interp in 54977841947SLeila Ghaffari const CeedScalar rho = q[0][i]; 55077841947SLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 55177841947SLeila Ghaffari q[2][i] / rho, 55277841947SLeila Ghaffari q[3][i] / rho 55377841947SLeila Ghaffari }; 55477841947SLeila Ghaffari const CeedScalar E = q[4][i]; 55577841947SLeila Ghaffari // -- Grad in 55677841947SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 55777841947SLeila Ghaffari dq[1][0][i], 55877841947SLeila Ghaffari dq[2][0][i] 55977841947SLeila Ghaffari }; 56077841947SLeila Ghaffari // *INDENT-OFF* 56177841947SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 56277841947SLeila Ghaffari dq[1][1][i], 56377841947SLeila Ghaffari dq[2][1][i]}, 56477841947SLeila Ghaffari {dq[0][2][i], 56577841947SLeila Ghaffari dq[1][2][i], 56677841947SLeila Ghaffari dq[2][2][i]}, 56777841947SLeila Ghaffari {dq[0][3][i], 56877841947SLeila Ghaffari dq[1][3][i], 56977841947SLeila Ghaffari dq[2][3][i]} 57077841947SLeila Ghaffari }; 57177841947SLeila Ghaffari // *INDENT-ON* 57277841947SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 57377841947SLeila Ghaffari dq[1][4][i], 57477841947SLeila Ghaffari dq[2][4][i] 57577841947SLeila Ghaffari }; 57677841947SLeila Ghaffari // -- Interp-to-Interp q_data 57777841947SLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 57877841947SLeila Ghaffari // -- Interp-to-Grad q_data 57977841947SLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 58077841947SLeila Ghaffari // *INDENT-OFF* 58177841947SLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 58277841947SLeila Ghaffari q_data[2][i], 58377841947SLeila Ghaffari q_data[3][i]}, 58477841947SLeila Ghaffari {q_data[4][i], 58577841947SLeila Ghaffari q_data[5][i], 58677841947SLeila Ghaffari q_data[6][i]}, 58777841947SLeila Ghaffari {q_data[7][i], 58877841947SLeila Ghaffari q_data[8][i], 58977841947SLeila Ghaffari q_data[9][i]} 59077841947SLeila Ghaffari }; 59177841947SLeila Ghaffari // *INDENT-ON* 59277841947SLeila Ghaffari // -- Grad-to-Grad q_data 59377841947SLeila Ghaffari // dU/dx 59477841947SLeila Ghaffari CeedScalar du[3][3] = {{0}}; 59577841947SLeila Ghaffari CeedScalar drhodx[3] = {0}; 59677841947SLeila Ghaffari CeedScalar dEdx[3] = {0}; 59777841947SLeila Ghaffari CeedScalar dUdx[3][3] = {{0}}; 59877841947SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0}}; 59977841947SLeila Ghaffari for (int j=0; j<3; j++) { 60077841947SLeila Ghaffari for (int k=0; k<3; k++) { 60177841947SLeila Ghaffari du[j][k] = (dU[j][k] - drho[k]*u[j]) / rho; 60277841947SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 60377841947SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 60477841947SLeila Ghaffari for (int l=0; l<3; l++) { 60577841947SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 60677841947SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 60777841947SLeila Ghaffari } 60877841947SLeila Ghaffari } 60977841947SLeila Ghaffari } 61077841947SLeila Ghaffari CeedScalar dudx[3][3] = {{0}}; 61177841947SLeila Ghaffari for (int j=0; j<3; j++) 61277841947SLeila Ghaffari for (int k=0; k<3; k++) 61377841947SLeila Ghaffari for (int l=0; l<3; l++) 61477841947SLeila Ghaffari dudx[j][k] += du[j][l] * dXdx[l][k]; 61577841947SLeila Ghaffari // -- grad_T 61677841947SLeila Ghaffari const CeedScalar grad_T[3] = {(dEdx[0]/rho - E*drhodx[0]/(rho*rho) - /* *NOPAD* */ 61777841947SLeila Ghaffari (u[0]*dudx[0][0] + u[1]*dudx[1][0] + u[2]*dudx[2][0]))/cv, 61877841947SLeila Ghaffari (dEdx[1]/rho - E*drhodx[1]/(rho*rho) - /* *NOPAD* */ 61977841947SLeila Ghaffari (u[0]*dudx[0][1] + u[1]*dudx[1][1] + u[2]*dudx[2][1]))/cv, 62077841947SLeila Ghaffari (dEdx[2]/rho - E*drhodx[2]/(rho*rho) - /* *NOPAD* */ 62177841947SLeila Ghaffari (u[0]*dudx[0][2] + u[1]*dudx[1][2] + u[2]*dudx[2][2]) - g)/cv 62277841947SLeila Ghaffari }; 62377841947SLeila Ghaffari // -- Fuvisc 62477841947SLeila Ghaffari // ---- Symmetric 3x3 matrix 62577841947SLeila Ghaffari const CeedScalar Fu[6] = {mu*(dudx[0][0] * (2 + lambda) + /* *NOPAD* */ 62677841947SLeila Ghaffari lambda * (dudx[1][1] + dudx[2][2])), 62777841947SLeila Ghaffari mu*(dudx[0][1] + dudx[1][0]), /* *NOPAD* */ 62877841947SLeila Ghaffari mu*(dudx[0][2] + dudx[2][0]), /* *NOPAD* */ 62977841947SLeila Ghaffari mu*(dudx[1][1] * (2 + lambda) + /* *NOPAD* */ 63077841947SLeila Ghaffari lambda * (dudx[0][0] + dudx[2][2])), 63177841947SLeila Ghaffari mu*(dudx[1][2] + dudx[2][1]), /* *NOPAD* */ 63277841947SLeila Ghaffari mu*(dudx[2][2] * (2 + lambda) + /* *NOPAD* */ 63377841947SLeila Ghaffari lambda * (dudx[0][0] + dudx[1][1])) 63477841947SLeila Ghaffari }; 63577841947SLeila Ghaffari // -- Fevisc 63677841947SLeila Ghaffari const CeedScalar Fe[3] = {u[0]*Fu[0] + u[1]*Fu[1] + u[2]*Fu[2] + /* *NOPAD* */ 63777841947SLeila Ghaffari k*grad_T[0], /* *NOPAD* */ 63877841947SLeila Ghaffari u[0]*Fu[1] + u[1]*Fu[3] + u[2]*Fu[4] + /* *NOPAD* */ 63977841947SLeila Ghaffari k*grad_T[1], /* *NOPAD* */ 64077841947SLeila Ghaffari u[0]*Fu[2] + u[1]*Fu[4] + u[2]*Fu[5] + /* *NOPAD* */ 64177841947SLeila Ghaffari k*grad_T[2] /* *NOPAD* */ 64277841947SLeila Ghaffari }; 643*e6225c47SLeila Ghaffari // Pressure 644*e6225c47SLeila Ghaffari const CeedScalar 645*e6225c47SLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 646*e6225c47SLeila Ghaffari E_potential = rho*g*x[2][i], 647*e6225c47SLeila Ghaffari E_internal = E - E_kinetic - E_potential, 648*e6225c47SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 64977841947SLeila Ghaffari 65077841947SLeila Ghaffari // jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 651*e6225c47SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 652*e6225c47SLeila Ghaffari computeFluxJacobian_NS(jacob_F_conv, rho, u, E, gamma, g, x[2][i]); 65377841947SLeila Ghaffari 65477841947SLeila Ghaffari // jacob_F_conv_T = jacob_F_conv^T 65577841947SLeila Ghaffari CeedScalar jacob_F_conv_T[3][5][5]; 65677841947SLeila Ghaffari for (int j=0; j<3; j++) 65777841947SLeila Ghaffari for (int k=0; k<5; k++) 65877841947SLeila Ghaffari for (int l=0; l<5; l++) 65977841947SLeila Ghaffari jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k]; 66077841947SLeila Ghaffari // dqdx collects drhodx, dUdx and dEdx in one vector 66177841947SLeila Ghaffari CeedScalar dqdx[5][3]; 66277841947SLeila Ghaffari for (int j=0; j<3; j++) { 66377841947SLeila Ghaffari dqdx[0][j] = drhodx[j]; 66477841947SLeila Ghaffari dqdx[4][j] = dEdx[j]; 66577841947SLeila Ghaffari for (int k=0; k<3; k++) 66677841947SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 66777841947SLeila Ghaffari } 66877841947SLeila Ghaffari // strong_conv = dF/dq * dq/dx (Strong convection) 66977841947SLeila Ghaffari CeedScalar strong_conv[5] = {0}; 67077841947SLeila Ghaffari for (int j=0; j<3; j++) 67177841947SLeila Ghaffari for (int k=0; k<5; k++) 67277841947SLeila Ghaffari for (int l=0; l<5; l++) 67377841947SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 67477841947SLeila Ghaffari 67577841947SLeila Ghaffari // Body force 67677841947SLeila Ghaffari const CeedScalar body_force[5] = {0, 0, 0, -rho*g, 0}; 67777841947SLeila Ghaffari 67877841947SLeila Ghaffari // Strong residual 67977841947SLeila Ghaffari CeedScalar strong_res[5]; 68077841947SLeila Ghaffari for (int j=0; j<5; j++) 68177841947SLeila Ghaffari strong_res[j] = q_dot[j][i] + strong_conv[j] - body_force[j]; 68277841947SLeila Ghaffari 68377841947SLeila Ghaffari // The Physics 68477841947SLeila Ghaffari //-----mass matrix 68577841947SLeila Ghaffari for (int j=0; j<5; j++) 68677841947SLeila Ghaffari v[j][i] = wdetJ*q_dot[j][i]; 68777841947SLeila Ghaffari 68877841947SLeila Ghaffari // Zero dv so all future terms can safely sum into it 68977841947SLeila Ghaffari for (int j=0; j<5; j++) 69077841947SLeila Ghaffari for (int k=0; k<3; k++) 69177841947SLeila Ghaffari dv[k][j][i] = 0; 69277841947SLeila Ghaffari 69377841947SLeila Ghaffari // -- Density 69477841947SLeila Ghaffari // ---- u rho 69577841947SLeila Ghaffari for (int j=0; j<3; j++) 69677841947SLeila Ghaffari dv[j][0][i] -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 69777841947SLeila Ghaffari rho*u[2]*dXdx[j][2]); 69877841947SLeila Ghaffari // -- Momentum 69977841947SLeila Ghaffari // ---- rho (u x u) + P I3 70077841947SLeila Ghaffari for (int j=0; j<3; j++) 70177841947SLeila Ghaffari for (int k=0; k<3; k++) 70277841947SLeila Ghaffari dv[k][j+1][i] -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0))*dXdx[k][0] + 70377841947SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0))*dXdx[k][1] + 70477841947SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0))*dXdx[k][2]); 70577841947SLeila Ghaffari // ---- Fuvisc 70677841947SLeila Ghaffari const CeedInt Fuviscidx[3][3] = {{0, 1, 2}, {1, 3, 4}, {2, 4, 5}}; // symmetric matrix indices 70777841947SLeila Ghaffari for (int j=0; j<3; j++) 70877841947SLeila Ghaffari for (int k=0; k<3; k++) 70977841947SLeila Ghaffari dv[k][j+1][i] += wdetJ*(Fu[Fuviscidx[j][0]]*dXdx[k][0] + 71077841947SLeila Ghaffari Fu[Fuviscidx[j][1]]*dXdx[k][1] + 71177841947SLeila Ghaffari Fu[Fuviscidx[j][2]]*dXdx[k][2]); 71277841947SLeila Ghaffari // -- Total Energy Density 71377841947SLeila Ghaffari // ---- (E + P) u 71477841947SLeila Ghaffari for (int j=0; j<3; j++) 71577841947SLeila Ghaffari dv[j][4][i] -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 71677841947SLeila Ghaffari u[2]*dXdx[j][2]); 71777841947SLeila Ghaffari // ---- Fevisc 71877841947SLeila Ghaffari for (int j=0; j<3; j++) 71977841947SLeila Ghaffari dv[j][4][i] += wdetJ * (Fe[0]*dXdx[j][0] + Fe[1]*dXdx[j][1] + 72077841947SLeila Ghaffari Fe[2]*dXdx[j][2]); 72177841947SLeila Ghaffari // Body Force 72277841947SLeila Ghaffari for (int j=0; j<5; j++) 72377841947SLeila Ghaffari v[j][i] -= wdetJ*body_force[j]; 72477841947SLeila Ghaffari 72577841947SLeila Ghaffari //Stabilization 72677841947SLeila Ghaffari CeedScalar uX[3]; 727*e6225c47SLeila Ghaffari for (int j=0; j<3; j++) 728*e6225c47SLeila Ghaffari uX[j] = dXdx[j][0]*u[0] + dXdx[j][1]*u[1] + dXdx[j][2]*u[2]; 72977841947SLeila Ghaffari const CeedScalar uiujgij = uX[0]*uX[0] + uX[1]*uX[1] + uX[2]*uX[2]; 73077841947SLeila Ghaffari const CeedScalar Cc = 1.; 73177841947SLeila Ghaffari const CeedScalar Ce = 1.; 73277841947SLeila Ghaffari const CeedScalar f1 = rho * sqrt(uiujgij); 73377841947SLeila Ghaffari const CeedScalar TauC = (Cc * f1) / 73477841947SLeila Ghaffari (8 * (dXdxdXdxT[0][0] + dXdxdXdxT[1][1] + dXdxdXdxT[2][2])); 73577841947SLeila Ghaffari const CeedScalar TauM = 1. / (f1>1. ? f1 : 1.); 73677841947SLeila Ghaffari const CeedScalar TauE = TauM / (Ce * cv); 73777841947SLeila Ghaffari const CeedScalar Tau[5] = {TauC, TauM, TauM, TauM, TauE}; 738*e6225c47SLeila Ghaffari 73977841947SLeila Ghaffari CeedScalar stab[5][3]; 74077841947SLeila Ghaffari switch (context->stabilization) { 74177841947SLeila Ghaffari case 0: // Galerkin 74277841947SLeila Ghaffari break; 74377841947SLeila Ghaffari case 1: // SU 74477841947SLeila Ghaffari for (int j=0; j<3; j++) 74577841947SLeila Ghaffari for (int k=0; k<5; k++) 74677841947SLeila Ghaffari for (int l=0; l<5; l++) 74777841947SLeila Ghaffari stab[k][j] = jacob_F_conv_T[j][k][l] * Tau[l] * strong_conv[l]; 74877841947SLeila Ghaffari 74977841947SLeila Ghaffari for (int j=0; j<5; j++) 75077841947SLeila Ghaffari for (int k=0; k<3; k++) 75177841947SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 75277841947SLeila Ghaffari stab[j][1] * dXdx[k][1] + 75377841947SLeila Ghaffari stab[j][2] * dXdx[k][2]); 75477841947SLeila Ghaffari break; 75577841947SLeila Ghaffari case 2: // SUPG 75677841947SLeila Ghaffari for (int j=0; j<3; j++) 75777841947SLeila Ghaffari for (int k=0; k<5; k++) 75877841947SLeila Ghaffari for (int l=0; l<5; l++) 75977841947SLeila Ghaffari stab[k][j] = jacob_F_conv_T[j][k][l] * Tau[l] * strong_res[l]; 76077841947SLeila Ghaffari 76177841947SLeila Ghaffari for (int j=0; j<5; j++) 76277841947SLeila Ghaffari for (int k=0; k<3; k++) 76377841947SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 76477841947SLeila Ghaffari stab[j][1] * dXdx[k][1] + 76577841947SLeila Ghaffari stab[j][2] * dXdx[k][2]); 76677841947SLeila Ghaffari break; 76777841947SLeila Ghaffari } 76877841947SLeila Ghaffari 76977841947SLeila Ghaffari } // End Quadrature Point Loop 77077841947SLeila Ghaffari 77177841947SLeila Ghaffari // Return 77277841947SLeila Ghaffari return 0; 77377841947SLeila Ghaffari } 77477841947SLeila Ghaffari // ***************************************************************************** 77577841947SLeila Ghaffari 77677841947SLeila Ghaffari #endif // densitycurrent_h 777