1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3a515125bSLeila Ghaffari 4a515125bSLeila Ghaffari /// @file 5a515125bSLeila Ghaffari /// Euler traveling vortex initial condition and operator for Navier-Stokes 6a515125bSLeila Ghaffari /// example using PETSc 7a515125bSLeila Ghaffari 8a515125bSLeila Ghaffari // Model from: 904e40bb6SJeremy L Thompson // On the Order of Accuracy and Numerical Performance of Two Classes of Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011). 103a8779fbSJames Wright #include <ceed.h> 11d0cce58aSJeremy L Thompson #include <math.h> 122b916ea7SJeremy L Thompson 13704b8bbeSJames Wright #include "utils.h" 14a515125bSLeila Ghaffari 15a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext; 16a515125bSLeila Ghaffari struct EulerContext_ { 17a515125bSLeila Ghaffari CeedScalar center[3]; 18a515125bSLeila Ghaffari CeedScalar curr_time; 19a515125bSLeila Ghaffari CeedScalar vortex_strength; 20d8a22b9eSJed Brown CeedScalar c_tau; 21a515125bSLeila Ghaffari CeedScalar mean_velocity[3]; 22a515125bSLeila Ghaffari bool implicit; 23139613f2SLeila Ghaffari int euler_test; 24139613f2SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 25a515125bSLeila Ghaffari }; 26a515125bSLeila Ghaffari 27a515125bSLeila Ghaffari // ***************************************************************************** 28a515125bSLeila Ghaffari // This function sets the initial conditions 29a515125bSLeila Ghaffari // 30a515125bSLeila Ghaffari // Temperature: 31a515125bSLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 32a515125bSLeila Ghaffari // Density: 33a515125bSLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1)) 34a515125bSLeila Ghaffari // Pressure: 35a515125bSLeila Ghaffari // P = rho * T 36a515125bSLeila Ghaffari // Velocity: 37a515125bSLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 38a515125bSLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 39a515125bSLeila Ghaffari // Velocity/Momentum Density: 40a515125bSLeila Ghaffari // Ui = rho ui 41a515125bSLeila Ghaffari // Total Energy: 42a515125bSLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2 43a515125bSLeila Ghaffari // 44a515125bSLeila Ghaffari // Constants: 45a515125bSLeila Ghaffari // cv , Specific heat, constant volume 46a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 47a515125bSLeila Ghaffari // vortex_strength , Strength of vortex 48a515125bSLeila Ghaffari // center , Location of bubble center 49a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 50a515125bSLeila Ghaffari // 51a515125bSLeila Ghaffari // ***************************************************************************** 52a515125bSLeila Ghaffari 53a515125bSLeila Ghaffari // ***************************************************************************** 5404e40bb6SJeremy L Thompson // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling 5504e40bb6SJeremy L Thompson // vortex 56a515125bSLeila Ghaffari // ***************************************************************************** 572b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) { 58a515125bSLeila Ghaffari // Context 59a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 60a515125bSLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength; 61a515125bSLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain 62a515125bSLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity; 63a515125bSLeila Ghaffari 64a515125bSLeila Ghaffari // Setup 65a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 66a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 67a515125bSLeila Ghaffari const CeedScalar R = 1.; 68a515125bSLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates 69a515125bSLeila Ghaffari // Vortex center 70a515125bSLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time; 71a515125bSLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time; 72a515125bSLeila Ghaffari 73a515125bSLeila Ghaffari const CeedScalar x0 = x - xc; 74a515125bSLeila Ghaffari const CeedScalar y0 = y - yc; 75a515125bSLeila Ghaffari const CeedScalar r = sqrt(x0 * x0 + y0 * y0); 76a515125bSLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI); 772b916ea7SJeremy L Thompson const CeedScalar delta_T = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI); 78a515125bSLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 792b916ea7SJeremy L Thompson const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI); 80a515125bSLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.}; 81a515125bSLeila Ghaffari 82a515125bSLeila Ghaffari // Initial Conditions 83a515125bSLeila Ghaffari switch (context->euler_test) { 84a515125bSLeila Ghaffari case 0: // Traveling vortex 85a515125bSLeila Ghaffari T = 1 + delta_T; 86a515125bSLeila Ghaffari // P = rho * T 87a515125bSLeila Ghaffari // P = S * rho^gamma 88a515125bSLeila Ghaffari // Solve for rho, then substitute for P 89139613f2SLeila Ghaffari rho = pow(T / S_vortex, 1 / (gamma - 1.)); 90a515125bSLeila Ghaffari P = rho * T; 91a515125bSLeila Ghaffari u[0] = mean_velocity[0] - C * y0; 92a515125bSLeila Ghaffari u[1] = mean_velocity[1] + C * x0; 93a515125bSLeila Ghaffari 94a515125bSLeila Ghaffari // Assign exact solution 95a515125bSLeila Ghaffari q[0] = rho; 96a515125bSLeila Ghaffari q[1] = rho * u[0]; 97a515125bSLeila Ghaffari q[2] = rho * u[1]; 98a515125bSLeila Ghaffari q[3] = rho * u[2]; 99a515125bSLeila Ghaffari q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.; 100a515125bSLeila Ghaffari break; 101a515125bSLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant 102a515125bSLeila Ghaffari rho = 1.; 103a515125bSLeila Ghaffari E = 2.; 104a515125bSLeila Ghaffari 105a515125bSLeila Ghaffari // Assign exact solution 106a515125bSLeila Ghaffari q[0] = rho; 107a515125bSLeila Ghaffari q[1] = rho * u[0]; 108a515125bSLeila Ghaffari q[2] = rho * u[1]; 109a515125bSLeila Ghaffari q[3] = rho * u[2]; 110a515125bSLeila Ghaffari q[4] = E; 111a515125bSLeila Ghaffari break; 112a515125bSLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant 113a515125bSLeila Ghaffari rho = 1.; 114a515125bSLeila Ghaffari E = 2.; 115a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 116a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 117a515125bSLeila Ghaffari 118a515125bSLeila Ghaffari // Assign exact solution 119a515125bSLeila Ghaffari q[0] = rho; 120a515125bSLeila Ghaffari q[1] = rho * u[0]; 121a515125bSLeila Ghaffari q[2] = rho * u[1]; 122a515125bSLeila Ghaffari q[3] = rho * u[2]; 123a515125bSLeila Ghaffari q[4] = E; 124a515125bSLeila Ghaffari break; 12504e40bb6SJeremy 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 12604e40bb6SJeremy L Thompson // bubble won't diffuse 127a515125bSLeila Ghaffari // (for Euler, where there is no thermal conductivity) 128a515125bSLeila Ghaffari P = 1.; 129a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r); 130a515125bSLeila Ghaffari rho = P / (R * T); 131a515125bSLeila Ghaffari 132a515125bSLeila Ghaffari // Assign exact solution 133a515125bSLeila Ghaffari q[0] = rho; 134a515125bSLeila Ghaffari q[1] = rho * u[0]; 135a515125bSLeila Ghaffari q[2] = rho * u[1]; 136a515125bSLeila Ghaffari q[3] = rho * u[2]; 137a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 138a515125bSLeila Ghaffari break; 13904e40bb6SJeremy L Thompson case 4: // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant), 14004e40bb6SJeremy L Thompson // It should be transported across the domain, but velocity stays constant 141a515125bSLeila Ghaffari P = 1.; 142a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r); 143a515125bSLeila Ghaffari rho = P / (R * T); 144a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 145a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 146a515125bSLeila Ghaffari 147a515125bSLeila Ghaffari // Assign exact solution 148a515125bSLeila Ghaffari q[0] = rho; 149a515125bSLeila Ghaffari q[1] = rho * u[0]; 150a515125bSLeila Ghaffari q[2] = rho * u[1]; 151a515125bSLeila Ghaffari q[3] = rho * u[2]; 152a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 153a515125bSLeila Ghaffari break; 1540df2634dSLeila Ghaffari case 5: // non-smooth thermal bubble - cylinder 1550df2634dSLeila Ghaffari P = 1.; 1560df2634dSLeila Ghaffari T = 1. - (r < 1. ? S_bubble : 0.); 1570df2634dSLeila Ghaffari rho = P / (R * T); 1580df2634dSLeila Ghaffari u[0] = mean_velocity[0]; 1590df2634dSLeila Ghaffari u[1] = mean_velocity[1]; 1600df2634dSLeila Ghaffari 1610df2634dSLeila Ghaffari // Assign exact solution 1620df2634dSLeila Ghaffari q[0] = rho; 1630df2634dSLeila Ghaffari q[1] = rho * u[0]; 1640df2634dSLeila Ghaffari q[2] = rho * u[1]; 1650df2634dSLeila Ghaffari q[3] = rho * u[2]; 1660df2634dSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 1670df2634dSLeila Ghaffari break; 168a515125bSLeila Ghaffari } 169a515125bSLeila Ghaffari return 0; 170a515125bSLeila Ghaffari } 171a515125bSLeila Ghaffari 172a515125bSLeila Ghaffari // ***************************************************************************** 173139613f2SLeila Ghaffari // Helper function for computing flux Jacobian 174139613f2SLeila Ghaffari // ***************************************************************************** 1752b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 176139613f2SLeila Ghaffari const CeedScalar gamma) { 177139613f2SLeila Ghaffari CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2]; // Velocity square 178139613f2SLeila Ghaffari for (CeedInt i = 0; i < 3; i++) { // Jacobian matrices for 3 directions 179139613f2SLeila Ghaffari for (CeedInt j = 0; j < 3; j++) { // Rows of each Jacobian matrix 180139613f2SLeila Ghaffari dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j]; 181139613f2SLeila Ghaffari for (CeedInt k = 0; k < 3; k++) { // Columns of each Jacobian matrix 182139613f2SLeila Ghaffari dF[i][0][k + 1] = ((i == k) ? 1. : 0.); 1832b916ea7SJeremy 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.); 1842b916ea7SJeremy L Thompson dF[i][4][k + 1] = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k]; 185139613f2SLeila Ghaffari } 186139613f2SLeila Ghaffari dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.); 187139613f2SLeila Ghaffari } 188139613f2SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho); 189139613f2SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 190139613f2SLeila Ghaffari } 191139613f2SLeila Ghaffari } 192139613f2SLeila Ghaffari 193139613f2SLeila Ghaffari // ***************************************************************************** 194d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant) 195d8a22b9eSJed Brown // Model from: 196d8a22b9eSJed Brown // Stabilized Methods for Compressible Flows, Hughes et al 2010 197d8a22b9eSJed Brown // 198d8a22b9eSJed Brown // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 199d8a22b9eSJed Brown // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 200d8a22b9eSJed Brown // 201d8a22b9eSJed Brown // Where 202d8a22b9eSJed Brown // c_tau = stabilization constant (0.5 is reported as "optimal") 203d8a22b9eSJed Brown // h[i] = 2 length(dxdX[i]) 204d8a22b9eSJed Brown // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 205d8a22b9eSJed Brown // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 20604e40bb6SJeremy L Thompson // rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i 207d8a22b9eSJed Brown // ***************************************************************************** 2082b916ea7SJeremy 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, 2092b916ea7SJeremy L Thompson const CeedScalar c_tau) { 210493642f1SJames Wright for (CeedInt i = 0; i < 3; i++) { 211d8a22b9eSJed Brown // length of element in direction i 2122b916ea7SJeremy L Thompson CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]); 213d8a22b9eSJed Brown // fastest wave in direction i 214d8a22b9eSJed Brown CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 215d8a22b9eSJed Brown Tau_x[i] = c_tau * h / fastest_wave; 216d8a22b9eSJed Brown } 217d8a22b9eSJed Brown } 218d8a22b9eSJed Brown 219d8a22b9eSJed Brown // ***************************************************************************** 220a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex 221a515125bSLeila Ghaffari // ***************************************************************************** 2222b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 223a515125bSLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 224a515125bSLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 225b193fadcSJames Wright 226a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 227a515125bSLeila Ghaffari 2283d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 229a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 230139613f2SLeila Ghaffari CeedScalar q[5] = {0.}; 231a515125bSLeila Ghaffari 232a515125bSLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx); 2332b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 234b193fadcSJames Wright } 235a515125bSLeila Ghaffari return 0; 236a515125bSLeila Ghaffari } 237a515125bSLeila Ghaffari 238a515125bSLeila Ghaffari // ***************************************************************************** 23904e40bb6SJeremy L Thompson // This QFunction implements the following formulation of Euler equations with explicit time stepping method 240a515125bSLeila Ghaffari // 24104e40bb6SJeremy L Thompson // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density. 242a515125bSLeila Ghaffari // 243a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 244a515125bSLeila Ghaffari // rho - Mass Density 245a515125bSLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 246a515125bSLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 247a515125bSLeila Ghaffari // 248a515125bSLeila Ghaffari // Euler Equations: 249a515125bSLeila Ghaffari // drho/dt + div( U ) = 0 250a515125bSLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0 251a515125bSLeila Ghaffari // dE/dt + div( (E + P) u ) = 0 252a515125bSLeila Ghaffari // 253a515125bSLeila Ghaffari // Equation of State: 254a515125bSLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2) 255a515125bSLeila Ghaffari // 256a515125bSLeila Ghaffari // Constants: 257a515125bSLeila Ghaffari // cv , Specific heat, constant volume 258a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 259a515125bSLeila Ghaffari // g , Gravity 260a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 261a515125bSLeila Ghaffari // ***************************************************************************** 2622b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 2633d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 2643d65b166SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 265ade49511SJames Wright const CeedScalar(*q_data) = in[2]; 2663d65b166SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 2673d65b166SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 268a515125bSLeila Ghaffari 269139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 270d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 271a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 272a515125bSLeila Ghaffari 2733d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 274a515125bSLeila Ghaffari // Setup 275a515125bSLeila Ghaffari // -- Interp in 276a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 2772b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 278a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 2792b916ea7SJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 2802b916ea7SJeremy L Thompson const CeedScalar dU[3][3] = { 2812b916ea7SJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 2822b916ea7SJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 2832b916ea7SJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 284139613f2SLeila Ghaffari }; 2852b916ea7SJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 286ade49511SJames Wright CeedScalar wdetJ, dXdx[3][3]; 287ade49511SJames Wright QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx); 288139613f2SLeila Ghaffari // dU/dx 289139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 290139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 291139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 292139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 293493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 294493642f1SJames Wright for (CeedInt k = 0; k < 3; k++) { 295139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 296139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 297493642f1SJames Wright for (CeedInt l = 0; l < 3; l++) { 298139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 299139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 300139613f2SLeila Ghaffari } 301139613f2SLeila Ghaffari } 302139613f2SLeila Ghaffari } 303139613f2SLeila Ghaffari // Pressure 3042b916ea7SJeremy 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, 305139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 306a515125bSLeila Ghaffari 307a515125bSLeila Ghaffari // The Physics 308a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 309493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) { 310139613f2SLeila Ghaffari v[j][i] = 0.; 3112b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 312a515125bSLeila Ghaffari } 313a515125bSLeila Ghaffari 314a515125bSLeila Ghaffari // -- Density 315a515125bSLeila Ghaffari // ---- u rho 3162b916ea7SJeremy 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]); 317a515125bSLeila Ghaffari // -- Momentum 318a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 3192b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3202b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { 3212b916ea7SJeremy 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] + 322139613f2SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 3232b916ea7SJeremy L Thompson } 3242b916ea7SJeremy L Thompson } 325a515125bSLeila Ghaffari // -- Total Energy Density 326a515125bSLeila Ghaffari // ---- (E + P) u 3272b916ea7SJeremy 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]); 328139613f2SLeila Ghaffari 329139613f2SLeila Ghaffari // --Stabilization terms 330139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 331139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 332d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 333139613f2SLeila Ghaffari 334139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 335139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 336493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 337139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 338139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 3392b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 340139613f2SLeila Ghaffari } 341139613f2SLeila Ghaffari 342139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 343139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 3442b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3452b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 3462b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 3472b916ea7SJeremy L Thompson } 3482b916ea7SJeremy L Thompson } 349139613f2SLeila Ghaffari 350d8a22b9eSJed Brown // Stabilization 351d8a22b9eSJed Brown // -- Tau elements 352d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 353d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 354d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 355139613f2SLeila Ghaffari 356d8a22b9eSJed Brown // -- Stabilization method: none or SU 357bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 358139613f2SLeila Ghaffari switch (context->stabilization) { 359139613f2SLeila Ghaffari case 0: // Galerkin 360139613f2SLeila Ghaffari break; 361139613f2SLeila Ghaffari case 1: // SU 3622b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3632b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 3642b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 3652b916ea7SJeremy L Thompson } 3662b916ea7SJeremy L Thompson } 367139613f2SLeila Ghaffari 3682b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 3692b916ea7SJeremy 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]); 3702b916ea7SJeremy L Thompson } 371139613f2SLeila Ghaffari break; 372139613f2SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 373139613f2SLeila Ghaffari break; 374139613f2SLeila Ghaffari } 375b193fadcSJames Wright } 376a515125bSLeila Ghaffari return 0; 377a515125bSLeila Ghaffari } 378a515125bSLeila Ghaffari 379a515125bSLeila Ghaffari // ***************************************************************************** 38004e40bb6SJeremy L Thompson // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method 381a515125bSLeila Ghaffari // ***************************************************************************** 3822b916ea7SJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 3833d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 3843d65b166SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 3853d65b166SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 386ade49511SJames Wright const CeedScalar(*q_data) = in[3]; 3873d65b166SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 3883d65b166SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 38980f5d3cbSJames Wright CeedScalar *jac_data = out[2]; 390a515125bSLeila Ghaffari 391139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 392d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 393a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 39480f5d3cbSJames Wright const CeedScalar zeros[14] = {0.}; 395a515125bSLeila Ghaffari 3963d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 397a515125bSLeila Ghaffari // Setup 398a515125bSLeila Ghaffari // -- Interp in 399a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 4002b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 401a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 4022b916ea7SJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 4032b916ea7SJeremy L Thompson const CeedScalar dU[3][3] = { 4042b916ea7SJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 4052b916ea7SJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 4062b916ea7SJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 407139613f2SLeila Ghaffari }; 4082b916ea7SJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 409ade49511SJames Wright CeedScalar wdetJ, dXdx[3][3]; 410ade49511SJames Wright QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx); 411139613f2SLeila Ghaffari // dU/dx 412139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 413139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 414139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 415139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 416493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 417493642f1SJames Wright for (CeedInt k = 0; k < 3; k++) { 418139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 419139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 420493642f1SJames Wright for (CeedInt l = 0; l < 3; l++) { 421139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 422139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 423139613f2SLeila Ghaffari } 424139613f2SLeila Ghaffari } 425139613f2SLeila Ghaffari } 4262b916ea7SJeremy 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, 427139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 428a515125bSLeila Ghaffari 429a515125bSLeila Ghaffari // The Physics 430a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 431493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) { 432139613f2SLeila Ghaffari v[j][i] = 0.; 4332b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 434a515125bSLeila Ghaffari } 435a515125bSLeila Ghaffari //-----mass matrix 4362b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i]; 437a515125bSLeila Ghaffari 438a515125bSLeila Ghaffari // -- Density 439a515125bSLeila Ghaffari // ---- u rho 4402b916ea7SJeremy 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]); 441a515125bSLeila Ghaffari // -- Momentum 442a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 4432b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 4442b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { 4452b916ea7SJeremy 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] + 446139613f2SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 4472b916ea7SJeremy L Thompson } 4482b916ea7SJeremy L Thompson } 449a515125bSLeila Ghaffari // -- Total Energy Density 450a515125bSLeila Ghaffari // ---- (E + P) u 4512b916ea7SJeremy 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]); 452139613f2SLeila Ghaffari 453139613f2SLeila Ghaffari // -- Stabilization terms 454139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 455139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 456d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 457139613f2SLeila Ghaffari 458139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 459139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 460493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 461139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 462139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 4632b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 464139613f2SLeila Ghaffari } 465139613f2SLeila Ghaffari 466139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 467139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 4682b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 4692b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 4702b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 4712b916ea7SJeremy L Thompson } 4722b916ea7SJeremy L Thompson } 473139613f2SLeila Ghaffari 474139613f2SLeila Ghaffari // ---- Strong residual 475139613f2SLeila Ghaffari CeedScalar strong_res[5]; 4762b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j]; 477139613f2SLeila Ghaffari 478d8a22b9eSJed Brown // Stabilization 479d8a22b9eSJed Brown // -- Tau elements 480d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 481d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 482d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 483139613f2SLeila Ghaffari 484d8a22b9eSJed Brown // -- Stabilization method: none, SU, or SUPG 485bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 486139613f2SLeila Ghaffari switch (context->stabilization) { 487139613f2SLeila Ghaffari case 0: // Galerkin 488139613f2SLeila Ghaffari break; 489139613f2SLeila Ghaffari case 1: // SU 4902b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 4912b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 4922b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 4932b916ea7SJeremy L Thompson } 4942b916ea7SJeremy L Thompson } 495139613f2SLeila Ghaffari 4962b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 4972b916ea7SJeremy 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]); 4982b916ea7SJeremy L Thompson } 499139613f2SLeila Ghaffari break; 500139613f2SLeila Ghaffari case 2: // SUPG 5012b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 5022b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 5032b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l]; 5042b916ea7SJeremy L Thompson } 5052b916ea7SJeremy L Thompson } 506139613f2SLeila Ghaffari 5072b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 5082b916ea7SJeremy 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]); 5092b916ea7SJeremy L Thompson } 510139613f2SLeila Ghaffari break; 511139613f2SLeila Ghaffari } 51280f5d3cbSJames Wright StoredValuesPack(Q, i, 0, 14, zeros, jac_data); 513b193fadcSJames Wright } 514a515125bSLeila Ghaffari return 0; 515a515125bSLeila Ghaffari } 516a515125bSLeila Ghaffari // ***************************************************************************** 51704e40bb6SJeremy L Thompson // This QFunction sets the inflow boundary conditions for the traveling vortex problem. 518a515125bSLeila Ghaffari // 51904e40bb6SJeremy L Thompson // Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly. 520a515125bSLeila Ghaffari // ***************************************************************************** 5212b916ea7SJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 522ade49511SJames Wright const CeedScalar(*q_data_sur) = in[2]; 523a515125bSLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 524b193fadcSJames Wright 525a515125bSLeila Ghaffari EulerContext context = (EulerContext)ctx; 526a515125bSLeila Ghaffari const int euler_test = context->euler_test; 527ade49511SJames Wright const bool is_implicit = context->implicit; 528a515125bSLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 529a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 530a515125bSLeila Ghaffari const CeedScalar R = 1.; 531a515125bSLeila Ghaffari CeedScalar T_inlet; 532a515125bSLeila Ghaffari CeedScalar P_inlet; 533a515125bSLeila Ghaffari 534a515125bSLeila Ghaffari // For test cases 1 and 3 the background velocity is zero 5352b916ea7SJeremy L Thompson if (euler_test == 1 || euler_test == 3) { 536a515125bSLeila Ghaffari for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.; 5372b916ea7SJeremy L Thompson } 538a515125bSLeila Ghaffari 539a515125bSLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4 540a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 541a515125bSLeila Ghaffari else T_inlet = P_inlet = 1.; 542a515125bSLeila Ghaffari 5433d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 544*78e8b7daSJames Wright CeedScalar wdetJb, normal[3]; 545*78e8b7daSJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, normal); 546ade49511SJames Wright wdetJb *= is_implicit ? -1. : 1.; 547a515125bSLeila Ghaffari 548a515125bSLeila Ghaffari // face_normal = Normal vector of the face 549*78e8b7daSJames Wright const CeedScalar face_normal = Dot3(normal, mean_velocity); 550a515125bSLeila Ghaffari // The Physics 551a515125bSLeila Ghaffari // Zero v so all future terms can safely sum into it 552493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 553a515125bSLeila Ghaffari 554a515125bSLeila Ghaffari // Implementing in/outflow BCs 555002797a3SLeila Ghaffari if (face_normal > 0) { 556a515125bSLeila Ghaffari } else { // inflow 557a515125bSLeila Ghaffari const CeedScalar rho_inlet = P_inlet / (R * T_inlet); 5582b916ea7SJeremy L Thompson const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.; 559a515125bSLeila Ghaffari // incoming total energy 560a515125bSLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 561a515125bSLeila Ghaffari 562a515125bSLeila Ghaffari // The Physics 563a515125bSLeila Ghaffari // -- Density 564a515125bSLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal; 565a515125bSLeila Ghaffari 566a515125bSLeila Ghaffari // -- Momentum 567*78e8b7daSJames Wright for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho_inlet * face_normal * mean_velocity[j] + normal[j] * P_inlet); 568a515125bSLeila Ghaffari 569a515125bSLeila Ghaffari // -- Total Energy Density 570a515125bSLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 571a515125bSLeila Ghaffari } 572b193fadcSJames Wright } 573a515125bSLeila Ghaffari return 0; 574a515125bSLeila Ghaffari } 575a515125bSLeila Ghaffari 576a515125bSLeila Ghaffari // ***************************************************************************** 57704e40bb6SJeremy L Thompson // This QFunction sets the outflow boundary conditions for the Euler solver. 57868ef3d20SLeila Ghaffari // 57968ef3d20SLeila Ghaffari // Outflow BCs: 58004e40bb6SJeremy L Thompson // The validity of the weak form of the governing equations is extended to the outflow. 58168ef3d20SLeila Ghaffari // ***************************************************************************** 5822b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 5833d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 584ade49511SJames Wright const CeedScalar(*q_data_sur) = in[2]; 58568ef3d20SLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 586b193fadcSJames Wright 58768ef3d20SLeila Ghaffari EulerContext context = (EulerContext)ctx; 588ade49511SJames Wright const bool is_implicit = context->implicit; 58968ef3d20SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 59068ef3d20SLeila Ghaffari 59168ef3d20SLeila Ghaffari const CeedScalar gamma = 1.4; 59268ef3d20SLeila Ghaffari 5933d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 59468ef3d20SLeila Ghaffari // Setup 59568ef3d20SLeila Ghaffari // -- Interp in 59668ef3d20SLeila Ghaffari const CeedScalar rho = q[0][i]; 5972b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 59868ef3d20SLeila Ghaffari const CeedScalar E = q[4][i]; 59968ef3d20SLeila Ghaffari 600*78e8b7daSJames Wright CeedScalar wdetJb, normal[3]; 601*78e8b7daSJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, normal); 602ade49511SJames Wright wdetJb *= is_implicit ? -1. : 1.; 60368ef3d20SLeila Ghaffari 60468ef3d20SLeila Ghaffari // face_normal = Normal vector of the face 605*78e8b7daSJames Wright const CeedScalar face_normal = Dot3(normal, mean_velocity); 60668ef3d20SLeila Ghaffari // The Physics 60768ef3d20SLeila Ghaffari // Zero v so all future terms can safely sum into it 608493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0; 60968ef3d20SLeila Ghaffari 61068ef3d20SLeila Ghaffari // Implementing in/outflow BCs 61168ef3d20SLeila Ghaffari if (face_normal > 0) { // outflow 61268ef3d20SLeila Ghaffari const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.; 61368ef3d20SLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 614*78e8b7daSJames Wright const CeedScalar u_normal = Dot3(normal, u); // Normal velocity 61568ef3d20SLeila Ghaffari // The Physics 61668ef3d20SLeila Ghaffari // -- Density 61768ef3d20SLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal; 61868ef3d20SLeila Ghaffari 61968ef3d20SLeila Ghaffari // -- Momentum 620*78e8b7daSJames Wright for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + normal[j] * P); 62168ef3d20SLeila Ghaffari 62268ef3d20SLeila Ghaffari // -- Total Energy Density 62368ef3d20SLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P); 62468ef3d20SLeila Ghaffari } 625b193fadcSJames Wright } 62668ef3d20SLeila Ghaffari return 0; 62768ef3d20SLeila Ghaffari } 628