1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3a515125bSLeila Ghaffari // 4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5a515125bSLeila Ghaffari // 6727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7a515125bSLeila Ghaffari 8a515125bSLeila Ghaffari /// @file 9a515125bSLeila Ghaffari /// Euler traveling vortex initial condition and operator for Navier-Stokes 10a515125bSLeila Ghaffari /// example using PETSc 11a515125bSLeila Ghaffari 12a515125bSLeila Ghaffari // Model from: 1304e40bb6SJeremy L Thompson // On the Order of Accuracy and Numerical Performance of Two Classes of Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011). 143a8779fbSJames Wright #include <ceed.h> 15d0cce58aSJeremy L Thompson #include <math.h> 162b916ea7SJeremy L Thompson 17704b8bbeSJames Wright #include "utils.h" 18a515125bSLeila Ghaffari 19a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext; 20a515125bSLeila Ghaffari struct EulerContext_ { 21a515125bSLeila Ghaffari CeedScalar center[3]; 22a515125bSLeila Ghaffari CeedScalar curr_time; 23a515125bSLeila Ghaffari CeedScalar vortex_strength; 24d8a22b9eSJed Brown CeedScalar c_tau; 25a515125bSLeila Ghaffari CeedScalar mean_velocity[3]; 26a515125bSLeila Ghaffari bool implicit; 27139613f2SLeila Ghaffari int euler_test; 28139613f2SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 29a515125bSLeila Ghaffari }; 30a515125bSLeila Ghaffari 31a515125bSLeila Ghaffari // ***************************************************************************** 32a515125bSLeila Ghaffari // This function sets the initial conditions 33a515125bSLeila Ghaffari // 34a515125bSLeila Ghaffari // Temperature: 35a515125bSLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 36a515125bSLeila Ghaffari // Density: 37a515125bSLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1)) 38a515125bSLeila Ghaffari // Pressure: 39a515125bSLeila Ghaffari // P = rho * T 40a515125bSLeila Ghaffari // Velocity: 41a515125bSLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 42a515125bSLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 43a515125bSLeila Ghaffari // Velocity/Momentum Density: 44a515125bSLeila Ghaffari // Ui = rho ui 45a515125bSLeila Ghaffari // Total Energy: 46a515125bSLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2 47a515125bSLeila Ghaffari // 48a515125bSLeila Ghaffari // Constants: 49a515125bSLeila Ghaffari // cv , Specific heat, constant volume 50a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 51a515125bSLeila Ghaffari // vortex_strength , Strength of vortex 52a515125bSLeila Ghaffari // center , Location of bubble center 53a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 54a515125bSLeila Ghaffari // 55a515125bSLeila Ghaffari // ***************************************************************************** 56a515125bSLeila Ghaffari 57a515125bSLeila Ghaffari // ***************************************************************************** 5804e40bb6SJeremy L Thompson // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling 5904e40bb6SJeremy L Thompson // vortex 60a515125bSLeila Ghaffari // ***************************************************************************** 612b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) { 62a515125bSLeila Ghaffari // Context 63a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 64a515125bSLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength; 65a515125bSLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain 66a515125bSLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity; 67a515125bSLeila Ghaffari 68a515125bSLeila Ghaffari // Setup 69a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 70a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 71a515125bSLeila Ghaffari const CeedScalar R = 1.; 72a515125bSLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates 73a515125bSLeila Ghaffari // Vortex center 74a515125bSLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time; 75a515125bSLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time; 76a515125bSLeila Ghaffari 77a515125bSLeila Ghaffari const CeedScalar x0 = x - xc; 78a515125bSLeila Ghaffari const CeedScalar y0 = y - yc; 79a515125bSLeila Ghaffari const CeedScalar r = sqrt(x0 * x0 + y0 * y0); 80a515125bSLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI); 812b916ea7SJeremy L Thompson const CeedScalar delta_T = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI); 82a515125bSLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 832b916ea7SJeremy L Thompson const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI); 84a515125bSLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.}; 85a515125bSLeila Ghaffari 86a515125bSLeila Ghaffari // Initial Conditions 87a515125bSLeila Ghaffari switch (context->euler_test) { 88a515125bSLeila Ghaffari case 0: // Traveling vortex 89a515125bSLeila Ghaffari T = 1 + delta_T; 90a515125bSLeila Ghaffari // P = rho * T 91a515125bSLeila Ghaffari // P = S * rho^gamma 92a515125bSLeila Ghaffari // Solve for rho, then substitute for P 93139613f2SLeila Ghaffari rho = pow(T / S_vortex, 1 / (gamma - 1.)); 94a515125bSLeila Ghaffari P = rho * T; 95a515125bSLeila Ghaffari u[0] = mean_velocity[0] - C * y0; 96a515125bSLeila Ghaffari u[1] = mean_velocity[1] + C * x0; 97a515125bSLeila Ghaffari 98a515125bSLeila Ghaffari // Assign exact solution 99a515125bSLeila Ghaffari q[0] = rho; 100a515125bSLeila Ghaffari q[1] = rho * u[0]; 101a515125bSLeila Ghaffari q[2] = rho * u[1]; 102a515125bSLeila Ghaffari q[3] = rho * u[2]; 103a515125bSLeila Ghaffari q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.; 104a515125bSLeila Ghaffari break; 105a515125bSLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant 106a515125bSLeila Ghaffari rho = 1.; 107a515125bSLeila Ghaffari E = 2.; 108a515125bSLeila Ghaffari 109a515125bSLeila Ghaffari // Assign exact solution 110a515125bSLeila Ghaffari q[0] = rho; 111a515125bSLeila Ghaffari q[1] = rho * u[0]; 112a515125bSLeila Ghaffari q[2] = rho * u[1]; 113a515125bSLeila Ghaffari q[3] = rho * u[2]; 114a515125bSLeila Ghaffari q[4] = E; 115a515125bSLeila Ghaffari break; 116a515125bSLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant 117a515125bSLeila Ghaffari rho = 1.; 118a515125bSLeila Ghaffari E = 2.; 119a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 120a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 121a515125bSLeila Ghaffari 122a515125bSLeila Ghaffari // Assign exact solution 123a515125bSLeila Ghaffari q[0] = rho; 124a515125bSLeila Ghaffari q[1] = rho * u[0]; 125a515125bSLeila Ghaffari q[2] = rho * u[1]; 126a515125bSLeila Ghaffari q[3] = rho * u[2]; 127a515125bSLeila Ghaffari q[4] = E; 128a515125bSLeila Ghaffari break; 12904e40bb6SJeremy 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 13004e40bb6SJeremy L Thompson // bubble won't diffuse 131a515125bSLeila Ghaffari // (for Euler, where there is no thermal conductivity) 132a515125bSLeila Ghaffari P = 1.; 133a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r); 134a515125bSLeila Ghaffari rho = P / (R * T); 135a515125bSLeila Ghaffari 136a515125bSLeila Ghaffari // Assign exact solution 137a515125bSLeila Ghaffari q[0] = rho; 138a515125bSLeila Ghaffari q[1] = rho * u[0]; 139a515125bSLeila Ghaffari q[2] = rho * u[1]; 140a515125bSLeila Ghaffari q[3] = rho * u[2]; 141a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 142a515125bSLeila Ghaffari break; 14304e40bb6SJeremy L Thompson case 4: // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant), 14404e40bb6SJeremy L Thompson // It should be transported across the domain, but velocity stays constant 145a515125bSLeila Ghaffari P = 1.; 146a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r); 147a515125bSLeila Ghaffari rho = P / (R * T); 148a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 149a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 150a515125bSLeila Ghaffari 151a515125bSLeila Ghaffari // Assign exact solution 152a515125bSLeila Ghaffari q[0] = rho; 153a515125bSLeila Ghaffari q[1] = rho * u[0]; 154a515125bSLeila Ghaffari q[2] = rho * u[1]; 155a515125bSLeila Ghaffari q[3] = rho * u[2]; 156a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 157a515125bSLeila Ghaffari break; 1580df2634dSLeila Ghaffari case 5: // non-smooth thermal bubble - cylinder 1590df2634dSLeila Ghaffari P = 1.; 1600df2634dSLeila Ghaffari T = 1. - (r < 1. ? S_bubble : 0.); 1610df2634dSLeila Ghaffari rho = P / (R * T); 1620df2634dSLeila Ghaffari u[0] = mean_velocity[0]; 1630df2634dSLeila Ghaffari u[1] = mean_velocity[1]; 1640df2634dSLeila Ghaffari 1650df2634dSLeila Ghaffari // Assign exact solution 1660df2634dSLeila Ghaffari q[0] = rho; 1670df2634dSLeila Ghaffari q[1] = rho * u[0]; 1680df2634dSLeila Ghaffari q[2] = rho * u[1]; 1690df2634dSLeila Ghaffari q[3] = rho * u[2]; 1700df2634dSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 1710df2634dSLeila Ghaffari break; 172a515125bSLeila Ghaffari } 173a515125bSLeila Ghaffari return 0; 174a515125bSLeila Ghaffari } 175a515125bSLeila Ghaffari 176a515125bSLeila Ghaffari // ***************************************************************************** 177139613f2SLeila Ghaffari // Helper function for computing flux Jacobian 178139613f2SLeila Ghaffari // ***************************************************************************** 1792b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 180139613f2SLeila Ghaffari const CeedScalar gamma) { 181139613f2SLeila Ghaffari CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2]; // Velocity square 182139613f2SLeila Ghaffari for (CeedInt i = 0; i < 3; i++) { // Jacobian matrices for 3 directions 183139613f2SLeila Ghaffari for (CeedInt j = 0; j < 3; j++) { // Rows of each Jacobian matrix 184139613f2SLeila Ghaffari dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j]; 185139613f2SLeila Ghaffari for (CeedInt k = 0; k < 3; k++) { // Columns of each Jacobian matrix 186139613f2SLeila Ghaffari dF[i][0][k + 1] = ((i == k) ? 1. : 0.); 1872b916ea7SJeremy 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.); 1882b916ea7SJeremy L Thompson dF[i][4][k + 1] = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k]; 189139613f2SLeila Ghaffari } 190139613f2SLeila Ghaffari dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.); 191139613f2SLeila Ghaffari } 192139613f2SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho); 193139613f2SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 194139613f2SLeila Ghaffari } 195139613f2SLeila Ghaffari } 196139613f2SLeila Ghaffari 197139613f2SLeila Ghaffari // ***************************************************************************** 198d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant) 199d8a22b9eSJed Brown // Model from: 200d8a22b9eSJed Brown // Stabilized Methods for Compressible Flows, Hughes et al 2010 201d8a22b9eSJed Brown // 202d8a22b9eSJed Brown // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 203d8a22b9eSJed Brown // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 204d8a22b9eSJed Brown // 205d8a22b9eSJed Brown // Where 206d8a22b9eSJed Brown // c_tau = stabilization constant (0.5 is reported as "optimal") 207d8a22b9eSJed Brown // h[i] = 2 length(dxdX[i]) 208d8a22b9eSJed Brown // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 209d8a22b9eSJed Brown // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 21004e40bb6SJeremy L Thompson // rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i 211d8a22b9eSJed Brown // ***************************************************************************** 2122b916ea7SJeremy 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, 2132b916ea7SJeremy L Thompson const CeedScalar c_tau) { 214493642f1SJames Wright for (CeedInt i = 0; i < 3; i++) { 215d8a22b9eSJed Brown // length of element in direction i 2162b916ea7SJeremy L Thompson CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]); 217d8a22b9eSJed Brown // fastest wave in direction i 218d8a22b9eSJed Brown CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 219d8a22b9eSJed Brown Tau_x[i] = c_tau * h / fastest_wave; 220d8a22b9eSJed Brown } 221d8a22b9eSJed Brown } 222d8a22b9eSJed Brown 223d8a22b9eSJed Brown // ***************************************************************************** 224a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex 225a515125bSLeila Ghaffari // ***************************************************************************** 2262b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 227a515125bSLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 228a515125bSLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 229*b193fadcSJames Wright 230a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 231a515125bSLeila Ghaffari 2323d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 233a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 234139613f2SLeila Ghaffari CeedScalar q[5] = {0.}; 235a515125bSLeila Ghaffari 236a515125bSLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx); 2372b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 238*b193fadcSJames Wright } 239a515125bSLeila Ghaffari return 0; 240a515125bSLeila Ghaffari } 241a515125bSLeila Ghaffari 242a515125bSLeila Ghaffari // ***************************************************************************** 24304e40bb6SJeremy L Thompson // This QFunction implements the following formulation of Euler equations with explicit time stepping method 244a515125bSLeila Ghaffari // 24504e40bb6SJeremy L Thompson // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density. 246a515125bSLeila Ghaffari // 247a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 248a515125bSLeila Ghaffari // rho - Mass Density 249a515125bSLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 250a515125bSLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 251a515125bSLeila Ghaffari // 252a515125bSLeila Ghaffari // Euler Equations: 253a515125bSLeila Ghaffari // drho/dt + div( U ) = 0 254a515125bSLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0 255a515125bSLeila Ghaffari // dE/dt + div( (E + P) u ) = 0 256a515125bSLeila Ghaffari // 257a515125bSLeila Ghaffari // Equation of State: 258a515125bSLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2) 259a515125bSLeila Ghaffari // 260a515125bSLeila Ghaffari // Constants: 261a515125bSLeila Ghaffari // cv , Specific heat, constant volume 262a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 263a515125bSLeila Ghaffari // g , Gravity 264a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 265a515125bSLeila Ghaffari // ***************************************************************************** 2662b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 2673d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 2683d65b166SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 269ade49511SJames Wright const CeedScalar(*q_data) = in[2]; 2703d65b166SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 2713d65b166SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 272a515125bSLeila Ghaffari 273139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 274d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 275a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 276a515125bSLeila Ghaffari 2773d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 278a515125bSLeila Ghaffari // Setup 279a515125bSLeila Ghaffari // -- Interp in 280a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 2812b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 282a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 2832b916ea7SJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 2842b916ea7SJeremy L Thompson const CeedScalar dU[3][3] = { 2852b916ea7SJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 2862b916ea7SJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 2872b916ea7SJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 288139613f2SLeila Ghaffari }; 2892b916ea7SJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 290ade49511SJames Wright CeedScalar wdetJ, dXdx[3][3]; 291ade49511SJames Wright QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx); 292139613f2SLeila Ghaffari // dU/dx 293139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 294139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 295139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 296139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 297493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 298493642f1SJames Wright for (CeedInt k = 0; k < 3; k++) { 299139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 300139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 301493642f1SJames Wright for (CeedInt l = 0; l < 3; l++) { 302139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 303139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 304139613f2SLeila Ghaffari } 305139613f2SLeila Ghaffari } 306139613f2SLeila Ghaffari } 307139613f2SLeila Ghaffari // Pressure 3082b916ea7SJeremy 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, 309139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 310a515125bSLeila Ghaffari 311a515125bSLeila Ghaffari // The Physics 312a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 313493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) { 314139613f2SLeila Ghaffari v[j][i] = 0.; 3152b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 316a515125bSLeila Ghaffari } 317a515125bSLeila Ghaffari 318a515125bSLeila Ghaffari // -- Density 319a515125bSLeila Ghaffari // ---- u rho 3202b916ea7SJeremy 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]); 321a515125bSLeila Ghaffari // -- Momentum 322a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 3232b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3242b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { 3252b916ea7SJeremy 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] + 326139613f2SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 3272b916ea7SJeremy L Thompson } 3282b916ea7SJeremy L Thompson } 329a515125bSLeila Ghaffari // -- Total Energy Density 330a515125bSLeila Ghaffari // ---- (E + P) u 3312b916ea7SJeremy 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]); 332139613f2SLeila Ghaffari 333139613f2SLeila Ghaffari // --Stabilization terms 334139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 335139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 336d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 337139613f2SLeila Ghaffari 338139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 339139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 340493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 341139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 342139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 3432b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 344139613f2SLeila Ghaffari } 345139613f2SLeila Ghaffari 346139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 347139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 3482b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3492b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 3502b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 3512b916ea7SJeremy L Thompson } 3522b916ea7SJeremy L Thompson } 353139613f2SLeila Ghaffari 354d8a22b9eSJed Brown // Stabilization 355d8a22b9eSJed Brown // -- Tau elements 356d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 357d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 358d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 359139613f2SLeila Ghaffari 360d8a22b9eSJed Brown // -- Stabilization method: none or SU 361bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 362139613f2SLeila Ghaffari switch (context->stabilization) { 363139613f2SLeila Ghaffari case 0: // Galerkin 364139613f2SLeila Ghaffari break; 365139613f2SLeila Ghaffari case 1: // SU 3662b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3672b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 3682b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 3692b916ea7SJeremy L Thompson } 3702b916ea7SJeremy L Thompson } 371139613f2SLeila Ghaffari 3722b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 3732b916ea7SJeremy 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]); 3742b916ea7SJeremy L Thompson } 375139613f2SLeila Ghaffari break; 376139613f2SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 377139613f2SLeila Ghaffari break; 378139613f2SLeila Ghaffari } 379*b193fadcSJames Wright } 380a515125bSLeila Ghaffari return 0; 381a515125bSLeila Ghaffari } 382a515125bSLeila Ghaffari 383a515125bSLeila Ghaffari // ***************************************************************************** 38404e40bb6SJeremy L Thompson // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method 385a515125bSLeila Ghaffari // ***************************************************************************** 3862b916ea7SJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 3873d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 3883d65b166SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 3893d65b166SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 390ade49511SJames Wright const CeedScalar(*q_data) = in[3]; 3913d65b166SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 3923d65b166SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 39380f5d3cbSJames Wright CeedScalar *jac_data = out[2]; 394a515125bSLeila Ghaffari 395139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 396d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 397a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 39880f5d3cbSJames Wright const CeedScalar zeros[14] = {0.}; 399a515125bSLeila Ghaffari 4003d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 401a515125bSLeila Ghaffari // Setup 402a515125bSLeila Ghaffari // -- Interp in 403a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 4042b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 405a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 4062b916ea7SJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 4072b916ea7SJeremy L Thompson const CeedScalar dU[3][3] = { 4082b916ea7SJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 4092b916ea7SJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 4102b916ea7SJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 411139613f2SLeila Ghaffari }; 4122b916ea7SJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 413ade49511SJames Wright CeedScalar wdetJ, dXdx[3][3]; 414ade49511SJames Wright QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx); 415139613f2SLeila Ghaffari // dU/dx 416139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 417139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 418139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 419139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 420493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 421493642f1SJames Wright for (CeedInt k = 0; k < 3; k++) { 422139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 423139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 424493642f1SJames Wright for (CeedInt l = 0; l < 3; l++) { 425139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 426139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 427139613f2SLeila Ghaffari } 428139613f2SLeila Ghaffari } 429139613f2SLeila Ghaffari } 4302b916ea7SJeremy 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, 431139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 432a515125bSLeila Ghaffari 433a515125bSLeila Ghaffari // The Physics 434a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 435493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) { 436139613f2SLeila Ghaffari v[j][i] = 0.; 4372b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 438a515125bSLeila Ghaffari } 439a515125bSLeila Ghaffari //-----mass matrix 4402b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i]; 441a515125bSLeila Ghaffari 442a515125bSLeila Ghaffari // -- Density 443a515125bSLeila Ghaffari // ---- u rho 4442b916ea7SJeremy 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]); 445a515125bSLeila Ghaffari // -- Momentum 446a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 4472b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 4482b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { 4492b916ea7SJeremy 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] + 450139613f2SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 4512b916ea7SJeremy L Thompson } 4522b916ea7SJeremy L Thompson } 453a515125bSLeila Ghaffari // -- Total Energy Density 454a515125bSLeila Ghaffari // ---- (E + P) u 4552b916ea7SJeremy 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]); 456139613f2SLeila Ghaffari 457139613f2SLeila Ghaffari // -- Stabilization terms 458139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 459139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 460d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 461139613f2SLeila Ghaffari 462139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 463139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 464493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 465139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 466139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 4672b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 468139613f2SLeila Ghaffari } 469139613f2SLeila Ghaffari 470139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 471139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 4722b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 4732b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 4742b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 4752b916ea7SJeremy L Thompson } 4762b916ea7SJeremy L Thompson } 477139613f2SLeila Ghaffari 478139613f2SLeila Ghaffari // ---- Strong residual 479139613f2SLeila Ghaffari CeedScalar strong_res[5]; 4802b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j]; 481139613f2SLeila Ghaffari 482d8a22b9eSJed Brown // Stabilization 483d8a22b9eSJed Brown // -- Tau elements 484d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 485d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 486d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 487139613f2SLeila Ghaffari 488d8a22b9eSJed Brown // -- Stabilization method: none, SU, or SUPG 489bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 490139613f2SLeila Ghaffari switch (context->stabilization) { 491139613f2SLeila Ghaffari case 0: // Galerkin 492139613f2SLeila Ghaffari break; 493139613f2SLeila Ghaffari case 1: // SU 4942b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 4952b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 4962b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 4972b916ea7SJeremy L Thompson } 4982b916ea7SJeremy L Thompson } 499139613f2SLeila Ghaffari 5002b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 5012b916ea7SJeremy 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]); 5022b916ea7SJeremy L Thompson } 503139613f2SLeila Ghaffari break; 504139613f2SLeila Ghaffari case 2: // SUPG 5052b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 5062b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 5072b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l]; 5082b916ea7SJeremy L Thompson } 5092b916ea7SJeremy L Thompson } 510139613f2SLeila Ghaffari 5112b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 5122b916ea7SJeremy 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]); 5132b916ea7SJeremy L Thompson } 514139613f2SLeila Ghaffari break; 515139613f2SLeila Ghaffari } 51680f5d3cbSJames Wright StoredValuesPack(Q, i, 0, 14, zeros, jac_data); 517*b193fadcSJames Wright } 518a515125bSLeila Ghaffari return 0; 519a515125bSLeila Ghaffari } 520a515125bSLeila Ghaffari // ***************************************************************************** 52104e40bb6SJeremy L Thompson // This QFunction sets the inflow boundary conditions for the traveling vortex problem. 522a515125bSLeila Ghaffari // 52304e40bb6SJeremy L Thompson // Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly. 524a515125bSLeila Ghaffari // ***************************************************************************** 5252b916ea7SJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 526ade49511SJames Wright const CeedScalar(*q_data_sur) = in[2]; 527a515125bSLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 528*b193fadcSJames Wright 529a515125bSLeila Ghaffari EulerContext context = (EulerContext)ctx; 530a515125bSLeila Ghaffari const int euler_test = context->euler_test; 531ade49511SJames Wright const bool is_implicit = context->implicit; 532a515125bSLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 533a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 534a515125bSLeila Ghaffari const CeedScalar R = 1.; 535a515125bSLeila Ghaffari CeedScalar T_inlet; 536a515125bSLeila Ghaffari CeedScalar P_inlet; 537a515125bSLeila Ghaffari 538a515125bSLeila Ghaffari // For test cases 1 and 3 the background velocity is zero 5392b916ea7SJeremy L Thompson if (euler_test == 1 || euler_test == 3) { 540a515125bSLeila Ghaffari for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.; 5412b916ea7SJeremy L Thompson } 542a515125bSLeila Ghaffari 543a515125bSLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4 544a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 545a515125bSLeila Ghaffari else T_inlet = P_inlet = 1.; 546a515125bSLeila Ghaffari 5473d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 548ade49511SJames Wright CeedScalar wdetJb, norm[3]; 549ade49511SJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm); 550ade49511SJames Wright wdetJb *= is_implicit ? -1. : 1.; 551a515125bSLeila Ghaffari 552a515125bSLeila Ghaffari // face_normal = Normal vector of the face 5532b916ea7SJeremy L Thompson const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2]; 554a515125bSLeila Ghaffari // The Physics 555a515125bSLeila Ghaffari // Zero v so all future terms can safely sum into it 556493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 557a515125bSLeila Ghaffari 558a515125bSLeila Ghaffari // Implementing in/outflow BCs 559002797a3SLeila Ghaffari if (face_normal > 0) { 560a515125bSLeila Ghaffari } else { // inflow 561a515125bSLeila Ghaffari const CeedScalar rho_inlet = P_inlet / (R * T_inlet); 5622b916ea7SJeremy L Thompson const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.; 563a515125bSLeila Ghaffari // incoming total energy 564a515125bSLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 565a515125bSLeila Ghaffari 566a515125bSLeila Ghaffari // The Physics 567a515125bSLeila Ghaffari // -- Density 568a515125bSLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal; 569a515125bSLeila Ghaffari 570a515125bSLeila Ghaffari // -- Momentum 5712b916ea7SJeremy 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); 572a515125bSLeila Ghaffari 573a515125bSLeila Ghaffari // -- Total Energy Density 574a515125bSLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 575a515125bSLeila Ghaffari } 576*b193fadcSJames Wright } 577a515125bSLeila Ghaffari return 0; 578a515125bSLeila Ghaffari } 579a515125bSLeila Ghaffari 580a515125bSLeila Ghaffari // ***************************************************************************** 58104e40bb6SJeremy L Thompson // This QFunction sets the outflow boundary conditions for the Euler solver. 58268ef3d20SLeila Ghaffari // 58368ef3d20SLeila Ghaffari // Outflow BCs: 58404e40bb6SJeremy L Thompson // The validity of the weak form of the governing equations is extended to the outflow. 58568ef3d20SLeila Ghaffari // ***************************************************************************** 5862b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 5873d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 588ade49511SJames Wright const CeedScalar(*q_data_sur) = in[2]; 58968ef3d20SLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 590*b193fadcSJames Wright 59168ef3d20SLeila Ghaffari EulerContext context = (EulerContext)ctx; 592ade49511SJames Wright const bool is_implicit = context->implicit; 59368ef3d20SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 59468ef3d20SLeila Ghaffari 59568ef3d20SLeila Ghaffari const CeedScalar gamma = 1.4; 59668ef3d20SLeila Ghaffari 5973d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 59868ef3d20SLeila Ghaffari // Setup 59968ef3d20SLeila Ghaffari // -- Interp in 60068ef3d20SLeila Ghaffari const CeedScalar rho = q[0][i]; 6012b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 60268ef3d20SLeila Ghaffari const CeedScalar E = q[4][i]; 60368ef3d20SLeila Ghaffari 604ade49511SJames Wright CeedScalar wdetJb, norm[3]; 605ade49511SJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm); 606ade49511SJames Wright wdetJb *= is_implicit ? -1. : 1.; 60768ef3d20SLeila Ghaffari 60868ef3d20SLeila Ghaffari // face_normal = Normal vector of the face 6092b916ea7SJeremy L Thompson const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2]; 61068ef3d20SLeila Ghaffari // The Physics 61168ef3d20SLeila Ghaffari // Zero v so all future terms can safely sum into it 612493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0; 61368ef3d20SLeila Ghaffari 61468ef3d20SLeila Ghaffari // Implementing in/outflow BCs 61568ef3d20SLeila Ghaffari if (face_normal > 0) { // outflow 61668ef3d20SLeila Ghaffari const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.; 61768ef3d20SLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 6182b916ea7SJeremy L Thompson const CeedScalar u_normal = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2]; // Normal velocity 61968ef3d20SLeila Ghaffari // The Physics 62068ef3d20SLeila Ghaffari // -- Density 62168ef3d20SLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal; 62268ef3d20SLeila Ghaffari 62368ef3d20SLeila Ghaffari // -- Momentum 6242b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 62568ef3d20SLeila Ghaffari 62668ef3d20SLeila Ghaffari // -- Total Energy Density 62768ef3d20SLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P); 62868ef3d20SLeila Ghaffari } 629*b193fadcSJames Wright } 63068ef3d20SLeila Ghaffari return 0; 63168ef3d20SLeila Ghaffari } 632