1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, 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: 13a515125bSLeila Ghaffari // On the Order of Accuracy and Numerical Performance of Two Classes of 14a515125bSLeila Ghaffari // Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011). 15a515125bSLeila Ghaffari 16a515125bSLeila Ghaffari #ifndef eulervortex_h 17a515125bSLeila Ghaffari #define eulervortex_h 18a515125bSLeila Ghaffari 19a515125bSLeila Ghaffari #include <math.h> 203a8779fbSJames Wright #include <ceed.h> 21*704b8bbeSJames Wright #include "utils.h" 22a515125bSLeila Ghaffari 23a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext; 24a515125bSLeila Ghaffari struct EulerContext_ { 25a515125bSLeila Ghaffari CeedScalar center[3]; 26a515125bSLeila Ghaffari CeedScalar curr_time; 27a515125bSLeila Ghaffari CeedScalar vortex_strength; 28d8a22b9eSJed Brown CeedScalar c_tau; 29a515125bSLeila Ghaffari CeedScalar mean_velocity[3]; 30a515125bSLeila Ghaffari bool implicit; 31139613f2SLeila Ghaffari int euler_test; 32139613f2SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 33a515125bSLeila Ghaffari }; 34a515125bSLeila Ghaffari 35a515125bSLeila Ghaffari // ***************************************************************************** 36a515125bSLeila Ghaffari // This function sets the initial conditions 37a515125bSLeila Ghaffari // 38a515125bSLeila Ghaffari // Temperature: 39a515125bSLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 40a515125bSLeila Ghaffari // Density: 41a515125bSLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1)) 42a515125bSLeila Ghaffari // Pressure: 43a515125bSLeila Ghaffari // P = rho * T 44a515125bSLeila Ghaffari // Velocity: 45a515125bSLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 46a515125bSLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 47a515125bSLeila Ghaffari // Velocity/Momentum Density: 48a515125bSLeila Ghaffari // Ui = rho ui 49a515125bSLeila Ghaffari // Total Energy: 50a515125bSLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2 51a515125bSLeila Ghaffari // 52a515125bSLeila Ghaffari // Constants: 53a515125bSLeila Ghaffari // cv , Specific heat, constant volume 54a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 55a515125bSLeila Ghaffari // vortex_strength , Strength of vortex 56a515125bSLeila Ghaffari // center , Location of bubble center 57a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 58a515125bSLeila Ghaffari // 59a515125bSLeila Ghaffari // ***************************************************************************** 60a515125bSLeila Ghaffari 61a515125bSLeila Ghaffari // ***************************************************************************** 62a515125bSLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 63a515125bSLeila Ghaffari // (currently not implemented) and IC formulation for Euler traveling vortex 64a515125bSLeila Ghaffari // ***************************************************************************** 65a515125bSLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, 66a515125bSLeila Ghaffari const CeedScalar X[], CeedInt Nf, CeedScalar q[], 67a515125bSLeila Ghaffari void *ctx) { 68a515125bSLeila Ghaffari // Context 69a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 70a515125bSLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength; 71a515125bSLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain 72a515125bSLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity; 73a515125bSLeila Ghaffari 74a515125bSLeila Ghaffari // Setup 75a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 76a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 77a515125bSLeila Ghaffari const CeedScalar R = 1.; 78a515125bSLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates 79a515125bSLeila Ghaffari // Vortex center 80a515125bSLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time; 81a515125bSLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time; 82a515125bSLeila Ghaffari 83a515125bSLeila Ghaffari const CeedScalar x0 = x - xc; 84a515125bSLeila Ghaffari const CeedScalar y0 = y - yc; 85a515125bSLeila Ghaffari const CeedScalar r = sqrt( x0*x0 + y0*y0 ); 86a515125bSLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r*r)/2.) / (2. * M_PI); 87139613f2SLeila Ghaffari const CeedScalar delta_T = - (gamma - 1.) * vortex_strength * vortex_strength * 88139613f2SLeila Ghaffari exp(1 - r*r) / (8. * gamma * M_PI * M_PI); 89a515125bSLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 90a515125bSLeila Ghaffari const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / 91a515125bSLeila Ghaffari (8.*gamma*M_PI*M_PI); 92a515125bSLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.}; 93a515125bSLeila Ghaffari 94a515125bSLeila Ghaffari // Initial Conditions 95a515125bSLeila Ghaffari switch (context->euler_test) { 96a515125bSLeila Ghaffari case 0: // Traveling vortex 97a515125bSLeila Ghaffari T = 1 + delta_T; 98a515125bSLeila Ghaffari // P = rho * T 99a515125bSLeila Ghaffari // P = S * rho^gamma 100a515125bSLeila Ghaffari // Solve for rho, then substitute for P 101139613f2SLeila Ghaffari rho = pow(T/S_vortex, 1 / (gamma - 1.)); 102a515125bSLeila Ghaffari P = rho * T; 103a515125bSLeila Ghaffari u[0] = mean_velocity[0] - C*y0; 104a515125bSLeila Ghaffari u[1] = mean_velocity[1] + C*x0; 105a515125bSLeila Ghaffari 106a515125bSLeila Ghaffari // Assign exact solution 107a515125bSLeila Ghaffari q[0] = rho; 108a515125bSLeila Ghaffari q[1] = rho * u[0]; 109a515125bSLeila Ghaffari q[2] = rho * u[1]; 110a515125bSLeila Ghaffari q[3] = rho * u[2]; 111a515125bSLeila Ghaffari q[4] = P / (gamma - 1.) + rho * (u[0]*u[0] + u[1]*u[1]) / 2.; 112a515125bSLeila Ghaffari break; 113a515125bSLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant 114a515125bSLeila Ghaffari rho = 1.; 115a515125bSLeila Ghaffari E = 2.; 116a515125bSLeila Ghaffari 117a515125bSLeila Ghaffari // Assign exact solution 118a515125bSLeila Ghaffari q[0] = rho; 119a515125bSLeila Ghaffari q[1] = rho * u[0]; 120a515125bSLeila Ghaffari q[2] = rho * u[1]; 121a515125bSLeila Ghaffari q[3] = rho * u[2]; 122a515125bSLeila Ghaffari q[4] = E; 123a515125bSLeila Ghaffari break; 124a515125bSLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant 125a515125bSLeila Ghaffari rho = 1.; 126a515125bSLeila Ghaffari E = 2.; 127a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 128a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 129a515125bSLeila Ghaffari 130a515125bSLeila Ghaffari // Assign exact solution 131a515125bSLeila Ghaffari q[0] = rho; 132a515125bSLeila Ghaffari q[1] = rho * u[0]; 133a515125bSLeila Ghaffari q[2] = rho * u[1]; 134a515125bSLeila Ghaffari q[3] = rho * u[2]; 135a515125bSLeila Ghaffari q[4] = E; 136a515125bSLeila Ghaffari break; 137a515125bSLeila Ghaffari case 3: // Velocity zero, pressure constant 138a515125bSLeila Ghaffari // (so density and internal energy will be non-constant), 139a515125bSLeila Ghaffari // but the velocity should stay zero and the bubble won't diffuse 140a515125bSLeila Ghaffari // (for Euler, where there is no thermal conductivity) 141a515125bSLeila Ghaffari P = 1.; 142a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 143a515125bSLeila Ghaffari rho = P / (R*T); 144a515125bSLeila Ghaffari 145a515125bSLeila Ghaffari // Assign exact solution 146a515125bSLeila Ghaffari q[0] = rho; 147a515125bSLeila Ghaffari q[1] = rho * u[0]; 148a515125bSLeila Ghaffari q[2] = rho * u[1]; 149a515125bSLeila Ghaffari q[3] = rho * u[2]; 150a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 151a515125bSLeila Ghaffari break; 152a515125bSLeila Ghaffari case 4: // Constant nonzero velocity, pressure constant 153a515125bSLeila Ghaffari // (so density and internal energy will be non-constant), 154a515125bSLeila Ghaffari // it should be transported across the domain, but velocity stays constant 155a515125bSLeila Ghaffari P = 1.; 156a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 157a515125bSLeila Ghaffari rho = P / (R*T); 158a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 159a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 160a515125bSLeila Ghaffari 161a515125bSLeila Ghaffari // Assign exact solution 162a515125bSLeila Ghaffari q[0] = rho; 163a515125bSLeila Ghaffari q[1] = rho * u[0]; 164a515125bSLeila Ghaffari q[2] = rho * u[1]; 165a515125bSLeila Ghaffari q[3] = rho * u[2]; 166a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 167a515125bSLeila Ghaffari break; 1680df2634dSLeila Ghaffari case 5: // non-smooth thermal bubble - cylinder 1690df2634dSLeila Ghaffari P = 1.; 1700df2634dSLeila Ghaffari T = 1. - (r < 1. ? S_bubble : 0.); 1710df2634dSLeila Ghaffari rho = P / (R*T); 1720df2634dSLeila Ghaffari u[0] = mean_velocity[0]; 1730df2634dSLeila Ghaffari u[1] = mean_velocity[1]; 1740df2634dSLeila Ghaffari 1750df2634dSLeila Ghaffari // Assign exact solution 1760df2634dSLeila Ghaffari q[0] = rho; 1770df2634dSLeila Ghaffari q[1] = rho * u[0]; 1780df2634dSLeila Ghaffari q[2] = rho * u[1]; 1790df2634dSLeila Ghaffari q[3] = rho * u[2]; 1800df2634dSLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 1810df2634dSLeila Ghaffari break; 182a515125bSLeila Ghaffari } 183a515125bSLeila Ghaffari // Return 184a515125bSLeila Ghaffari return 0; 185a515125bSLeila Ghaffari } 186a515125bSLeila Ghaffari 187a515125bSLeila Ghaffari // ***************************************************************************** 188139613f2SLeila Ghaffari // Helper function for computing flux Jacobian 189139613f2SLeila Ghaffari // ***************************************************************************** 190d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], 191139613f2SLeila Ghaffari const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 192139613f2SLeila Ghaffari const CeedScalar gamma) { 193139613f2SLeila Ghaffari CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 194139613f2SLeila Ghaffari for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 195139613f2SLeila Ghaffari for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix 196139613f2SLeila Ghaffari dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j]; 197139613f2SLeila Ghaffari for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix 198139613f2SLeila Ghaffari dF[i][0][k+1] = ((i==k) ? 1. : 0.); 199139613f2SLeila Ghaffari dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) + 200139613f2SLeila Ghaffari ((i==k) ? u[j] : 0.) - 201139613f2SLeila Ghaffari ((i==j) ? u[k] : 0.) * (gamma-1.); 202139613f2SLeila Ghaffari dF[i][4][k+1] = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) - 203139613f2SLeila Ghaffari (gamma-1.)*u[i]*u[k]; 204139613f2SLeila Ghaffari } 205139613f2SLeila Ghaffari dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.); 206139613f2SLeila Ghaffari } 207139613f2SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho); 208139613f2SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 209139613f2SLeila Ghaffari } 210139613f2SLeila Ghaffari } 211139613f2SLeila Ghaffari 212139613f2SLeila Ghaffari // ***************************************************************************** 213d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant) 214d8a22b9eSJed Brown // Model from: 215d8a22b9eSJed Brown // Stabilized Methods for Compressible Flows, Hughes et al 2010 216d8a22b9eSJed Brown // 217d8a22b9eSJed Brown // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 218d8a22b9eSJed Brown // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 219d8a22b9eSJed Brown // 220d8a22b9eSJed Brown // Where 221d8a22b9eSJed Brown // c_tau = stabilization constant (0.5 is reported as "optimal") 222d8a22b9eSJed Brown // h[i] = 2 length(dxdX[i]) 223d8a22b9eSJed Brown // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 224d8a22b9eSJed Brown // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 225d8a22b9eSJed Brown // rho(A[i]) = spectral radius of the convective flux Jacobian i, 226d8a22b9eSJed Brown // wave speed in direction i 227d8a22b9eSJed Brown // ***************************************************************************** 228d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], 229d8a22b9eSJed Brown const CeedScalar dXdx[3][3], const CeedScalar u[3], 230d8a22b9eSJed Brown const CeedScalar sound_speed, const CeedScalar c_tau) { 231493642f1SJames Wright for (CeedInt i=0; i<3; i++) { 232d8a22b9eSJed Brown // length of element in direction i 233d8a22b9eSJed Brown CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] + 234d8a22b9eSJed Brown dXdx[2][i]*dXdx[2][i]); 235d8a22b9eSJed Brown // fastest wave in direction i 236d8a22b9eSJed Brown CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 237d8a22b9eSJed Brown Tau_x[i] = c_tau * h / fastest_wave; 238d8a22b9eSJed Brown } 239d8a22b9eSJed Brown } 240d8a22b9eSJed Brown 241d8a22b9eSJed Brown // ***************************************************************************** 242a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex 243a515125bSLeila Ghaffari // ***************************************************************************** 244a515125bSLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, 245a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 246a515125bSLeila Ghaffari // Inputs 247a515125bSLeila Ghaffari const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 248a515125bSLeila Ghaffari 249a515125bSLeila Ghaffari // Outputs 250a515125bSLeila Ghaffari CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 251a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 252a515125bSLeila Ghaffari 253a515125bSLeila Ghaffari CeedPragmaSIMD 254a515125bSLeila Ghaffari // Quadrature Point Loop 255a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 256a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 257139613f2SLeila Ghaffari CeedScalar q[5] = {0.}; 258a515125bSLeila Ghaffari 259a515125bSLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx); 260a515125bSLeila Ghaffari 261a515125bSLeila Ghaffari for (CeedInt j=0; j<5; j++) 262a515125bSLeila Ghaffari q0[j][i] = q[j]; 263a515125bSLeila Ghaffari } // End of Quadrature Point Loop 264a515125bSLeila Ghaffari 265a515125bSLeila Ghaffari // Return 266a515125bSLeila Ghaffari return 0; 267a515125bSLeila Ghaffari } 268a515125bSLeila Ghaffari 269a515125bSLeila Ghaffari // ***************************************************************************** 270a515125bSLeila Ghaffari // This QFunction implements the following formulation of Euler equations 271a515125bSLeila Ghaffari // with explicit time stepping method 272a515125bSLeila Ghaffari // 273a515125bSLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation 274a515125bSLeila Ghaffari // form with state variables of density, momentum density, and total 275a515125bSLeila Ghaffari // energy density. 276a515125bSLeila Ghaffari // 277a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 278a515125bSLeila Ghaffari // rho - Mass Density 279a515125bSLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 280a515125bSLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 281a515125bSLeila Ghaffari // 282a515125bSLeila Ghaffari // Euler Equations: 283a515125bSLeila Ghaffari // drho/dt + div( U ) = 0 284a515125bSLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0 285a515125bSLeila Ghaffari // dE/dt + div( (E + P) u ) = 0 286a515125bSLeila Ghaffari // 287a515125bSLeila Ghaffari // Equation of State: 288a515125bSLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2) 289a515125bSLeila Ghaffari // 290a515125bSLeila Ghaffari // Constants: 291a515125bSLeila Ghaffari // cv , Specific heat, constant volume 292a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 293a515125bSLeila Ghaffari // g , Gravity 294a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 295a515125bSLeila Ghaffari // ***************************************************************************** 296a515125bSLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, 297a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 298a515125bSLeila Ghaffari // *INDENT-OFF* 299a515125bSLeila Ghaffari // Inputs 300a515125bSLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 301139613f2SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 302a515125bSLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 303a515125bSLeila Ghaffari // Outputs 304a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 305a515125bSLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 306a515125bSLeila Ghaffari 307139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 308d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 309a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 310a515125bSLeila Ghaffari 311a515125bSLeila Ghaffari CeedPragmaSIMD 312a515125bSLeila Ghaffari // Quadrature Point Loop 313a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 314a515125bSLeila Ghaffari // *INDENT-OFF* 315a515125bSLeila Ghaffari // Setup 316a515125bSLeila Ghaffari // -- Interp in 317a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 318a515125bSLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 319a515125bSLeila Ghaffari q[2][i] / rho, 320a515125bSLeila Ghaffari q[3][i] / rho 321a515125bSLeila Ghaffari }; 322a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 323139613f2SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 324139613f2SLeila Ghaffari dq[1][0][i], 325139613f2SLeila Ghaffari dq[2][0][i] 326139613f2SLeila Ghaffari }; 327139613f2SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 328139613f2SLeila Ghaffari dq[1][1][i], 329139613f2SLeila Ghaffari dq[2][1][i]}, 330139613f2SLeila Ghaffari {dq[0][2][i], 331139613f2SLeila Ghaffari dq[1][2][i], 332139613f2SLeila Ghaffari dq[2][2][i]}, 333139613f2SLeila Ghaffari {dq[0][3][i], 334139613f2SLeila Ghaffari dq[1][3][i], 335139613f2SLeila Ghaffari dq[2][3][i]} 336139613f2SLeila Ghaffari }; 337139613f2SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 338139613f2SLeila Ghaffari dq[1][4][i], 339139613f2SLeila Ghaffari dq[2][4][i] 340139613f2SLeila Ghaffari }; 341a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 342a515125bSLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 343a515125bSLeila Ghaffari // -- Interp-to-Grad q_data 344a515125bSLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 345a515125bSLeila Ghaffari // *INDENT-OFF* 346a515125bSLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 347a515125bSLeila Ghaffari q_data[2][i], 348a515125bSLeila Ghaffari q_data[3][i]}, 349a515125bSLeila Ghaffari {q_data[4][i], 350a515125bSLeila Ghaffari q_data[5][i], 351a515125bSLeila Ghaffari q_data[6][i]}, 352a515125bSLeila Ghaffari {q_data[7][i], 353a515125bSLeila Ghaffari q_data[8][i], 354a515125bSLeila Ghaffari q_data[9][i]} 355a515125bSLeila Ghaffari }; 356a515125bSLeila Ghaffari // *INDENT-ON* 357139613f2SLeila Ghaffari // dU/dx 358139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 359139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 360139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 361139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 362493642f1SJames Wright for (CeedInt j=0; j<3; j++) { 363493642f1SJames Wright for (CeedInt k=0; k<3; k++) { 364139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 365139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 366493642f1SJames Wright for (CeedInt l=0; l<3; l++) { 367139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 368139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 369139613f2SLeila Ghaffari } 370139613f2SLeila Ghaffari } 371139613f2SLeila Ghaffari } 372139613f2SLeila Ghaffari // Pressure 373a515125bSLeila Ghaffari const CeedScalar 374a515125bSLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 375a515125bSLeila Ghaffari E_internal = E - E_kinetic, 376139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 377a515125bSLeila Ghaffari 378a515125bSLeila Ghaffari // The Physics 379a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 380493642f1SJames Wright for (CeedInt j=0; j<5; j++) { 381139613f2SLeila Ghaffari v[j][i] = 0.; 382493642f1SJames Wright for (CeedInt k=0; k<3; k++) 383139613f2SLeila Ghaffari dv[k][j][i] = 0.; 384a515125bSLeila Ghaffari } 385a515125bSLeila Ghaffari 386a515125bSLeila Ghaffari // -- Density 387a515125bSLeila Ghaffari // ---- u rho 388493642f1SJames Wright for (CeedInt j=0; j<3; j++) 389a515125bSLeila Ghaffari dv[j][0][i] += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 390a515125bSLeila Ghaffari rho*u[2]*dXdx[j][2]); 391a515125bSLeila Ghaffari // -- Momentum 392a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 393493642f1SJames Wright for (CeedInt j=0; j<3; j++) 394493642f1SJames Wright for (CeedInt k=0; k<3; k++) 395139613f2SLeila Ghaffari dv[k][j+1][i] += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 396139613f2SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 397139613f2SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 398a515125bSLeila Ghaffari // -- Total Energy Density 399a515125bSLeila Ghaffari // ---- (E + P) u 400493642f1SJames Wright for (CeedInt j=0; j<3; j++) 401a515125bSLeila Ghaffari dv[j][4][i] += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 402a515125bSLeila Ghaffari u[2]*dXdx[j][2]); 403139613f2SLeila Ghaffari 404139613f2SLeila Ghaffari // --Stabilization terms 405139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 406139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 407d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 408139613f2SLeila Ghaffari 409139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 410139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 411493642f1SJames Wright for (CeedInt j=0; j<3; j++) { 412139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 413139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 414493642f1SJames Wright for (CeedInt k=0; k<3; k++) 415139613f2SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 416139613f2SLeila Ghaffari } 417139613f2SLeila Ghaffari 418139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 419139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 420493642f1SJames Wright for (CeedInt j=0; j<3; j++) 421493642f1SJames Wright for (CeedInt k=0; k<5; k++) 422493642f1SJames Wright for (CeedInt l=0; l<5; l++) 423139613f2SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 424139613f2SLeila Ghaffari 425d8a22b9eSJed Brown // Stabilization 426d8a22b9eSJed Brown // -- Tau elements 427d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 428d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 429d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 430139613f2SLeila Ghaffari 431d8a22b9eSJed Brown // -- Stabilization method: none or SU 432bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 433139613f2SLeila Ghaffari switch (context->stabilization) { 434139613f2SLeila Ghaffari case 0: // Galerkin 435139613f2SLeila Ghaffari break; 436139613f2SLeila Ghaffari case 1: // SU 437493642f1SJames Wright for (CeedInt j=0; j<3; j++) 438493642f1SJames Wright for (CeedInt k=0; k<5; k++) 439493642f1SJames Wright for (CeedInt l=0; l<5; l++) 440bb8a0c61SJames Wright stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 441139613f2SLeila Ghaffari 442493642f1SJames Wright for (CeedInt j=0; j<5; j++) 443493642f1SJames Wright for (CeedInt k=0; k<3; k++) 444139613f2SLeila Ghaffari dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] + 445139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 446139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 447139613f2SLeila Ghaffari break; 448139613f2SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 449139613f2SLeila Ghaffari break; 450139613f2SLeila Ghaffari } 451139613f2SLeila Ghaffari 452a515125bSLeila Ghaffari } // End Quadrature Point Loop 453a515125bSLeila Ghaffari 454a515125bSLeila Ghaffari // Return 455a515125bSLeila Ghaffari return 0; 456a515125bSLeila Ghaffari } 457a515125bSLeila Ghaffari 458a515125bSLeila Ghaffari // ***************************************************************************** 459a515125bSLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above) 460a515125bSLeila Ghaffari // with implicit time stepping method 461a515125bSLeila Ghaffari // 462a515125bSLeila Ghaffari // ***************************************************************************** 463a515125bSLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, 464a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 465a515125bSLeila Ghaffari // *INDENT-OFF* 466a515125bSLeila Ghaffari // Inputs 467a515125bSLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 468139613f2SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 469a515125bSLeila Ghaffari (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 470a515125bSLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 471a515125bSLeila Ghaffari // Outputs 472a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 473a515125bSLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 474a515125bSLeila Ghaffari 475139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 476d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 477a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 478a515125bSLeila Ghaffari 479a515125bSLeila Ghaffari CeedPragmaSIMD 480a515125bSLeila Ghaffari // Quadrature Point Loop 481a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 482a515125bSLeila Ghaffari // *INDENT-OFF* 483a515125bSLeila Ghaffari // Setup 484a515125bSLeila Ghaffari // -- Interp in 485a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 486a515125bSLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 487a515125bSLeila Ghaffari q[2][i] / rho, 488a515125bSLeila Ghaffari q[3][i] / rho 489a515125bSLeila Ghaffari }; 490a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 491139613f2SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 492139613f2SLeila Ghaffari dq[1][0][i], 493139613f2SLeila Ghaffari dq[2][0][i] 494139613f2SLeila Ghaffari }; 495139613f2SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 496139613f2SLeila Ghaffari dq[1][1][i], 497139613f2SLeila Ghaffari dq[2][1][i]}, 498139613f2SLeila Ghaffari {dq[0][2][i], 499139613f2SLeila Ghaffari dq[1][2][i], 500139613f2SLeila Ghaffari dq[2][2][i]}, 501139613f2SLeila Ghaffari {dq[0][3][i], 502139613f2SLeila Ghaffari dq[1][3][i], 503139613f2SLeila Ghaffari dq[2][3][i]} 504139613f2SLeila Ghaffari }; 505139613f2SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 506139613f2SLeila Ghaffari dq[1][4][i], 507139613f2SLeila Ghaffari dq[2][4][i] 508139613f2SLeila Ghaffari }; 509a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 510a515125bSLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 511a515125bSLeila Ghaffari // -- Interp-to-Grad q_data 512a515125bSLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 513a515125bSLeila Ghaffari // *INDENT-OFF* 514a515125bSLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 515a515125bSLeila Ghaffari q_data[2][i], 516a515125bSLeila Ghaffari q_data[3][i]}, 517a515125bSLeila Ghaffari {q_data[4][i], 518a515125bSLeila Ghaffari q_data[5][i], 519a515125bSLeila Ghaffari q_data[6][i]}, 520a515125bSLeila Ghaffari {q_data[7][i], 521a515125bSLeila Ghaffari q_data[8][i], 522a515125bSLeila Ghaffari q_data[9][i]} 523a515125bSLeila Ghaffari }; 524a515125bSLeila Ghaffari // *INDENT-ON* 525139613f2SLeila Ghaffari // dU/dx 526139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 527139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 528139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 529139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 530493642f1SJames Wright for (CeedInt j=0; j<3; j++) { 531493642f1SJames Wright for (CeedInt k=0; k<3; k++) { 532139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 533139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 534493642f1SJames Wright for (CeedInt l=0; l<3; l++) { 535139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 536139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 537139613f2SLeila Ghaffari } 538139613f2SLeila Ghaffari } 539139613f2SLeila Ghaffari } 540a515125bSLeila Ghaffari const CeedScalar 541a515125bSLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 542a515125bSLeila Ghaffari E_internal = E - E_kinetic, 543139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 544a515125bSLeila Ghaffari 545a515125bSLeila Ghaffari // The Physics 546a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 547493642f1SJames Wright for (CeedInt j=0; j<5; j++) { 548139613f2SLeila Ghaffari v[j][i] = 0.; 549493642f1SJames Wright for (CeedInt k=0; k<3; k++) 550139613f2SLeila Ghaffari dv[k][j][i] = 0.; 551a515125bSLeila Ghaffari } 552a515125bSLeila Ghaffari //-----mass matrix 553493642f1SJames Wright for (CeedInt j=0; j<5; j++) 554a515125bSLeila Ghaffari v[j][i] += wdetJ*q_dot[j][i]; 555a515125bSLeila Ghaffari 556a515125bSLeila Ghaffari // -- Density 557a515125bSLeila Ghaffari // ---- u rho 558493642f1SJames Wright for (CeedInt j=0; j<3; j++) 559a515125bSLeila Ghaffari dv[j][0][i] -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 560a515125bSLeila Ghaffari rho*u[2]*dXdx[j][2]); 561a515125bSLeila Ghaffari // -- Momentum 562a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 563493642f1SJames Wright for (CeedInt j=0; j<3; j++) 564493642f1SJames Wright for (CeedInt k=0; k<3; k++) 565139613f2SLeila Ghaffari dv[k][j+1][i] -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 566139613f2SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 567139613f2SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 568a515125bSLeila Ghaffari // -- Total Energy Density 569a515125bSLeila Ghaffari // ---- (E + P) u 570493642f1SJames Wright for (CeedInt j=0; j<3; j++) 571a515125bSLeila Ghaffari dv[j][4][i] -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 572a515125bSLeila Ghaffari u[2]*dXdx[j][2]); 573139613f2SLeila Ghaffari 574139613f2SLeila Ghaffari // -- Stabilization terms 575139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 576139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 577d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 578139613f2SLeila Ghaffari 579139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 580139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 581493642f1SJames Wright for (CeedInt j=0; j<3; j++) { 582139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 583139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 584493642f1SJames Wright for (CeedInt k=0; k<3; k++) 585139613f2SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 586139613f2SLeila Ghaffari } 587139613f2SLeila Ghaffari 588139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 589139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 590493642f1SJames Wright for (CeedInt j=0; j<3; j++) 591493642f1SJames Wright for (CeedInt k=0; k<5; k++) 592493642f1SJames Wright for (CeedInt l=0; l<5; l++) 593139613f2SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 594139613f2SLeila Ghaffari 595139613f2SLeila Ghaffari // ---- Strong residual 596139613f2SLeila Ghaffari CeedScalar strong_res[5]; 597493642f1SJames Wright for (CeedInt j=0; j<5; j++) 598139613f2SLeila Ghaffari strong_res[j] = q_dot[j][i] + strong_conv[j]; 599139613f2SLeila Ghaffari 600d8a22b9eSJed Brown // Stabilization 601d8a22b9eSJed Brown // -- Tau elements 602d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 603d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 604d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 605139613f2SLeila Ghaffari 606d8a22b9eSJed Brown // -- Stabilization method: none, SU, or SUPG 607bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 608139613f2SLeila Ghaffari switch (context->stabilization) { 609139613f2SLeila Ghaffari case 0: // Galerkin 610139613f2SLeila Ghaffari break; 611139613f2SLeila Ghaffari case 1: // SU 612493642f1SJames Wright for (CeedInt j=0; j<3; j++) 613493642f1SJames Wright for (CeedInt k=0; k<5; k++) 614493642f1SJames Wright for (CeedInt l=0; l<5; l++) 615bb8a0c61SJames Wright stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 616139613f2SLeila Ghaffari 617493642f1SJames Wright for (CeedInt j=0; j<5; j++) 618493642f1SJames Wright for (CeedInt k=0; k<3; k++) 619139613f2SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 620139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 621139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 622139613f2SLeila Ghaffari break; 623139613f2SLeila Ghaffari case 2: // SUPG 624493642f1SJames Wright for (CeedInt j=0; j<3; j++) 625493642f1SJames Wright for (CeedInt k=0; k<5; k++) 626493642f1SJames Wright for (CeedInt l=0; l<5; l++) 627bb8a0c61SJames Wright stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l]; 628139613f2SLeila Ghaffari 629493642f1SJames Wright for (CeedInt j=0; j<5; j++) 630493642f1SJames Wright for (CeedInt k=0; k<3; k++) 631139613f2SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 632139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 633139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 634139613f2SLeila Ghaffari break; 635139613f2SLeila Ghaffari } 636a515125bSLeila Ghaffari } // End Quadrature Point Loop 637a515125bSLeila Ghaffari 638a515125bSLeila Ghaffari // Return 639a515125bSLeila Ghaffari return 0; 640a515125bSLeila Ghaffari } 641a515125bSLeila Ghaffari // ***************************************************************************** 642002797a3SLeila Ghaffari // This QFunction sets the inflow boundary conditions for 643002797a3SLeila Ghaffari // the traveling vortex problem. 644a515125bSLeila Ghaffari // 645a515125bSLeila Ghaffari // Prescribed T_inlet and P_inlet are converted to conservative variables 646a515125bSLeila Ghaffari // and applied weakly. 647a515125bSLeila Ghaffari // 648a515125bSLeila Ghaffari // ***************************************************************************** 649002797a3SLeila Ghaffari CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, 650a515125bSLeila Ghaffari const CeedScalar *const *in, 651a515125bSLeila Ghaffari CeedScalar *const *out) { 652a515125bSLeila Ghaffari // *INDENT-OFF* 653a515125bSLeila Ghaffari // Inputs 654dd64951cSJames Wright const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 655a515125bSLeila Ghaffari // Outputs 656a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 657a515125bSLeila Ghaffari // *INDENT-ON* 658a515125bSLeila Ghaffari EulerContext context = (EulerContext)ctx; 659a515125bSLeila Ghaffari const int euler_test = context->euler_test; 660a515125bSLeila Ghaffari const bool implicit = context->implicit; 661a515125bSLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 662a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 663a515125bSLeila Ghaffari const CeedScalar R = 1.; 664a515125bSLeila Ghaffari CeedScalar T_inlet; 665a515125bSLeila Ghaffari CeedScalar P_inlet; 666a515125bSLeila Ghaffari 667a515125bSLeila Ghaffari // For test cases 1 and 3 the background velocity is zero 668a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 3) 669a515125bSLeila Ghaffari for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.; 670a515125bSLeila Ghaffari 671a515125bSLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4 672a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 673a515125bSLeila Ghaffari else T_inlet = P_inlet = 1.; 674a515125bSLeila Ghaffari 675a515125bSLeila Ghaffari CeedPragmaSIMD 676a515125bSLeila Ghaffari // Quadrature Point Loop 677a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 678a515125bSLeila Ghaffari // Setup 679a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 680a515125bSLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 681a515125bSLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 682a515125bSLeila Ghaffari // We can effect this by swapping the sign on this weight 683a515125bSLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 684002797a3SLeila Ghaffari // ---- Normal vect 685a515125bSLeila Ghaffari const CeedScalar norm[3] = {q_data_sur[1][i], 686a515125bSLeila Ghaffari q_data_sur[2][i], 687a515125bSLeila Ghaffari q_data_sur[3][i] 688a515125bSLeila Ghaffari }; 689a515125bSLeila Ghaffari 690a515125bSLeila Ghaffari // face_normal = Normal vector of the face 691a515125bSLeila Ghaffari const CeedScalar face_normal = norm[0]*mean_velocity[0] + 692a515125bSLeila Ghaffari norm[1]*mean_velocity[1] + 693a515125bSLeila Ghaffari norm[2]*mean_velocity[2]; 694a515125bSLeila Ghaffari // The Physics 695a515125bSLeila Ghaffari // Zero v so all future terms can safely sum into it 696493642f1SJames Wright for (CeedInt j=0; j<5; j++) v[j][i] = 0.; 697a515125bSLeila Ghaffari 698a515125bSLeila Ghaffari // Implementing in/outflow BCs 699002797a3SLeila Ghaffari if (face_normal > 0) { 700a515125bSLeila Ghaffari } else { // inflow 701a515125bSLeila Ghaffari const CeedScalar rho_inlet = P_inlet/(R*T_inlet); 702a515125bSLeila Ghaffari const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] + 703a515125bSLeila Ghaffari mean_velocity[1]*mean_velocity[1]) / 2.; 704a515125bSLeila Ghaffari // incoming total energy 705a515125bSLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 706a515125bSLeila Ghaffari 707a515125bSLeila Ghaffari // The Physics 708a515125bSLeila Ghaffari // -- Density 709a515125bSLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal; 710a515125bSLeila Ghaffari 711a515125bSLeila Ghaffari // -- Momentum 712493642f1SJames Wright for (CeedInt j=0; j<3; j++) 713a515125bSLeila Ghaffari v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] + 714a515125bSLeila Ghaffari norm[j] * P_inlet); 715a515125bSLeila Ghaffari 716a515125bSLeila Ghaffari // -- Total Energy Density 717a515125bSLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 718a515125bSLeila Ghaffari } 719a515125bSLeila Ghaffari 720a515125bSLeila Ghaffari } // End Quadrature Point Loop 721a515125bSLeila Ghaffari return 0; 722a515125bSLeila Ghaffari } 723a515125bSLeila Ghaffari 724a515125bSLeila Ghaffari // ***************************************************************************** 72568ef3d20SLeila Ghaffari // This QFunction sets the outflow boundary conditions for 72668ef3d20SLeila Ghaffari // the Euler solver. 72768ef3d20SLeila Ghaffari // 72868ef3d20SLeila Ghaffari // Outflow BCs: 72968ef3d20SLeila Ghaffari // The validity of the weak form of the governing equations is 73068ef3d20SLeila Ghaffari // extended to the outflow. 73168ef3d20SLeila Ghaffari // 73268ef3d20SLeila Ghaffari // ***************************************************************************** 73368ef3d20SLeila Ghaffari CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, 73468ef3d20SLeila Ghaffari const CeedScalar *const *in, 73568ef3d20SLeila Ghaffari CeedScalar *const *out) { 73668ef3d20SLeila Ghaffari // *INDENT-OFF* 73768ef3d20SLeila Ghaffari // Inputs 73868ef3d20SLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 739dd64951cSJames Wright (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 74068ef3d20SLeila Ghaffari // Outputs 74168ef3d20SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 74268ef3d20SLeila Ghaffari // *INDENT-ON* 74368ef3d20SLeila Ghaffari EulerContext context = (EulerContext)ctx; 74468ef3d20SLeila Ghaffari const bool implicit = context->implicit; 74568ef3d20SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 74668ef3d20SLeila Ghaffari 74768ef3d20SLeila Ghaffari const CeedScalar gamma = 1.4; 74868ef3d20SLeila Ghaffari 74968ef3d20SLeila Ghaffari CeedPragmaSIMD 75068ef3d20SLeila Ghaffari // Quadrature Point Loop 75168ef3d20SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 75268ef3d20SLeila Ghaffari // Setup 75368ef3d20SLeila Ghaffari // -- Interp in 75468ef3d20SLeila Ghaffari const CeedScalar rho = q[0][i]; 75568ef3d20SLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 75668ef3d20SLeila Ghaffari q[2][i] / rho, 75768ef3d20SLeila Ghaffari q[3][i] / rho 75868ef3d20SLeila Ghaffari }; 75968ef3d20SLeila Ghaffari const CeedScalar E = q[4][i]; 76068ef3d20SLeila Ghaffari 76168ef3d20SLeila Ghaffari // -- Interp-to-Interp q_data 76268ef3d20SLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 76368ef3d20SLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 76468ef3d20SLeila Ghaffari // We can effect this by swapping the sign on this weight 76568ef3d20SLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 76668ef3d20SLeila Ghaffari // ---- Normal vectors 76768ef3d20SLeila Ghaffari const CeedScalar norm[3] = {q_data_sur[1][i], 76868ef3d20SLeila Ghaffari q_data_sur[2][i], 76968ef3d20SLeila Ghaffari q_data_sur[3][i] 77068ef3d20SLeila Ghaffari }; 77168ef3d20SLeila Ghaffari 77268ef3d20SLeila Ghaffari // face_normal = Normal vector of the face 77368ef3d20SLeila Ghaffari const CeedScalar face_normal = norm[0]*mean_velocity[0] + 77468ef3d20SLeila Ghaffari norm[1]*mean_velocity[1] + 77568ef3d20SLeila Ghaffari norm[2]*mean_velocity[2]; 77668ef3d20SLeila Ghaffari // The Physics 77768ef3d20SLeila Ghaffari // Zero v so all future terms can safely sum into it 778493642f1SJames Wright for (CeedInt j=0; j<5; j++) v[j][i] = 0; 77968ef3d20SLeila Ghaffari 78068ef3d20SLeila Ghaffari // Implementing in/outflow BCs 78168ef3d20SLeila Ghaffari if (face_normal > 0) { // outflow 78268ef3d20SLeila Ghaffari const CeedScalar E_kinetic = (u[0]*u[0] + u[1]*u[1]) / 2.; 78368ef3d20SLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 78468ef3d20SLeila Ghaffari const CeedScalar u_normal = norm[0]*u[0] + norm[1]*u[1] + 78568ef3d20SLeila Ghaffari norm[2]*u[2]; // Normal velocity 78668ef3d20SLeila Ghaffari // The Physics 78768ef3d20SLeila Ghaffari // -- Density 78868ef3d20SLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal; 78968ef3d20SLeila Ghaffari 79068ef3d20SLeila Ghaffari // -- Momentum 791493642f1SJames Wright for (CeedInt j=0; j<3; j++) 79268ef3d20SLeila Ghaffari v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + norm[j] * P); 79368ef3d20SLeila Ghaffari 79468ef3d20SLeila Ghaffari // -- Total Energy Density 79568ef3d20SLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P); 79668ef3d20SLeila Ghaffari } 79768ef3d20SLeila Ghaffari } // End Quadrature Point Loop 79868ef3d20SLeila Ghaffari return 0; 79968ef3d20SLeila Ghaffari } 80068ef3d20SLeila Ghaffari 80168ef3d20SLeila Ghaffari // ***************************************************************************** 802a515125bSLeila Ghaffari 803a515125bSLeila Ghaffari #endif // eulervortex_h 804