13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 377841947SLeila Ghaffari // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 577841947SLeila Ghaffari // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 777841947SLeila Ghaffari 877841947SLeila Ghaffari /// @file 977841947SLeila Ghaffari /// Euler traveling vortex initial condition and operator for Navier-Stokes 1077841947SLeila Ghaffari /// example using PETSc 1177841947SLeila Ghaffari 1277841947SLeila Ghaffari // Model from: 1377841947SLeila Ghaffari // On the Order of Accuracy and Numerical Performance of Two Classes of 1477841947SLeila Ghaffari // Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011). 1577841947SLeila Ghaffari 1677841947SLeila Ghaffari #ifndef eulervortex_h 1777841947SLeila Ghaffari #define eulervortex_h 1877841947SLeila Ghaffari 1977841947SLeila Ghaffari #include <math.h> 2088b783a1SJames Wright #include <ceed.h> 2177841947SLeila Ghaffari 2277841947SLeila Ghaffari #ifndef M_PI 2377841947SLeila Ghaffari #define M_PI 3.14159265358979323846 2477841947SLeila Ghaffari #endif 2577841947SLeila Ghaffari 2677841947SLeila Ghaffari typedef struct EulerContext_ *EulerContext; 2777841947SLeila Ghaffari struct EulerContext_ { 2877841947SLeila Ghaffari CeedScalar center[3]; 2977841947SLeila Ghaffari CeedScalar curr_time; 3077841947SLeila Ghaffari CeedScalar vortex_strength; 31932417b3SJed Brown CeedScalar c_tau; 3277841947SLeila Ghaffari CeedScalar mean_velocity[3]; 3377841947SLeila Ghaffari bool implicit; 34e6225c47SLeila Ghaffari int euler_test; 35e6225c47SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 3677841947SLeila Ghaffari }; 3777841947SLeila Ghaffari 3877841947SLeila Ghaffari // ***************************************************************************** 3977841947SLeila Ghaffari // This function sets the initial conditions 4077841947SLeila Ghaffari // 4177841947SLeila Ghaffari // Temperature: 4277841947SLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 4377841947SLeila Ghaffari // Density: 4477841947SLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1)) 4577841947SLeila Ghaffari // Pressure: 4677841947SLeila Ghaffari // P = rho * T 4777841947SLeila Ghaffari // Velocity: 4877841947SLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 4977841947SLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 5077841947SLeila Ghaffari // Velocity/Momentum Density: 5177841947SLeila Ghaffari // Ui = rho ui 5277841947SLeila Ghaffari // Total Energy: 5377841947SLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2 5477841947SLeila Ghaffari // 5577841947SLeila Ghaffari // Constants: 5677841947SLeila Ghaffari // cv , Specific heat, constant volume 5777841947SLeila Ghaffari // cp , Specific heat, constant pressure 5877841947SLeila Ghaffari // vortex_strength , Strength of vortex 5977841947SLeila Ghaffari // center , Location of bubble center 6077841947SLeila Ghaffari // gamma = cp / cv, Specific heat ratio 6177841947SLeila Ghaffari // 6277841947SLeila Ghaffari // ***************************************************************************** 6377841947SLeila Ghaffari 6477841947SLeila Ghaffari // ***************************************************************************** 6577841947SLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 6677841947SLeila Ghaffari // (currently not implemented) and IC formulation for Euler traveling vortex 6777841947SLeila Ghaffari // ***************************************************************************** 6877841947SLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, 6977841947SLeila Ghaffari const CeedScalar X[], CeedInt Nf, CeedScalar q[], 7077841947SLeila Ghaffari void *ctx) { 7177841947SLeila Ghaffari // Context 7277841947SLeila Ghaffari const EulerContext context = (EulerContext)ctx; 7377841947SLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength; 7477841947SLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain 7577841947SLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity; 7677841947SLeila Ghaffari 7777841947SLeila Ghaffari // Setup 7877841947SLeila Ghaffari const CeedScalar gamma = 1.4; 7977841947SLeila Ghaffari const CeedScalar cv = 2.5; 8077841947SLeila Ghaffari const CeedScalar R = 1.; 8177841947SLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates 8277841947SLeila Ghaffari // Vortex center 8377841947SLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time; 8477841947SLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time; 8577841947SLeila Ghaffari 8677841947SLeila Ghaffari const CeedScalar x0 = x - xc; 8777841947SLeila Ghaffari const CeedScalar y0 = y - yc; 8877841947SLeila Ghaffari const CeedScalar r = sqrt( x0*x0 + y0*y0 ); 8977841947SLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r*r)/2.) / (2. * M_PI); 90e6225c47SLeila Ghaffari const CeedScalar delta_T = - (gamma - 1.) * vortex_strength * vortex_strength * 91e6225c47SLeila Ghaffari exp(1 - r*r) / (8. * gamma * M_PI * M_PI); 9277841947SLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 9377841947SLeila Ghaffari const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / 9477841947SLeila Ghaffari (8.*gamma*M_PI*M_PI); 9577841947SLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.}; 9677841947SLeila Ghaffari 9777841947SLeila Ghaffari // Initial Conditions 9877841947SLeila Ghaffari switch (context->euler_test) { 9977841947SLeila Ghaffari case 0: // Traveling vortex 10077841947SLeila Ghaffari T = 1 + delta_T; 10177841947SLeila Ghaffari // P = rho * T 10277841947SLeila Ghaffari // P = S * rho^gamma 10377841947SLeila Ghaffari // Solve for rho, then substitute for P 104e6225c47SLeila Ghaffari rho = pow(T/S_vortex, 1 / (gamma - 1.)); 10577841947SLeila Ghaffari P = rho * T; 10677841947SLeila Ghaffari u[0] = mean_velocity[0] - C*y0; 10777841947SLeila Ghaffari u[1] = mean_velocity[1] + C*x0; 10877841947SLeila Ghaffari 10977841947SLeila Ghaffari // Assign exact solution 11077841947SLeila Ghaffari q[0] = rho; 11177841947SLeila Ghaffari q[1] = rho * u[0]; 11277841947SLeila Ghaffari q[2] = rho * u[1]; 11377841947SLeila Ghaffari q[3] = rho * u[2]; 11477841947SLeila Ghaffari q[4] = P / (gamma - 1.) + rho * (u[0]*u[0] + u[1]*u[1]) / 2.; 11577841947SLeila Ghaffari break; 11677841947SLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant 11777841947SLeila Ghaffari rho = 1.; 11877841947SLeila Ghaffari E = 2.; 11977841947SLeila Ghaffari 12077841947SLeila Ghaffari // Assign exact solution 12177841947SLeila Ghaffari q[0] = rho; 12277841947SLeila Ghaffari q[1] = rho * u[0]; 12377841947SLeila Ghaffari q[2] = rho * u[1]; 12477841947SLeila Ghaffari q[3] = rho * u[2]; 12577841947SLeila Ghaffari q[4] = E; 12677841947SLeila Ghaffari break; 12777841947SLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant 12877841947SLeila Ghaffari rho = 1.; 12977841947SLeila Ghaffari E = 2.; 13077841947SLeila Ghaffari u[0] = mean_velocity[0]; 13177841947SLeila Ghaffari u[1] = mean_velocity[1]; 13277841947SLeila Ghaffari 13377841947SLeila Ghaffari // Assign exact solution 13477841947SLeila Ghaffari q[0] = rho; 13577841947SLeila Ghaffari q[1] = rho * u[0]; 13677841947SLeila Ghaffari q[2] = rho * u[1]; 13777841947SLeila Ghaffari q[3] = rho * u[2]; 13877841947SLeila Ghaffari q[4] = E; 13977841947SLeila Ghaffari break; 14077841947SLeila Ghaffari case 3: // Velocity zero, pressure constant 14177841947SLeila Ghaffari // (so density and internal energy will be non-constant), 14277841947SLeila Ghaffari // but the velocity should stay zero and the bubble won't diffuse 14377841947SLeila Ghaffari // (for Euler, where there is no thermal conductivity) 14477841947SLeila Ghaffari P = 1.; 14577841947SLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 14677841947SLeila Ghaffari rho = P / (R*T); 14777841947SLeila Ghaffari 14877841947SLeila Ghaffari // Assign exact solution 14977841947SLeila Ghaffari q[0] = rho; 15077841947SLeila Ghaffari q[1] = rho * u[0]; 15177841947SLeila Ghaffari q[2] = rho * u[1]; 15277841947SLeila Ghaffari q[3] = rho * u[2]; 15377841947SLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 15477841947SLeila Ghaffari break; 15577841947SLeila Ghaffari case 4: // Constant nonzero velocity, pressure constant 15677841947SLeila Ghaffari // (so density and internal energy will be non-constant), 15777841947SLeila Ghaffari // it should be transported across the domain, but velocity stays constant 15877841947SLeila Ghaffari P = 1.; 15977841947SLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 16077841947SLeila Ghaffari rho = P / (R*T); 16177841947SLeila Ghaffari u[0] = mean_velocity[0]; 16277841947SLeila Ghaffari u[1] = mean_velocity[1]; 16377841947SLeila Ghaffari 16477841947SLeila Ghaffari // Assign exact solution 16577841947SLeila Ghaffari q[0] = rho; 16677841947SLeila Ghaffari q[1] = rho * u[0]; 16777841947SLeila Ghaffari q[2] = rho * u[1]; 16877841947SLeila Ghaffari q[3] = rho * u[2]; 16977841947SLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 17077841947SLeila Ghaffari break; 17132f166c6SLeila Ghaffari case 5: // non-smooth thermal bubble - cylinder 17232f166c6SLeila Ghaffari P = 1.; 17332f166c6SLeila Ghaffari T = 1. - (r < 1. ? S_bubble : 0.); 17432f166c6SLeila Ghaffari rho = P / (R*T); 17532f166c6SLeila Ghaffari u[0] = mean_velocity[0]; 17632f166c6SLeila Ghaffari u[1] = mean_velocity[1]; 17732f166c6SLeila Ghaffari 17832f166c6SLeila Ghaffari // Assign exact solution 17932f166c6SLeila Ghaffari q[0] = rho; 18032f166c6SLeila Ghaffari q[1] = rho * u[0]; 18132f166c6SLeila Ghaffari q[2] = rho * u[1]; 18232f166c6SLeila Ghaffari q[3] = rho * u[2]; 18332f166c6SLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 18432f166c6SLeila Ghaffari break; 18577841947SLeila Ghaffari } 18677841947SLeila Ghaffari // Return 18777841947SLeila Ghaffari return 0; 18877841947SLeila Ghaffari } 18977841947SLeila Ghaffari 19077841947SLeila Ghaffari // ***************************************************************************** 191e6225c47SLeila Ghaffari // Helper function for computing flux Jacobian 192e6225c47SLeila Ghaffari // ***************************************************************************** 193932417b3SJed Brown CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], 194e6225c47SLeila Ghaffari const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 195e6225c47SLeila Ghaffari const CeedScalar gamma) { 196e6225c47SLeila Ghaffari CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 197e6225c47SLeila Ghaffari for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 198e6225c47SLeila Ghaffari for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix 199e6225c47SLeila Ghaffari dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j]; 200e6225c47SLeila Ghaffari for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix 201e6225c47SLeila Ghaffari dF[i][0][k+1] = ((i==k) ? 1. : 0.); 202e6225c47SLeila Ghaffari dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) + 203e6225c47SLeila Ghaffari ((i==k) ? u[j] : 0.) - 204e6225c47SLeila Ghaffari ((i==j) ? u[k] : 0.) * (gamma-1.); 205e6225c47SLeila Ghaffari dF[i][4][k+1] = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) - 206e6225c47SLeila Ghaffari (gamma-1.)*u[i]*u[k]; 207e6225c47SLeila Ghaffari } 208e6225c47SLeila Ghaffari dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.); 209e6225c47SLeila Ghaffari } 210e6225c47SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho); 211e6225c47SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 212e6225c47SLeila Ghaffari } 213e6225c47SLeila Ghaffari } 214e6225c47SLeila Ghaffari 215e6225c47SLeila Ghaffari // ***************************************************************************** 216932417b3SJed Brown // Helper function for computing Tau elements (stabilization constant) 217932417b3SJed Brown // Model from: 218932417b3SJed Brown // Stabilized Methods for Compressible Flows, Hughes et al 2010 219932417b3SJed Brown // 220932417b3SJed Brown // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 221932417b3SJed Brown // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 222932417b3SJed Brown // 223932417b3SJed Brown // Where 224932417b3SJed Brown // c_tau = stabilization constant (0.5 is reported as "optimal") 225932417b3SJed Brown // h[i] = 2 length(dxdX[i]) 226932417b3SJed Brown // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 227932417b3SJed Brown // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 228932417b3SJed Brown // rho(A[i]) = spectral radius of the convective flux Jacobian i, 229932417b3SJed Brown // wave speed in direction i 230932417b3SJed Brown // ***************************************************************************** 231932417b3SJed Brown CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], 232932417b3SJed Brown const CeedScalar dXdx[3][3], const CeedScalar u[3], 233932417b3SJed Brown const CeedScalar sound_speed, const CeedScalar c_tau) { 234ba6664aeSJames Wright for (CeedInt i=0; i<3; i++) { 235932417b3SJed Brown // length of element in direction i 236932417b3SJed Brown CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] + 237932417b3SJed Brown dXdx[2][i]*dXdx[2][i]); 238932417b3SJed Brown // fastest wave in direction i 239932417b3SJed Brown CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 240932417b3SJed Brown Tau_x[i] = c_tau * h / fastest_wave; 241932417b3SJed Brown } 242932417b3SJed Brown } 243932417b3SJed Brown 244932417b3SJed Brown // ***************************************************************************** 24577841947SLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex 24677841947SLeila Ghaffari // ***************************************************************************** 24777841947SLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, 24877841947SLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 24977841947SLeila Ghaffari // Inputs 25077841947SLeila Ghaffari const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 25177841947SLeila Ghaffari 25277841947SLeila Ghaffari // Outputs 25377841947SLeila Ghaffari CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 25477841947SLeila Ghaffari const EulerContext context = (EulerContext)ctx; 25577841947SLeila Ghaffari 25677841947SLeila Ghaffari CeedPragmaSIMD 25777841947SLeila Ghaffari // Quadrature Point Loop 25877841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 25977841947SLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 260e6225c47SLeila Ghaffari CeedScalar q[5] = {0.}; 26177841947SLeila Ghaffari 26277841947SLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx); 26377841947SLeila Ghaffari 26477841947SLeila Ghaffari for (CeedInt j=0; j<5; j++) 26577841947SLeila Ghaffari q0[j][i] = q[j]; 26677841947SLeila Ghaffari } // End of Quadrature Point Loop 26777841947SLeila Ghaffari 26877841947SLeila Ghaffari // Return 26977841947SLeila Ghaffari return 0; 27077841947SLeila Ghaffari } 27177841947SLeila Ghaffari 27277841947SLeila Ghaffari // ***************************************************************************** 27377841947SLeila Ghaffari // This QFunction implements the following formulation of Euler equations 27477841947SLeila Ghaffari // with explicit time stepping method 27577841947SLeila Ghaffari // 27677841947SLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation 27777841947SLeila Ghaffari // form with state variables of density, momentum density, and total 27877841947SLeila Ghaffari // energy density. 27977841947SLeila Ghaffari // 28077841947SLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 28177841947SLeila Ghaffari // rho - Mass Density 28277841947SLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 28377841947SLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 28477841947SLeila Ghaffari // 28577841947SLeila Ghaffari // Euler Equations: 28677841947SLeila Ghaffari // drho/dt + div( U ) = 0 28777841947SLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0 28877841947SLeila Ghaffari // dE/dt + div( (E + P) u ) = 0 28977841947SLeila Ghaffari // 29077841947SLeila Ghaffari // Equation of State: 29177841947SLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2) 29277841947SLeila Ghaffari // 29377841947SLeila Ghaffari // Constants: 29477841947SLeila Ghaffari // cv , Specific heat, constant volume 29577841947SLeila Ghaffari // cp , Specific heat, constant pressure 29677841947SLeila Ghaffari // g , Gravity 29777841947SLeila Ghaffari // gamma = cp / cv, Specific heat ratio 29877841947SLeila Ghaffari // ***************************************************************************** 29977841947SLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, 30077841947SLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 30177841947SLeila Ghaffari // *INDENT-OFF* 30277841947SLeila Ghaffari // Inputs 30377841947SLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 304e6225c47SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 30577841947SLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 30677841947SLeila Ghaffari // Outputs 30777841947SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 30877841947SLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 30977841947SLeila Ghaffari 310e6225c47SLeila Ghaffari EulerContext context = (EulerContext)ctx; 311932417b3SJed Brown const CeedScalar c_tau = context->c_tau; 31277841947SLeila Ghaffari const CeedScalar gamma = 1.4; 31377841947SLeila Ghaffari 31477841947SLeila Ghaffari CeedPragmaSIMD 31577841947SLeila Ghaffari // Quadrature Point Loop 31677841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 31777841947SLeila Ghaffari // *INDENT-OFF* 31877841947SLeila Ghaffari // Setup 31977841947SLeila Ghaffari // -- Interp in 32077841947SLeila Ghaffari const CeedScalar rho = q[0][i]; 32177841947SLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 32277841947SLeila Ghaffari q[2][i] / rho, 32377841947SLeila Ghaffari q[3][i] / rho 32477841947SLeila Ghaffari }; 32577841947SLeila Ghaffari const CeedScalar E = q[4][i]; 326e6225c47SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 327e6225c47SLeila Ghaffari dq[1][0][i], 328e6225c47SLeila Ghaffari dq[2][0][i] 329e6225c47SLeila Ghaffari }; 330e6225c47SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 331e6225c47SLeila Ghaffari dq[1][1][i], 332e6225c47SLeila Ghaffari dq[2][1][i]}, 333e6225c47SLeila Ghaffari {dq[0][2][i], 334e6225c47SLeila Ghaffari dq[1][2][i], 335e6225c47SLeila Ghaffari dq[2][2][i]}, 336e6225c47SLeila Ghaffari {dq[0][3][i], 337e6225c47SLeila Ghaffari dq[1][3][i], 338e6225c47SLeila Ghaffari dq[2][3][i]} 339e6225c47SLeila Ghaffari }; 340e6225c47SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 341e6225c47SLeila Ghaffari dq[1][4][i], 342e6225c47SLeila Ghaffari dq[2][4][i] 343e6225c47SLeila Ghaffari }; 34477841947SLeila Ghaffari // -- Interp-to-Interp q_data 34577841947SLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 34677841947SLeila Ghaffari // -- Interp-to-Grad q_data 34777841947SLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 34877841947SLeila Ghaffari // *INDENT-OFF* 34977841947SLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 35077841947SLeila Ghaffari q_data[2][i], 35177841947SLeila Ghaffari q_data[3][i]}, 35277841947SLeila Ghaffari {q_data[4][i], 35377841947SLeila Ghaffari q_data[5][i], 35477841947SLeila Ghaffari q_data[6][i]}, 35577841947SLeila Ghaffari {q_data[7][i], 35677841947SLeila Ghaffari q_data[8][i], 35777841947SLeila Ghaffari q_data[9][i]} 35877841947SLeila Ghaffari }; 35977841947SLeila Ghaffari // *INDENT-ON* 360e6225c47SLeila Ghaffari // dU/dx 361e6225c47SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 362e6225c47SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 363e6225c47SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 364e6225c47SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 365ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) { 366ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) { 367e6225c47SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 368e6225c47SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 369ba6664aeSJames Wright for (CeedInt l=0; l<3; l++) { 370e6225c47SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 371e6225c47SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 372e6225c47SLeila Ghaffari } 373e6225c47SLeila Ghaffari } 374e6225c47SLeila Ghaffari } 375e6225c47SLeila Ghaffari // Pressure 37677841947SLeila Ghaffari const CeedScalar 37777841947SLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 37877841947SLeila Ghaffari E_internal = E - E_kinetic, 379e6225c47SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 38077841947SLeila Ghaffari 38177841947SLeila Ghaffari // The Physics 38277841947SLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 383ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) { 384e6225c47SLeila Ghaffari v[j][i] = 0.; 385ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 386e6225c47SLeila Ghaffari dv[k][j][i] = 0.; 38777841947SLeila Ghaffari } 38877841947SLeila Ghaffari 38977841947SLeila Ghaffari // -- Density 39077841947SLeila Ghaffari // ---- u rho 391ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 39277841947SLeila Ghaffari dv[j][0][i] += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 39377841947SLeila Ghaffari rho*u[2]*dXdx[j][2]); 39477841947SLeila Ghaffari // -- Momentum 39577841947SLeila Ghaffari // ---- rho (u x u) + P I3 396ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 397ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 398e6225c47SLeila Ghaffari dv[k][j+1][i] += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 399e6225c47SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 400e6225c47SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 40177841947SLeila Ghaffari // -- Total Energy Density 40277841947SLeila Ghaffari // ---- (E + P) u 403ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 40477841947SLeila Ghaffari dv[j][4][i] += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 40577841947SLeila Ghaffari u[2]*dXdx[j][2]); 406e6225c47SLeila Ghaffari 407e6225c47SLeila Ghaffari // --Stabilization terms 408e6225c47SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 409e6225c47SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 410932417b3SJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 411e6225c47SLeila Ghaffari 412e6225c47SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 413e6225c47SLeila Ghaffari CeedScalar dqdx[5][3]; 414ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) { 415e6225c47SLeila Ghaffari dqdx[0][j] = drhodx[j]; 416e6225c47SLeila Ghaffari dqdx[4][j] = dEdx[j]; 417ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 418e6225c47SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 419e6225c47SLeila Ghaffari } 420e6225c47SLeila Ghaffari 421e6225c47SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 422e6225c47SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 423ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 424ba6664aeSJames Wright for (CeedInt k=0; k<5; k++) 425ba6664aeSJames Wright for (CeedInt l=0; l<5; l++) 426e6225c47SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 427e6225c47SLeila Ghaffari 428932417b3SJed Brown // Stabilization 429932417b3SJed Brown // -- Tau elements 430932417b3SJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 431932417b3SJed Brown CeedScalar Tau_x[3] = {0.}; 432932417b3SJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 433e6225c47SLeila Ghaffari 434932417b3SJed Brown // -- Stabilization method: none or SU 43588626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 436e6225c47SLeila Ghaffari switch (context->stabilization) { 437e6225c47SLeila Ghaffari case 0: // Galerkin 438e6225c47SLeila Ghaffari break; 439e6225c47SLeila Ghaffari case 1: // SU 440ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 441ba6664aeSJames Wright for (CeedInt k=0; k<5; k++) 442ba6664aeSJames Wright for (CeedInt l=0; l<5; l++) 44388626eedSJames Wright stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 444e6225c47SLeila Ghaffari 445ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) 446ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 447e6225c47SLeila Ghaffari dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] + 448e6225c47SLeila Ghaffari stab[j][1] * dXdx[k][1] + 449e6225c47SLeila Ghaffari stab[j][2] * dXdx[k][2]); 450e6225c47SLeila Ghaffari break; 451e6225c47SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 452e6225c47SLeila Ghaffari break; 453e6225c47SLeila Ghaffari } 454e6225c47SLeila Ghaffari 45577841947SLeila Ghaffari } // End Quadrature Point Loop 45677841947SLeila Ghaffari 45777841947SLeila Ghaffari // Return 45877841947SLeila Ghaffari return 0; 45977841947SLeila Ghaffari } 46077841947SLeila Ghaffari 46177841947SLeila Ghaffari // ***************************************************************************** 46277841947SLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above) 46377841947SLeila Ghaffari // with implicit time stepping method 46477841947SLeila Ghaffari // 46577841947SLeila Ghaffari // ***************************************************************************** 46677841947SLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, 46777841947SLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 46877841947SLeila Ghaffari // *INDENT-OFF* 46977841947SLeila Ghaffari // Inputs 47077841947SLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 471e6225c47SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 47277841947SLeila Ghaffari (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 47377841947SLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 47477841947SLeila Ghaffari // Outputs 47577841947SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 47677841947SLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 47777841947SLeila Ghaffari 478e6225c47SLeila Ghaffari EulerContext context = (EulerContext)ctx; 479932417b3SJed Brown const CeedScalar c_tau = context->c_tau; 48077841947SLeila Ghaffari const CeedScalar gamma = 1.4; 48177841947SLeila Ghaffari 48277841947SLeila Ghaffari CeedPragmaSIMD 48377841947SLeila Ghaffari // Quadrature Point Loop 48477841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 48577841947SLeila Ghaffari // *INDENT-OFF* 48677841947SLeila Ghaffari // Setup 48777841947SLeila Ghaffari // -- Interp in 48877841947SLeila Ghaffari const CeedScalar rho = q[0][i]; 48977841947SLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 49077841947SLeila Ghaffari q[2][i] / rho, 49177841947SLeila Ghaffari q[3][i] / rho 49277841947SLeila Ghaffari }; 49377841947SLeila Ghaffari const CeedScalar E = q[4][i]; 494e6225c47SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 495e6225c47SLeila Ghaffari dq[1][0][i], 496e6225c47SLeila Ghaffari dq[2][0][i] 497e6225c47SLeila Ghaffari }; 498e6225c47SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 499e6225c47SLeila Ghaffari dq[1][1][i], 500e6225c47SLeila Ghaffari dq[2][1][i]}, 501e6225c47SLeila Ghaffari {dq[0][2][i], 502e6225c47SLeila Ghaffari dq[1][2][i], 503e6225c47SLeila Ghaffari dq[2][2][i]}, 504e6225c47SLeila Ghaffari {dq[0][3][i], 505e6225c47SLeila Ghaffari dq[1][3][i], 506e6225c47SLeila Ghaffari dq[2][3][i]} 507e6225c47SLeila Ghaffari }; 508e6225c47SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 509e6225c47SLeila Ghaffari dq[1][4][i], 510e6225c47SLeila Ghaffari dq[2][4][i] 511e6225c47SLeila Ghaffari }; 51277841947SLeila Ghaffari // -- Interp-to-Interp q_data 51377841947SLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 51477841947SLeila Ghaffari // -- Interp-to-Grad q_data 51577841947SLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 51677841947SLeila Ghaffari // *INDENT-OFF* 51777841947SLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 51877841947SLeila Ghaffari q_data[2][i], 51977841947SLeila Ghaffari q_data[3][i]}, 52077841947SLeila Ghaffari {q_data[4][i], 52177841947SLeila Ghaffari q_data[5][i], 52277841947SLeila Ghaffari q_data[6][i]}, 52377841947SLeila Ghaffari {q_data[7][i], 52477841947SLeila Ghaffari q_data[8][i], 52577841947SLeila Ghaffari q_data[9][i]} 52677841947SLeila Ghaffari }; 52777841947SLeila Ghaffari // *INDENT-ON* 528e6225c47SLeila Ghaffari // dU/dx 529e6225c47SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 530e6225c47SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 531e6225c47SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 532e6225c47SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 533ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) { 534ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) { 535e6225c47SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 536e6225c47SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 537ba6664aeSJames Wright for (CeedInt l=0; l<3; l++) { 538e6225c47SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 539e6225c47SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 540e6225c47SLeila Ghaffari } 541e6225c47SLeila Ghaffari } 542e6225c47SLeila Ghaffari } 54377841947SLeila Ghaffari const CeedScalar 54477841947SLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 54577841947SLeila Ghaffari E_internal = E - E_kinetic, 546e6225c47SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 54777841947SLeila Ghaffari 54877841947SLeila Ghaffari // The Physics 54977841947SLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 550ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) { 551e6225c47SLeila Ghaffari v[j][i] = 0.; 552ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 553e6225c47SLeila Ghaffari dv[k][j][i] = 0.; 55477841947SLeila Ghaffari } 55577841947SLeila Ghaffari //-----mass matrix 556ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) 55777841947SLeila Ghaffari v[j][i] += wdetJ*q_dot[j][i]; 55877841947SLeila Ghaffari 55977841947SLeila Ghaffari // -- Density 56077841947SLeila Ghaffari // ---- u rho 561ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 56277841947SLeila Ghaffari dv[j][0][i] -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 56377841947SLeila Ghaffari rho*u[2]*dXdx[j][2]); 56477841947SLeila Ghaffari // -- Momentum 56577841947SLeila Ghaffari // ---- rho (u x u) + P I3 566ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 567ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 568e6225c47SLeila Ghaffari dv[k][j+1][i] -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 569e6225c47SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 570e6225c47SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 57177841947SLeila Ghaffari // -- Total Energy Density 57277841947SLeila Ghaffari // ---- (E + P) u 573ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 57477841947SLeila Ghaffari dv[j][4][i] -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 57577841947SLeila Ghaffari u[2]*dXdx[j][2]); 576e6225c47SLeila Ghaffari 577e6225c47SLeila Ghaffari // -- Stabilization terms 578e6225c47SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 579e6225c47SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 580932417b3SJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 581e6225c47SLeila Ghaffari 582e6225c47SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 583e6225c47SLeila Ghaffari CeedScalar dqdx[5][3]; 584ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) { 585e6225c47SLeila Ghaffari dqdx[0][j] = drhodx[j]; 586e6225c47SLeila Ghaffari dqdx[4][j] = dEdx[j]; 587ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 588e6225c47SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 589e6225c47SLeila Ghaffari } 590e6225c47SLeila Ghaffari 591e6225c47SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 592e6225c47SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 593ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 594ba6664aeSJames Wright for (CeedInt k=0; k<5; k++) 595ba6664aeSJames Wright for (CeedInt l=0; l<5; l++) 596e6225c47SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 597e6225c47SLeila Ghaffari 598e6225c47SLeila Ghaffari // ---- Strong residual 599e6225c47SLeila Ghaffari CeedScalar strong_res[5]; 600ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) 601e6225c47SLeila Ghaffari strong_res[j] = q_dot[j][i] + strong_conv[j]; 602e6225c47SLeila Ghaffari 603932417b3SJed Brown // Stabilization 604932417b3SJed Brown // -- Tau elements 605932417b3SJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 606932417b3SJed Brown CeedScalar Tau_x[3] = {0.}; 607932417b3SJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 608e6225c47SLeila Ghaffari 609932417b3SJed Brown // -- Stabilization method: none, SU, or SUPG 61088626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 611e6225c47SLeila Ghaffari switch (context->stabilization) { 612e6225c47SLeila Ghaffari case 0: // Galerkin 613e6225c47SLeila Ghaffari break; 614e6225c47SLeila Ghaffari case 1: // SU 615ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 616ba6664aeSJames Wright for (CeedInt k=0; k<5; k++) 617ba6664aeSJames Wright for (CeedInt l=0; l<5; l++) 61888626eedSJames Wright stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 619e6225c47SLeila Ghaffari 620ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) 621ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 622e6225c47SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 623e6225c47SLeila Ghaffari stab[j][1] * dXdx[k][1] + 624e6225c47SLeila Ghaffari stab[j][2] * dXdx[k][2]); 625e6225c47SLeila Ghaffari break; 626e6225c47SLeila Ghaffari case 2: // SUPG 627ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 628ba6664aeSJames Wright for (CeedInt k=0; k<5; k++) 629ba6664aeSJames Wright for (CeedInt l=0; l<5; l++) 63088626eedSJames Wright stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l]; 631e6225c47SLeila Ghaffari 632ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) 633ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 634e6225c47SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 635e6225c47SLeila Ghaffari stab[j][1] * dXdx[k][1] + 636e6225c47SLeila Ghaffari stab[j][2] * dXdx[k][2]); 637e6225c47SLeila Ghaffari break; 638e6225c47SLeila Ghaffari } 63977841947SLeila Ghaffari } // End Quadrature Point Loop 64077841947SLeila Ghaffari 64177841947SLeila Ghaffari // Return 64277841947SLeila Ghaffari return 0; 64377841947SLeila Ghaffari } 64477841947SLeila Ghaffari // ***************************************************************************** 6452fe7aee7SLeila Ghaffari // This QFunction sets the inflow boundary conditions for 6462fe7aee7SLeila Ghaffari // the traveling vortex problem. 64777841947SLeila Ghaffari // 64877841947SLeila Ghaffari // Prescribed T_inlet and P_inlet are converted to conservative variables 64977841947SLeila Ghaffari // and applied weakly. 65077841947SLeila Ghaffari // 65177841947SLeila Ghaffari // ***************************************************************************** 6522fe7aee7SLeila Ghaffari CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, 65377841947SLeila Ghaffari const CeedScalar *const *in, 65477841947SLeila Ghaffari CeedScalar *const *out) { 65577841947SLeila Ghaffari // *INDENT-OFF* 65677841947SLeila Ghaffari // Inputs 657*e8b03feeSJames Wright const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 65877841947SLeila Ghaffari // Outputs 65977841947SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 66077841947SLeila Ghaffari // *INDENT-ON* 66177841947SLeila Ghaffari EulerContext context = (EulerContext)ctx; 66277841947SLeila Ghaffari const int euler_test = context->euler_test; 66377841947SLeila Ghaffari const bool implicit = context->implicit; 66477841947SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 66577841947SLeila Ghaffari const CeedScalar cv = 2.5; 66677841947SLeila Ghaffari const CeedScalar R = 1.; 66777841947SLeila Ghaffari CeedScalar T_inlet; 66877841947SLeila Ghaffari CeedScalar P_inlet; 66977841947SLeila Ghaffari 67077841947SLeila Ghaffari // For test cases 1 and 3 the background velocity is zero 67177841947SLeila Ghaffari if (euler_test == 1 || euler_test == 3) 67277841947SLeila Ghaffari for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.; 67377841947SLeila Ghaffari 67477841947SLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4 67577841947SLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 67677841947SLeila Ghaffari else T_inlet = P_inlet = 1.; 67777841947SLeila Ghaffari 67877841947SLeila Ghaffari CeedPragmaSIMD 67977841947SLeila Ghaffari // Quadrature Point Loop 68077841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 68177841947SLeila Ghaffari // Setup 68277841947SLeila Ghaffari // -- Interp-to-Interp q_data 68377841947SLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 68477841947SLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 68577841947SLeila Ghaffari // We can effect this by swapping the sign on this weight 68677841947SLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 6872fe7aee7SLeila Ghaffari // ---- Normal vect 68877841947SLeila Ghaffari const CeedScalar norm[3] = {q_data_sur[1][i], 68977841947SLeila Ghaffari q_data_sur[2][i], 69077841947SLeila Ghaffari q_data_sur[3][i] 69177841947SLeila Ghaffari }; 69277841947SLeila Ghaffari 69377841947SLeila Ghaffari // face_normal = Normal vector of the face 69477841947SLeila Ghaffari const CeedScalar face_normal = norm[0]*mean_velocity[0] + 69577841947SLeila Ghaffari norm[1]*mean_velocity[1] + 69677841947SLeila Ghaffari norm[2]*mean_velocity[2]; 69777841947SLeila Ghaffari // The Physics 69877841947SLeila Ghaffari // Zero v so all future terms can safely sum into it 699ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) v[j][i] = 0.; 70077841947SLeila Ghaffari 70177841947SLeila Ghaffari // Implementing in/outflow BCs 7022fe7aee7SLeila Ghaffari if (face_normal > 0) { 70377841947SLeila Ghaffari } else { // inflow 70477841947SLeila Ghaffari const CeedScalar rho_inlet = P_inlet/(R*T_inlet); 70577841947SLeila Ghaffari const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] + 70677841947SLeila Ghaffari mean_velocity[1]*mean_velocity[1]) / 2.; 70777841947SLeila Ghaffari // incoming total energy 70877841947SLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 70977841947SLeila Ghaffari 71077841947SLeila Ghaffari // The Physics 71177841947SLeila Ghaffari // -- Density 71277841947SLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal; 71377841947SLeila Ghaffari 71477841947SLeila Ghaffari // -- Momentum 715ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 71677841947SLeila Ghaffari v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] + 71777841947SLeila Ghaffari norm[j] * P_inlet); 71877841947SLeila Ghaffari 71977841947SLeila Ghaffari // -- Total Energy Density 72077841947SLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 72177841947SLeila Ghaffari } 72277841947SLeila Ghaffari 72377841947SLeila Ghaffari } // End Quadrature Point Loop 72477841947SLeila Ghaffari return 0; 72577841947SLeila Ghaffari } 72677841947SLeila Ghaffari 72777841947SLeila Ghaffari // ***************************************************************************** 72855e76554SLeila Ghaffari // This QFunction sets the outflow boundary conditions for 72955e76554SLeila Ghaffari // the Euler solver. 73055e76554SLeila Ghaffari // 73155e76554SLeila Ghaffari // Outflow BCs: 73255e76554SLeila Ghaffari // The validity of the weak form of the governing equations is 73355e76554SLeila Ghaffari // extended to the outflow. 73455e76554SLeila Ghaffari // 73555e76554SLeila Ghaffari // ***************************************************************************** 73655e76554SLeila Ghaffari CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, 73755e76554SLeila Ghaffari const CeedScalar *const *in, 73855e76554SLeila Ghaffari CeedScalar *const *out) { 73955e76554SLeila Ghaffari // *INDENT-OFF* 74055e76554SLeila Ghaffari // Inputs 74155e76554SLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 742*e8b03feeSJames Wright (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 74355e76554SLeila Ghaffari // Outputs 74455e76554SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 74555e76554SLeila Ghaffari // *INDENT-ON* 74655e76554SLeila Ghaffari EulerContext context = (EulerContext)ctx; 74755e76554SLeila Ghaffari const bool implicit = context->implicit; 74855e76554SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 74955e76554SLeila Ghaffari 75055e76554SLeila Ghaffari const CeedScalar gamma = 1.4; 75155e76554SLeila Ghaffari 75255e76554SLeila Ghaffari CeedPragmaSIMD 75355e76554SLeila Ghaffari // Quadrature Point Loop 75455e76554SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 75555e76554SLeila Ghaffari // Setup 75655e76554SLeila Ghaffari // -- Interp in 75755e76554SLeila Ghaffari const CeedScalar rho = q[0][i]; 75855e76554SLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 75955e76554SLeila Ghaffari q[2][i] / rho, 76055e76554SLeila Ghaffari q[3][i] / rho 76155e76554SLeila Ghaffari }; 76255e76554SLeila Ghaffari const CeedScalar E = q[4][i]; 76355e76554SLeila Ghaffari 76455e76554SLeila Ghaffari // -- Interp-to-Interp q_data 76555e76554SLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 76655e76554SLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 76755e76554SLeila Ghaffari // We can effect this by swapping the sign on this weight 76855e76554SLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 76955e76554SLeila Ghaffari // ---- Normal vectors 77055e76554SLeila Ghaffari const CeedScalar norm[3] = {q_data_sur[1][i], 77155e76554SLeila Ghaffari q_data_sur[2][i], 77255e76554SLeila Ghaffari q_data_sur[3][i] 77355e76554SLeila Ghaffari }; 77455e76554SLeila Ghaffari 77555e76554SLeila Ghaffari // face_normal = Normal vector of the face 77655e76554SLeila Ghaffari const CeedScalar face_normal = norm[0]*mean_velocity[0] + 77755e76554SLeila Ghaffari norm[1]*mean_velocity[1] + 77855e76554SLeila Ghaffari norm[2]*mean_velocity[2]; 77955e76554SLeila Ghaffari // The Physics 78055e76554SLeila Ghaffari // Zero v so all future terms can safely sum into it 781ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) v[j][i] = 0; 78255e76554SLeila Ghaffari 78355e76554SLeila Ghaffari // Implementing in/outflow BCs 78455e76554SLeila Ghaffari if (face_normal > 0) { // outflow 78555e76554SLeila Ghaffari const CeedScalar E_kinetic = (u[0]*u[0] + u[1]*u[1]) / 2.; 78655e76554SLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 78755e76554SLeila Ghaffari const CeedScalar u_normal = norm[0]*u[0] + norm[1]*u[1] + 78855e76554SLeila Ghaffari norm[2]*u[2]; // Normal velocity 78955e76554SLeila Ghaffari // The Physics 79055e76554SLeila Ghaffari // -- Density 79155e76554SLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal; 79255e76554SLeila Ghaffari 79355e76554SLeila Ghaffari // -- Momentum 794ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 79555e76554SLeila Ghaffari v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + norm[j] * P); 79655e76554SLeila Ghaffari 79755e76554SLeila Ghaffari // -- Total Energy Density 79855e76554SLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P); 79955e76554SLeila Ghaffari } 80055e76554SLeila Ghaffari } // End Quadrature Point Loop 80155e76554SLeila Ghaffari return 0; 80255e76554SLeila Ghaffari } 80355e76554SLeila Ghaffari 80455e76554SLeila Ghaffari // ***************************************************************************** 80577841947SLeila Ghaffari 80677841947SLeila Ghaffari #endif // eulervortex_h 807