15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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: 13ea61e9acSJeremy L Thompson // On the Order of Accuracy and Numerical Performance of Two Classes of Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011). 1488b783a1SJames Wright #include <ceed.h> 15c9c2c079SJeremy L Thompson #include <math.h> 162b730f8bSJeremy L Thompson 1713fa47b2SJames Wright #include "utils.h" 1877841947SLeila Ghaffari 1977841947SLeila Ghaffari typedef struct EulerContext_ *EulerContext; 2077841947SLeila Ghaffari struct EulerContext_ { 2177841947SLeila Ghaffari CeedScalar center[3]; 2277841947SLeila Ghaffari CeedScalar curr_time; 2377841947SLeila Ghaffari CeedScalar vortex_strength; 24932417b3SJed Brown CeedScalar c_tau; 2577841947SLeila Ghaffari CeedScalar mean_velocity[3]; 2677841947SLeila Ghaffari bool implicit; 27e6225c47SLeila Ghaffari int euler_test; 28e6225c47SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 2977841947SLeila Ghaffari }; 3077841947SLeila Ghaffari 3177841947SLeila Ghaffari // ***************************************************************************** 3277841947SLeila Ghaffari // This function sets the initial conditions 3377841947SLeila Ghaffari // 3477841947SLeila Ghaffari // Temperature: 3577841947SLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 3677841947SLeila Ghaffari // Density: 3777841947SLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1)) 3877841947SLeila Ghaffari // Pressure: 3977841947SLeila Ghaffari // P = rho * T 4077841947SLeila Ghaffari // Velocity: 4177841947SLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 4277841947SLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 4377841947SLeila Ghaffari // Velocity/Momentum Density: 4477841947SLeila Ghaffari // Ui = rho ui 4577841947SLeila Ghaffari // Total Energy: 4677841947SLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2 4777841947SLeila Ghaffari // 4877841947SLeila Ghaffari // Constants: 4977841947SLeila Ghaffari // cv , Specific heat, constant volume 5077841947SLeila Ghaffari // cp , Specific heat, constant pressure 5177841947SLeila Ghaffari // vortex_strength , Strength of vortex 5277841947SLeila Ghaffari // center , Location of bubble center 5377841947SLeila Ghaffari // gamma = cp / cv, Specific heat ratio 5477841947SLeila Ghaffari // 5577841947SLeila Ghaffari // ***************************************************************************** 5677841947SLeila Ghaffari 5777841947SLeila Ghaffari // ***************************************************************************** 58ea61e9acSJeremy L Thompson // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling 59ea61e9acSJeremy L Thompson // vortex 6077841947SLeila Ghaffari // ***************************************************************************** 612b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) { 6277841947SLeila Ghaffari // Context 6377841947SLeila Ghaffari const EulerContext context = (EulerContext)ctx; 6477841947SLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength; 6577841947SLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain 6677841947SLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity; 6777841947SLeila Ghaffari 6877841947SLeila Ghaffari // Setup 6977841947SLeila Ghaffari const CeedScalar gamma = 1.4; 7077841947SLeila Ghaffari const CeedScalar cv = 2.5; 7177841947SLeila Ghaffari const CeedScalar R = 1.; 7277841947SLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates 7377841947SLeila Ghaffari // Vortex center 7477841947SLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time; 7577841947SLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time; 7677841947SLeila Ghaffari 7777841947SLeila Ghaffari const CeedScalar x0 = x - xc; 7877841947SLeila Ghaffari const CeedScalar y0 = y - yc; 7977841947SLeila Ghaffari const CeedScalar r = sqrt(x0 * x0 + y0 * y0); 8077841947SLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI); 812b730f8bSJeremy L Thompson const CeedScalar delta_T = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI); 8277841947SLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 832b730f8bSJeremy L Thompson const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI); 8477841947SLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.}; 8577841947SLeila Ghaffari 8677841947SLeila Ghaffari // Initial Conditions 8777841947SLeila Ghaffari switch (context->euler_test) { 8877841947SLeila Ghaffari case 0: // Traveling vortex 8977841947SLeila Ghaffari T = 1 + delta_T; 9077841947SLeila Ghaffari // P = rho * T 9177841947SLeila Ghaffari // P = S * rho^gamma 9277841947SLeila Ghaffari // Solve for rho, then substitute for P 93e6225c47SLeila Ghaffari rho = pow(T / S_vortex, 1 / (gamma - 1.)); 9477841947SLeila Ghaffari P = rho * T; 9577841947SLeila Ghaffari u[0] = mean_velocity[0] - C * y0; 9677841947SLeila Ghaffari u[1] = mean_velocity[1] + C * x0; 9777841947SLeila Ghaffari 9877841947SLeila Ghaffari // Assign exact solution 9977841947SLeila Ghaffari q[0] = rho; 10077841947SLeila Ghaffari q[1] = rho * u[0]; 10177841947SLeila Ghaffari q[2] = rho * u[1]; 10277841947SLeila Ghaffari q[3] = rho * u[2]; 10377841947SLeila Ghaffari q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.; 10477841947SLeila Ghaffari break; 10577841947SLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant 10677841947SLeila Ghaffari rho = 1.; 10777841947SLeila Ghaffari E = 2.; 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] = E; 11577841947SLeila Ghaffari break; 11677841947SLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant 11777841947SLeila Ghaffari rho = 1.; 11877841947SLeila Ghaffari E = 2.; 11977841947SLeila Ghaffari u[0] = mean_velocity[0]; 12077841947SLeila Ghaffari u[1] = mean_velocity[1]; 12177841947SLeila Ghaffari 12277841947SLeila Ghaffari // Assign exact solution 12377841947SLeila Ghaffari q[0] = rho; 12477841947SLeila Ghaffari q[1] = rho * u[0]; 12577841947SLeila Ghaffari q[2] = rho * u[1]; 12677841947SLeila Ghaffari q[3] = rho * u[2]; 12777841947SLeila Ghaffari q[4] = E; 12877841947SLeila Ghaffari break; 129ea61e9acSJeremy L Thompson case 3: // Velocity zero, pressure constant (so density and internal energy will be non-constant), but the velocity should stay zero and the 130ea61e9acSJeremy L Thompson // bubble won't diffuse 13177841947SLeila Ghaffari // (for Euler, where there is no thermal conductivity) 13277841947SLeila Ghaffari P = 1.; 13377841947SLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r); 13477841947SLeila Ghaffari rho = P / (R * T); 13577841947SLeila Ghaffari 13677841947SLeila Ghaffari // Assign exact solution 13777841947SLeila Ghaffari q[0] = rho; 13877841947SLeila Ghaffari q[1] = rho * u[0]; 13977841947SLeila Ghaffari q[2] = rho * u[1]; 14077841947SLeila Ghaffari q[3] = rho * u[2]; 14177841947SLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 14277841947SLeila Ghaffari break; 143ea61e9acSJeremy L Thompson case 4: // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant), 144ea61e9acSJeremy L Thompson // It should be transported across the domain, but velocity stays constant 14577841947SLeila Ghaffari P = 1.; 14677841947SLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r); 14777841947SLeila Ghaffari rho = P / (R * T); 14877841947SLeila Ghaffari u[0] = mean_velocity[0]; 14977841947SLeila Ghaffari u[1] = mean_velocity[1]; 15077841947SLeila Ghaffari 15177841947SLeila Ghaffari // Assign exact solution 15277841947SLeila Ghaffari q[0] = rho; 15377841947SLeila Ghaffari q[1] = rho * u[0]; 15477841947SLeila Ghaffari q[2] = rho * u[1]; 15577841947SLeila Ghaffari q[3] = rho * u[2]; 15677841947SLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 15777841947SLeila Ghaffari break; 15832f166c6SLeila Ghaffari case 5: // non-smooth thermal bubble - cylinder 15932f166c6SLeila Ghaffari P = 1.; 16032f166c6SLeila Ghaffari T = 1. - (r < 1. ? S_bubble : 0.); 16132f166c6SLeila Ghaffari rho = P / (R * T); 16232f166c6SLeila Ghaffari u[0] = mean_velocity[0]; 16332f166c6SLeila Ghaffari u[1] = mean_velocity[1]; 16432f166c6SLeila Ghaffari 16532f166c6SLeila Ghaffari // Assign exact solution 16632f166c6SLeila Ghaffari q[0] = rho; 16732f166c6SLeila Ghaffari q[1] = rho * u[0]; 16832f166c6SLeila Ghaffari q[2] = rho * u[1]; 16932f166c6SLeila Ghaffari q[3] = rho * u[2]; 17032f166c6SLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 17132f166c6SLeila Ghaffari break; 17277841947SLeila Ghaffari } 17377841947SLeila Ghaffari return 0; 17477841947SLeila Ghaffari } 17577841947SLeila Ghaffari 17677841947SLeila Ghaffari // ***************************************************************************** 177e6225c47SLeila Ghaffari // Helper function for computing flux Jacobian 178e6225c47SLeila Ghaffari // ***************************************************************************** 1792b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 180e6225c47SLeila Ghaffari const CeedScalar gamma) { 181e6225c47SLeila Ghaffari CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2]; // Velocity square 182e6225c47SLeila Ghaffari for (CeedInt i = 0; i < 3; i++) { // Jacobian matrices for 3 directions 183e6225c47SLeila Ghaffari for (CeedInt j = 0; j < 3; j++) { // Rows of each Jacobian matrix 184e6225c47SLeila Ghaffari dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j]; 185e6225c47SLeila Ghaffari for (CeedInt k = 0; k < 3; k++) { // Columns of each Jacobian matrix 186e6225c47SLeila Ghaffari dF[i][0][k + 1] = ((i == k) ? 1. : 0.); 1872b730f8bSJeremy L Thompson dF[i][j + 1][k + 1] = ((j == k) ? u[i] : 0.) + ((i == k) ? u[j] : 0.) - ((i == j) ? u[k] : 0.) * (gamma - 1.); 1882b730f8bSJeremy L Thompson dF[i][4][k + 1] = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k]; 189e6225c47SLeila Ghaffari } 190e6225c47SLeila Ghaffari dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.); 191e6225c47SLeila Ghaffari } 192e6225c47SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho); 193e6225c47SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 194e6225c47SLeila Ghaffari } 195e6225c47SLeila Ghaffari } 196e6225c47SLeila Ghaffari 197e6225c47SLeila Ghaffari // ***************************************************************************** 198932417b3SJed Brown // Helper function for computing Tau elements (stabilization constant) 199932417b3SJed Brown // Model from: 200932417b3SJed Brown // Stabilized Methods for Compressible Flows, Hughes et al 2010 201932417b3SJed Brown // 202932417b3SJed Brown // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 203932417b3SJed Brown // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 204932417b3SJed Brown // 205932417b3SJed Brown // Where 206932417b3SJed Brown // c_tau = stabilization constant (0.5 is reported as "optimal") 207932417b3SJed Brown // h[i] = 2 length(dxdX[i]) 208932417b3SJed Brown // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 209932417b3SJed Brown // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 210ea61e9acSJeremy L Thompson // rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i 211932417b3SJed Brown // ***************************************************************************** 2122b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], const CeedScalar dXdx[3][3], const CeedScalar u[3], const CeedScalar sound_speed, 2132b730f8bSJeremy L Thompson const CeedScalar c_tau) { 214ba6664aeSJames Wright for (CeedInt i = 0; i < 3; i++) { 215932417b3SJed Brown // length of element in direction i 2162b730f8bSJeremy L Thompson CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]); 217932417b3SJed Brown // fastest wave in direction i 218932417b3SJed Brown CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 219932417b3SJed Brown Tau_x[i] = c_tau * h / fastest_wave; 220932417b3SJed Brown } 221932417b3SJed Brown } 222932417b3SJed Brown 223932417b3SJed Brown // ***************************************************************************** 22477841947SLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex 22577841947SLeila Ghaffari // ***************************************************************************** 2262b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 22777841947SLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 22877841947SLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 229*f0b01153SJames Wright 23077841947SLeila Ghaffari const EulerContext context = (EulerContext)ctx; 23177841947SLeila Ghaffari 23246603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 23377841947SLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 234e6225c47SLeila Ghaffari CeedScalar q[5] = {0.}; 23577841947SLeila Ghaffari 23677841947SLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx); 2372b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 238*f0b01153SJames Wright } 23977841947SLeila Ghaffari return 0; 24077841947SLeila Ghaffari } 24177841947SLeila Ghaffari 24277841947SLeila Ghaffari // ***************************************************************************** 243ea61e9acSJeremy L Thompson // This QFunction implements the following formulation of Euler equations with explicit time stepping method 24477841947SLeila Ghaffari // 245ea61e9acSJeremy L Thompson // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density. 24677841947SLeila Ghaffari // 24777841947SLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 24877841947SLeila Ghaffari // rho - Mass Density 24977841947SLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 25077841947SLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 25177841947SLeila Ghaffari // 25277841947SLeila Ghaffari // Euler Equations: 25377841947SLeila Ghaffari // drho/dt + div( U ) = 0 25477841947SLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0 25577841947SLeila Ghaffari // dE/dt + div( (E + P) u ) = 0 25677841947SLeila Ghaffari // 25777841947SLeila Ghaffari // Equation of State: 25877841947SLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2) 25977841947SLeila Ghaffari // 26077841947SLeila Ghaffari // Constants: 26177841947SLeila Ghaffari // cv , Specific heat, constant volume 26277841947SLeila Ghaffari // cp , Specific heat, constant pressure 26377841947SLeila Ghaffari // g , Gravity 26477841947SLeila Ghaffari // gamma = cp / cv, Specific heat ratio 26577841947SLeila Ghaffari // ***************************************************************************** 2662b730f8bSJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 26746603fc5SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 26846603fc5SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 269f3e15844SJames Wright const CeedScalar(*q_data) = in[2]; 27046603fc5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 27146603fc5SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 27277841947SLeila Ghaffari 273e6225c47SLeila Ghaffari EulerContext context = (EulerContext)ctx; 274932417b3SJed Brown const CeedScalar c_tau = context->c_tau; 27577841947SLeila Ghaffari const CeedScalar gamma = 1.4; 27677841947SLeila Ghaffari 27746603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 27877841947SLeila Ghaffari // Setup 27977841947SLeila Ghaffari // -- Interp in 28077841947SLeila Ghaffari const CeedScalar rho = q[0][i]; 2812b730f8bSJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 28277841947SLeila Ghaffari const CeedScalar E = q[4][i]; 2832b730f8bSJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 2842b730f8bSJeremy L Thompson const CeedScalar dU[3][3] = { 2852b730f8bSJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 2862b730f8bSJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 2872b730f8bSJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 288e6225c47SLeila Ghaffari }; 2892b730f8bSJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 290f3e15844SJames Wright CeedScalar wdetJ, dXdx[3][3]; 291f3e15844SJames Wright QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx); 292e6225c47SLeila Ghaffari // dU/dx 293e6225c47SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 294e6225c47SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 295e6225c47SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 296e6225c47SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 297ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) { 298ba6664aeSJames Wright for (CeedInt k = 0; k < 3; k++) { 299e6225c47SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 300e6225c47SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 301ba6664aeSJames Wright for (CeedInt l = 0; l < 3; l++) { 302e6225c47SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 303e6225c47SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 304e6225c47SLeila Ghaffari } 305e6225c47SLeila Ghaffari } 306e6225c47SLeila Ghaffari } 307e6225c47SLeila Ghaffari // Pressure 3082b730f8bSJeremy L Thompson const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic, 309e6225c47SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 31077841947SLeila Ghaffari 31177841947SLeila Ghaffari // The Physics 31277841947SLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 313ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) { 314e6225c47SLeila Ghaffari v[j][i] = 0.; 3152b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 31677841947SLeila Ghaffari } 31777841947SLeila Ghaffari 31877841947SLeila Ghaffari // -- Density 31977841947SLeila Ghaffari // ---- u rho 3202b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][0][i] += wdetJ * (rho * u[0] * dXdx[j][0] + rho * u[1] * dXdx[j][1] + rho * u[2] * dXdx[j][2]); 32177841947SLeila Ghaffari // -- Momentum 32277841947SLeila Ghaffari // ---- rho (u x u) + P I3 3232b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3242b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { 3252b730f8bSJeremy L Thompson dv[k][j + 1][i] += wdetJ * ((rho * u[j] * u[0] + (j == 0 ? P : 0.)) * dXdx[k][0] + (rho * u[j] * u[1] + (j == 1 ? P : 0.)) * dXdx[k][1] + 326e6225c47SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 3272b730f8bSJeremy L Thompson } 3282b730f8bSJeremy L Thompson } 32977841947SLeila Ghaffari // -- Total Energy Density 33077841947SLeila Ghaffari // ---- (E + P) u 3312b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][4][i] += wdetJ * (E + P) * (u[0] * dXdx[j][0] + u[1] * dXdx[j][1] + u[2] * dXdx[j][2]); 332e6225c47SLeila Ghaffari 333e6225c47SLeila Ghaffari // --Stabilization terms 334e6225c47SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 335e6225c47SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 336932417b3SJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 337e6225c47SLeila Ghaffari 338e6225c47SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 339e6225c47SLeila Ghaffari CeedScalar dqdx[5][3]; 340ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) { 341e6225c47SLeila Ghaffari dqdx[0][j] = drhodx[j]; 342e6225c47SLeila Ghaffari dqdx[4][j] = dEdx[j]; 3432b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 344e6225c47SLeila Ghaffari } 345e6225c47SLeila Ghaffari 346e6225c47SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 347e6225c47SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 3482b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3492b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 3502b730f8bSJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 3512b730f8bSJeremy L Thompson } 3522b730f8bSJeremy L Thompson } 353e6225c47SLeila Ghaffari 354932417b3SJed Brown // Stabilization 355932417b3SJed Brown // -- Tau elements 356932417b3SJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 357932417b3SJed Brown CeedScalar Tau_x[3] = {0.}; 358932417b3SJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 359e6225c47SLeila Ghaffari 360932417b3SJed Brown // -- Stabilization method: none or SU 36188626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 362e6225c47SLeila Ghaffari switch (context->stabilization) { 363e6225c47SLeila Ghaffari case 0: // Galerkin 364e6225c47SLeila Ghaffari break; 365e6225c47SLeila Ghaffari case 1: // SU 3662b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3672b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 3682b730f8bSJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 3692b730f8bSJeremy L Thompson } 3702b730f8bSJeremy L Thompson } 371e6225c47SLeila Ghaffari 3722b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 3732b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] -= wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]); 3742b730f8bSJeremy L Thompson } 375e6225c47SLeila Ghaffari break; 376e6225c47SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 377e6225c47SLeila Ghaffari break; 378e6225c47SLeila Ghaffari } 379*f0b01153SJames Wright } 38077841947SLeila Ghaffari return 0; 38177841947SLeila Ghaffari } 38277841947SLeila Ghaffari 38377841947SLeila Ghaffari // ***************************************************************************** 384ea61e9acSJeremy L Thompson // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method 38577841947SLeila Ghaffari // ***************************************************************************** 3862b730f8bSJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 38746603fc5SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 38846603fc5SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 38946603fc5SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 390f3e15844SJames Wright const CeedScalar(*q_data) = in[3]; 39146603fc5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 39246603fc5SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 39329ea4e10SJames Wright CeedScalar *jac_data = out[2]; 39477841947SLeila Ghaffari 395e6225c47SLeila Ghaffari EulerContext context = (EulerContext)ctx; 396932417b3SJed Brown const CeedScalar c_tau = context->c_tau; 39777841947SLeila Ghaffari const CeedScalar gamma = 1.4; 39829ea4e10SJames Wright const CeedScalar zeros[14] = {0.}; 39977841947SLeila Ghaffari 40046603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 40177841947SLeila Ghaffari // Setup 40277841947SLeila Ghaffari // -- Interp in 40377841947SLeila Ghaffari const CeedScalar rho = q[0][i]; 4042b730f8bSJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 40577841947SLeila Ghaffari const CeedScalar E = q[4][i]; 4062b730f8bSJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 4072b730f8bSJeremy L Thompson const CeedScalar dU[3][3] = { 4082b730f8bSJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 4092b730f8bSJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 4102b730f8bSJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 411e6225c47SLeila Ghaffari }; 4122b730f8bSJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 413f3e15844SJames Wright CeedScalar wdetJ, dXdx[3][3]; 414f3e15844SJames Wright QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx); 415e6225c47SLeila Ghaffari // dU/dx 416e6225c47SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 417e6225c47SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 418e6225c47SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 419e6225c47SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 420ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) { 421ba6664aeSJames Wright for (CeedInt k = 0; k < 3; k++) { 422e6225c47SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 423e6225c47SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 424ba6664aeSJames Wright for (CeedInt l = 0; l < 3; l++) { 425e6225c47SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 426e6225c47SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 427e6225c47SLeila Ghaffari } 428e6225c47SLeila Ghaffari } 429e6225c47SLeila Ghaffari } 4302b730f8bSJeremy L Thompson const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic, 431e6225c47SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 43277841947SLeila Ghaffari 43377841947SLeila Ghaffari // The Physics 43477841947SLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 435ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) { 436e6225c47SLeila Ghaffari v[j][i] = 0.; 4372b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 43877841947SLeila Ghaffari } 43977841947SLeila Ghaffari //-----mass matrix 4402b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i]; 44177841947SLeila Ghaffari 44277841947SLeila Ghaffari // -- Density 44377841947SLeila Ghaffari // ---- u rho 4442b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][0][i] -= wdetJ * (rho * u[0] * dXdx[j][0] + rho * u[1] * dXdx[j][1] + rho * u[2] * dXdx[j][2]); 44577841947SLeila Ghaffari // -- Momentum 44677841947SLeila Ghaffari // ---- rho (u x u) + P I3 4472b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 4482b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { 4492b730f8bSJeremy L Thompson dv[k][j + 1][i] -= wdetJ * ((rho * u[j] * u[0] + (j == 0 ? P : 0.)) * dXdx[k][0] + (rho * u[j] * u[1] + (j == 1 ? P : 0.)) * dXdx[k][1] + 450e6225c47SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 4512b730f8bSJeremy L Thompson } 4522b730f8bSJeremy L Thompson } 45377841947SLeila Ghaffari // -- Total Energy Density 45477841947SLeila Ghaffari // ---- (E + P) u 4552b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][4][i] -= wdetJ * (E + P) * (u[0] * dXdx[j][0] + u[1] * dXdx[j][1] + u[2] * dXdx[j][2]); 456e6225c47SLeila Ghaffari 457e6225c47SLeila Ghaffari // -- Stabilization terms 458e6225c47SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 459e6225c47SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 460932417b3SJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 461e6225c47SLeila Ghaffari 462e6225c47SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 463e6225c47SLeila Ghaffari CeedScalar dqdx[5][3]; 464ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) { 465e6225c47SLeila Ghaffari dqdx[0][j] = drhodx[j]; 466e6225c47SLeila Ghaffari dqdx[4][j] = dEdx[j]; 4672b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 468e6225c47SLeila Ghaffari } 469e6225c47SLeila Ghaffari 470e6225c47SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 471e6225c47SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 4722b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 4732b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 4742b730f8bSJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 4752b730f8bSJeremy L Thompson } 4762b730f8bSJeremy L Thompson } 477e6225c47SLeila Ghaffari 478e6225c47SLeila Ghaffari // ---- Strong residual 479e6225c47SLeila Ghaffari CeedScalar strong_res[5]; 4802b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j]; 481e6225c47SLeila Ghaffari 482932417b3SJed Brown // Stabilization 483932417b3SJed Brown // -- Tau elements 484932417b3SJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 485932417b3SJed Brown CeedScalar Tau_x[3] = {0.}; 486932417b3SJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 487e6225c47SLeila Ghaffari 488932417b3SJed Brown // -- Stabilization method: none, SU, or SUPG 48988626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 490e6225c47SLeila Ghaffari switch (context->stabilization) { 491e6225c47SLeila Ghaffari case 0: // Galerkin 492e6225c47SLeila Ghaffari break; 493e6225c47SLeila Ghaffari case 1: // SU 4942b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 4952b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 4962b730f8bSJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 4972b730f8bSJeremy L Thompson } 4982b730f8bSJeremy L Thompson } 499e6225c47SLeila Ghaffari 5002b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 5012b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] += wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]); 5022b730f8bSJeremy L Thompson } 503e6225c47SLeila Ghaffari break; 504e6225c47SLeila Ghaffari case 2: // SUPG 5052b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 5062b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 5072b730f8bSJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l]; 5082b730f8bSJeremy L Thompson } 5092b730f8bSJeremy L Thompson } 510e6225c47SLeila Ghaffari 5112b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 5122b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] += wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]); 5132b730f8bSJeremy L Thompson } 514e6225c47SLeila Ghaffari break; 515e6225c47SLeila Ghaffari } 51629ea4e10SJames Wright StoredValuesPack(Q, i, 0, 14, zeros, jac_data); 517*f0b01153SJames Wright } 51877841947SLeila Ghaffari return 0; 51977841947SLeila Ghaffari } 52077841947SLeila Ghaffari // ***************************************************************************** 521ea61e9acSJeremy L Thompson // This QFunction sets the inflow boundary conditions for the traveling vortex problem. 52277841947SLeila Ghaffari // 523ea61e9acSJeremy L Thompson // Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly. 52477841947SLeila Ghaffari // ***************************************************************************** 5252b730f8bSJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 526f3e15844SJames Wright const CeedScalar(*q_data_sur) = in[2]; 52777841947SLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 528*f0b01153SJames Wright 52977841947SLeila Ghaffari EulerContext context = (EulerContext)ctx; 53077841947SLeila Ghaffari const int euler_test = context->euler_test; 531f3e15844SJames Wright const bool is_implicit = context->implicit; 53277841947SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 53377841947SLeila Ghaffari const CeedScalar cv = 2.5; 53477841947SLeila Ghaffari const CeedScalar R = 1.; 53577841947SLeila Ghaffari CeedScalar T_inlet; 53677841947SLeila Ghaffari CeedScalar P_inlet; 53777841947SLeila Ghaffari 53877841947SLeila Ghaffari // For test cases 1 and 3 the background velocity is zero 5392b730f8bSJeremy L Thompson if (euler_test == 1 || euler_test == 3) { 54077841947SLeila Ghaffari for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.; 5412b730f8bSJeremy L Thompson } 54277841947SLeila Ghaffari 54377841947SLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4 54477841947SLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 54577841947SLeila Ghaffari else T_inlet = P_inlet = 1.; 54677841947SLeila Ghaffari 54746603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 548f3e15844SJames Wright CeedScalar wdetJb, norm[3]; 549f3e15844SJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm); 550f3e15844SJames Wright wdetJb *= is_implicit ? -1. : 1.; 55177841947SLeila Ghaffari 55277841947SLeila Ghaffari // face_normal = Normal vector of the face 5532b730f8bSJeremy L Thompson const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2]; 55477841947SLeila Ghaffari // The Physics 55577841947SLeila Ghaffari // Zero v so all future terms can safely sum into it 556ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 55777841947SLeila Ghaffari 55877841947SLeila Ghaffari // Implementing in/outflow BCs 5592fe7aee7SLeila Ghaffari if (face_normal > 0) { 56077841947SLeila Ghaffari } else { // inflow 56177841947SLeila Ghaffari const CeedScalar rho_inlet = P_inlet / (R * T_inlet); 5622b730f8bSJeremy L Thompson const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.; 56377841947SLeila Ghaffari // incoming total energy 56477841947SLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 56577841947SLeila Ghaffari 56677841947SLeila Ghaffari // The Physics 56777841947SLeila Ghaffari // -- Density 56877841947SLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal; 56977841947SLeila Ghaffari 57077841947SLeila Ghaffari // -- Momentum 5712b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho_inlet * face_normal * mean_velocity[j] + norm[j] * P_inlet); 57277841947SLeila Ghaffari 57377841947SLeila Ghaffari // -- Total Energy Density 57477841947SLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 57577841947SLeila Ghaffari } 576*f0b01153SJames Wright } 57777841947SLeila Ghaffari return 0; 57877841947SLeila Ghaffari } 57977841947SLeila Ghaffari 58077841947SLeila Ghaffari // ***************************************************************************** 581ea61e9acSJeremy L Thompson // This QFunction sets the outflow boundary conditions for the Euler solver. 58255e76554SLeila Ghaffari // 58355e76554SLeila Ghaffari // Outflow BCs: 584ea61e9acSJeremy L Thompson // The validity of the weak form of the governing equations is extended to the outflow. 58555e76554SLeila Ghaffari // ***************************************************************************** 5862b730f8bSJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 58746603fc5SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 588f3e15844SJames Wright const CeedScalar(*q_data_sur) = in[2]; 58955e76554SLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 590*f0b01153SJames Wright 59155e76554SLeila Ghaffari EulerContext context = (EulerContext)ctx; 592f3e15844SJames Wright const bool is_implicit = context->implicit; 59355e76554SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 59455e76554SLeila Ghaffari 59555e76554SLeila Ghaffari const CeedScalar gamma = 1.4; 59655e76554SLeila Ghaffari 59746603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 59855e76554SLeila Ghaffari // Setup 59955e76554SLeila Ghaffari // -- Interp in 60055e76554SLeila Ghaffari const CeedScalar rho = q[0][i]; 6012b730f8bSJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 60255e76554SLeila Ghaffari const CeedScalar E = q[4][i]; 60355e76554SLeila Ghaffari 604f3e15844SJames Wright CeedScalar wdetJb, norm[3]; 605f3e15844SJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm); 606f3e15844SJames Wright wdetJb *= is_implicit ? -1. : 1.; 60755e76554SLeila Ghaffari 60855e76554SLeila Ghaffari // face_normal = Normal vector of the face 6092b730f8bSJeremy L Thompson const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2]; 61055e76554SLeila Ghaffari // The Physics 61155e76554SLeila Ghaffari // Zero v so all future terms can safely sum into it 612ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0; 61355e76554SLeila Ghaffari 61455e76554SLeila Ghaffari // Implementing in/outflow BCs 61555e76554SLeila Ghaffari if (face_normal > 0) { // outflow 61655e76554SLeila Ghaffari const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.; 61755e76554SLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 6182b730f8bSJeremy L Thompson const CeedScalar u_normal = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2]; // Normal velocity 61955e76554SLeila Ghaffari // The Physics 62055e76554SLeila Ghaffari // -- Density 62155e76554SLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal; 62255e76554SLeila Ghaffari 62355e76554SLeila Ghaffari // -- Momentum 6242b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 62555e76554SLeila Ghaffari 62655e76554SLeila Ghaffari // -- Total Energy Density 62755e76554SLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P); 62855e76554SLeila Ghaffari } 629*f0b01153SJames Wright } 63055e76554SLeila Ghaffari return 0; 63155e76554SLeila Ghaffari } 632