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> 21a515125bSLeila Ghaffari 22a515125bSLeila Ghaffari #ifndef M_PI 23a515125bSLeila Ghaffari #define M_PI 3.14159265358979323846 24a515125bSLeila Ghaffari #endif 25a515125bSLeila Ghaffari 26a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext; 27a515125bSLeila Ghaffari struct EulerContext_ { 28a515125bSLeila Ghaffari CeedScalar center[3]; 29a515125bSLeila Ghaffari CeedScalar curr_time; 30a515125bSLeila Ghaffari CeedScalar vortex_strength; 31d8a22b9eSJed Brown CeedScalar c_tau; 32a515125bSLeila Ghaffari CeedScalar mean_velocity[3]; 33a515125bSLeila Ghaffari bool implicit; 34139613f2SLeila Ghaffari int euler_test; 35139613f2SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 36a515125bSLeila Ghaffari }; 37a515125bSLeila Ghaffari 38a515125bSLeila Ghaffari // ***************************************************************************** 39a515125bSLeila Ghaffari // This function sets the initial conditions 40a515125bSLeila Ghaffari // 41a515125bSLeila Ghaffari // Temperature: 42a515125bSLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 43a515125bSLeila Ghaffari // Density: 44a515125bSLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1)) 45a515125bSLeila Ghaffari // Pressure: 46a515125bSLeila Ghaffari // P = rho * T 47a515125bSLeila Ghaffari // Velocity: 48a515125bSLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 49a515125bSLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 50a515125bSLeila Ghaffari // Velocity/Momentum Density: 51a515125bSLeila Ghaffari // Ui = rho ui 52a515125bSLeila Ghaffari // Total Energy: 53a515125bSLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2 54a515125bSLeila Ghaffari // 55a515125bSLeila Ghaffari // Constants: 56a515125bSLeila Ghaffari // cv , Specific heat, constant volume 57a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 58a515125bSLeila Ghaffari // vortex_strength , Strength of vortex 59a515125bSLeila Ghaffari // center , Location of bubble center 60a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 61a515125bSLeila Ghaffari // 62a515125bSLeila Ghaffari // ***************************************************************************** 63a515125bSLeila Ghaffari 64a515125bSLeila Ghaffari // ***************************************************************************** 65a515125bSLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 66a515125bSLeila Ghaffari // (currently not implemented) and IC formulation for Euler traveling vortex 67a515125bSLeila Ghaffari // ***************************************************************************** 68a515125bSLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, 69a515125bSLeila Ghaffari const CeedScalar X[], CeedInt Nf, CeedScalar q[], 70a515125bSLeila Ghaffari void *ctx) { 71a515125bSLeila Ghaffari // Context 72a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 73a515125bSLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength; 74a515125bSLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain 75a515125bSLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity; 76a515125bSLeila Ghaffari 77a515125bSLeila Ghaffari // Setup 78a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 79a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 80a515125bSLeila Ghaffari const CeedScalar R = 1.; 81a515125bSLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates 82a515125bSLeila Ghaffari // Vortex center 83a515125bSLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time; 84a515125bSLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time; 85a515125bSLeila Ghaffari 86a515125bSLeila Ghaffari const CeedScalar x0 = x - xc; 87a515125bSLeila Ghaffari const CeedScalar y0 = y - yc; 88a515125bSLeila Ghaffari const CeedScalar r = sqrt( x0*x0 + y0*y0 ); 89a515125bSLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r*r)/2.) / (2. * M_PI); 90139613f2SLeila Ghaffari const CeedScalar delta_T = - (gamma - 1.) * vortex_strength * vortex_strength * 91139613f2SLeila Ghaffari exp(1 - r*r) / (8. * gamma * M_PI * M_PI); 92a515125bSLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 93a515125bSLeila Ghaffari const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / 94a515125bSLeila Ghaffari (8.*gamma*M_PI*M_PI); 95a515125bSLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.}; 96a515125bSLeila Ghaffari 97a515125bSLeila Ghaffari // Initial Conditions 98a515125bSLeila Ghaffari switch (context->euler_test) { 99a515125bSLeila Ghaffari case 0: // Traveling vortex 100a515125bSLeila Ghaffari T = 1 + delta_T; 101a515125bSLeila Ghaffari // P = rho * T 102a515125bSLeila Ghaffari // P = S * rho^gamma 103a515125bSLeila Ghaffari // Solve for rho, then substitute for P 104139613f2SLeila Ghaffari rho = pow(T/S_vortex, 1 / (gamma - 1.)); 105a515125bSLeila Ghaffari P = rho * T; 106a515125bSLeila Ghaffari u[0] = mean_velocity[0] - C*y0; 107a515125bSLeila Ghaffari u[1] = mean_velocity[1] + C*x0; 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] = P / (gamma - 1.) + rho * (u[0]*u[0] + u[1]*u[1]) / 2.; 115a515125bSLeila Ghaffari break; 116a515125bSLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant 117a515125bSLeila Ghaffari rho = 1.; 118a515125bSLeila Ghaffari E = 2.; 119a515125bSLeila Ghaffari 120a515125bSLeila Ghaffari // Assign exact solution 121a515125bSLeila Ghaffari q[0] = rho; 122a515125bSLeila Ghaffari q[1] = rho * u[0]; 123a515125bSLeila Ghaffari q[2] = rho * u[1]; 124a515125bSLeila Ghaffari q[3] = rho * u[2]; 125a515125bSLeila Ghaffari q[4] = E; 126a515125bSLeila Ghaffari break; 127a515125bSLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant 128a515125bSLeila Ghaffari rho = 1.; 129a515125bSLeila Ghaffari E = 2.; 130a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 131a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 132a515125bSLeila Ghaffari 133a515125bSLeila Ghaffari // Assign exact solution 134a515125bSLeila Ghaffari q[0] = rho; 135a515125bSLeila Ghaffari q[1] = rho * u[0]; 136a515125bSLeila Ghaffari q[2] = rho * u[1]; 137a515125bSLeila Ghaffari q[3] = rho * u[2]; 138a515125bSLeila Ghaffari q[4] = E; 139a515125bSLeila Ghaffari break; 140a515125bSLeila Ghaffari case 3: // Velocity zero, pressure constant 141a515125bSLeila Ghaffari // (so density and internal energy will be non-constant), 142a515125bSLeila Ghaffari // but the velocity should stay zero and the bubble won't diffuse 143a515125bSLeila Ghaffari // (for Euler, where there is no thermal conductivity) 144a515125bSLeila Ghaffari P = 1.; 145a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 146a515125bSLeila Ghaffari rho = P / (R*T); 147a515125bSLeila Ghaffari 148a515125bSLeila Ghaffari // Assign exact solution 149a515125bSLeila Ghaffari q[0] = rho; 150a515125bSLeila Ghaffari q[1] = rho * u[0]; 151a515125bSLeila Ghaffari q[2] = rho * u[1]; 152a515125bSLeila Ghaffari q[3] = rho * u[2]; 153a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 154a515125bSLeila Ghaffari break; 155a515125bSLeila Ghaffari case 4: // Constant nonzero velocity, pressure constant 156a515125bSLeila Ghaffari // (so density and internal energy will be non-constant), 157a515125bSLeila Ghaffari // it should be transported across the domain, but velocity stays constant 158a515125bSLeila Ghaffari P = 1.; 159a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 160a515125bSLeila Ghaffari rho = P / (R*T); 161a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 162a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 163a515125bSLeila Ghaffari 164a515125bSLeila Ghaffari // Assign exact solution 165a515125bSLeila Ghaffari q[0] = rho; 166a515125bSLeila Ghaffari q[1] = rho * u[0]; 167a515125bSLeila Ghaffari q[2] = rho * u[1]; 168a515125bSLeila Ghaffari q[3] = rho * u[2]; 169a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 170a515125bSLeila Ghaffari break; 1710df2634dSLeila Ghaffari case 5: // non-smooth thermal bubble - cylinder 1720df2634dSLeila Ghaffari P = 1.; 1730df2634dSLeila Ghaffari T = 1. - (r < 1. ? S_bubble : 0.); 1740df2634dSLeila Ghaffari rho = P / (R*T); 1750df2634dSLeila Ghaffari u[0] = mean_velocity[0]; 1760df2634dSLeila Ghaffari u[1] = mean_velocity[1]; 1770df2634dSLeila Ghaffari 1780df2634dSLeila Ghaffari // Assign exact solution 1790df2634dSLeila Ghaffari q[0] = rho; 1800df2634dSLeila Ghaffari q[1] = rho * u[0]; 1810df2634dSLeila Ghaffari q[2] = rho * u[1]; 1820df2634dSLeila Ghaffari q[3] = rho * u[2]; 1830df2634dSLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 1840df2634dSLeila Ghaffari break; 185a515125bSLeila Ghaffari } 186a515125bSLeila Ghaffari // Return 187a515125bSLeila Ghaffari return 0; 188a515125bSLeila Ghaffari } 189a515125bSLeila Ghaffari 190a515125bSLeila Ghaffari // ***************************************************************************** 191139613f2SLeila Ghaffari // Helper function for computing flux Jacobian 192139613f2SLeila Ghaffari // ***************************************************************************** 193d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], 194139613f2SLeila Ghaffari const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 195139613f2SLeila Ghaffari const CeedScalar gamma) { 196139613f2SLeila Ghaffari CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 197139613f2SLeila Ghaffari for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 198139613f2SLeila Ghaffari for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix 199139613f2SLeila Ghaffari dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j]; 200139613f2SLeila Ghaffari for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix 201139613f2SLeila Ghaffari dF[i][0][k+1] = ((i==k) ? 1. : 0.); 202139613f2SLeila Ghaffari dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) + 203139613f2SLeila Ghaffari ((i==k) ? u[j] : 0.) - 204139613f2SLeila Ghaffari ((i==j) ? u[k] : 0.) * (gamma-1.); 205139613f2SLeila Ghaffari dF[i][4][k+1] = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) - 206139613f2SLeila Ghaffari (gamma-1.)*u[i]*u[k]; 207139613f2SLeila Ghaffari } 208139613f2SLeila Ghaffari dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.); 209139613f2SLeila Ghaffari } 210139613f2SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho); 211139613f2SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 212139613f2SLeila Ghaffari } 213139613f2SLeila Ghaffari } 214139613f2SLeila Ghaffari 215139613f2SLeila Ghaffari // ***************************************************************************** 216d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant) 217d8a22b9eSJed Brown // Model from: 218d8a22b9eSJed Brown // Stabilized Methods for Compressible Flows, Hughes et al 2010 219d8a22b9eSJed Brown // 220d8a22b9eSJed Brown // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 221d8a22b9eSJed Brown // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 222d8a22b9eSJed Brown // 223d8a22b9eSJed Brown // Where 224d8a22b9eSJed Brown // c_tau = stabilization constant (0.5 is reported as "optimal") 225d8a22b9eSJed Brown // h[i] = 2 length(dxdX[i]) 226d8a22b9eSJed Brown // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 227d8a22b9eSJed Brown // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 228d8a22b9eSJed Brown // rho(A[i]) = spectral radius of the convective flux Jacobian i, 229d8a22b9eSJed Brown // wave speed in direction i 230d8a22b9eSJed Brown // ***************************************************************************** 231d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], 232d8a22b9eSJed Brown const CeedScalar dXdx[3][3], const CeedScalar u[3], 233d8a22b9eSJed Brown const CeedScalar sound_speed, const CeedScalar c_tau) { 234*493642f1SJames Wright for (CeedInt i=0; i<3; i++) { 235d8a22b9eSJed Brown // length of element in direction i 236d8a22b9eSJed Brown CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] + 237d8a22b9eSJed Brown dXdx[2][i]*dXdx[2][i]); 238d8a22b9eSJed Brown // fastest wave in direction i 239d8a22b9eSJed Brown CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 240d8a22b9eSJed Brown Tau_x[i] = c_tau * h / fastest_wave; 241d8a22b9eSJed Brown } 242d8a22b9eSJed Brown } 243d8a22b9eSJed Brown 244d8a22b9eSJed Brown // ***************************************************************************** 245a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex 246a515125bSLeila Ghaffari // ***************************************************************************** 247a515125bSLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, 248a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 249a515125bSLeila Ghaffari // Inputs 250a515125bSLeila Ghaffari const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 251a515125bSLeila Ghaffari 252a515125bSLeila Ghaffari // Outputs 253a515125bSLeila Ghaffari CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 254a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 255a515125bSLeila Ghaffari 256a515125bSLeila Ghaffari CeedPragmaSIMD 257a515125bSLeila Ghaffari // Quadrature Point Loop 258a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 259a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 260139613f2SLeila Ghaffari CeedScalar q[5] = {0.}; 261a515125bSLeila Ghaffari 262a515125bSLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx); 263a515125bSLeila Ghaffari 264a515125bSLeila Ghaffari for (CeedInt j=0; j<5; j++) 265a515125bSLeila Ghaffari q0[j][i] = q[j]; 266a515125bSLeila Ghaffari } // End of Quadrature Point Loop 267a515125bSLeila Ghaffari 268a515125bSLeila Ghaffari // Return 269a515125bSLeila Ghaffari return 0; 270a515125bSLeila Ghaffari } 271a515125bSLeila Ghaffari 272a515125bSLeila Ghaffari // ***************************************************************************** 273a515125bSLeila Ghaffari // This QFunction implements the following formulation of Euler equations 274a515125bSLeila Ghaffari // with explicit time stepping method 275a515125bSLeila Ghaffari // 276a515125bSLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation 277a515125bSLeila Ghaffari // form with state variables of density, momentum density, and total 278a515125bSLeila Ghaffari // energy density. 279a515125bSLeila Ghaffari // 280a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 281a515125bSLeila Ghaffari // rho - Mass Density 282a515125bSLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 283a515125bSLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 284a515125bSLeila Ghaffari // 285a515125bSLeila Ghaffari // Euler Equations: 286a515125bSLeila Ghaffari // drho/dt + div( U ) = 0 287a515125bSLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0 288a515125bSLeila Ghaffari // dE/dt + div( (E + P) u ) = 0 289a515125bSLeila Ghaffari // 290a515125bSLeila Ghaffari // Equation of State: 291a515125bSLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2) 292a515125bSLeila Ghaffari // 293a515125bSLeila Ghaffari // Constants: 294a515125bSLeila Ghaffari // cv , Specific heat, constant volume 295a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 296a515125bSLeila Ghaffari // g , Gravity 297a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 298a515125bSLeila Ghaffari // ***************************************************************************** 299a515125bSLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, 300a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 301a515125bSLeila Ghaffari // *INDENT-OFF* 302a515125bSLeila Ghaffari // Inputs 303a515125bSLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 304139613f2SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 305a515125bSLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 306a515125bSLeila Ghaffari // Outputs 307a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 308a515125bSLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 309a515125bSLeila Ghaffari 310139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 311d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 312a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 313a515125bSLeila Ghaffari 314a515125bSLeila Ghaffari CeedPragmaSIMD 315a515125bSLeila Ghaffari // Quadrature Point Loop 316a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 317a515125bSLeila Ghaffari // *INDENT-OFF* 318a515125bSLeila Ghaffari // Setup 319a515125bSLeila Ghaffari // -- Interp in 320a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 321a515125bSLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 322a515125bSLeila Ghaffari q[2][i] / rho, 323a515125bSLeila Ghaffari q[3][i] / rho 324a515125bSLeila Ghaffari }; 325a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 326139613f2SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 327139613f2SLeila Ghaffari dq[1][0][i], 328139613f2SLeila Ghaffari dq[2][0][i] 329139613f2SLeila Ghaffari }; 330139613f2SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 331139613f2SLeila Ghaffari dq[1][1][i], 332139613f2SLeila Ghaffari dq[2][1][i]}, 333139613f2SLeila Ghaffari {dq[0][2][i], 334139613f2SLeila Ghaffari dq[1][2][i], 335139613f2SLeila Ghaffari dq[2][2][i]}, 336139613f2SLeila Ghaffari {dq[0][3][i], 337139613f2SLeila Ghaffari dq[1][3][i], 338139613f2SLeila Ghaffari dq[2][3][i]} 339139613f2SLeila Ghaffari }; 340139613f2SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 341139613f2SLeila Ghaffari dq[1][4][i], 342139613f2SLeila Ghaffari dq[2][4][i] 343139613f2SLeila Ghaffari }; 344a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 345a515125bSLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 346a515125bSLeila Ghaffari // -- Interp-to-Grad q_data 347a515125bSLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 348a515125bSLeila Ghaffari // *INDENT-OFF* 349a515125bSLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 350a515125bSLeila Ghaffari q_data[2][i], 351a515125bSLeila Ghaffari q_data[3][i]}, 352a515125bSLeila Ghaffari {q_data[4][i], 353a515125bSLeila Ghaffari q_data[5][i], 354a515125bSLeila Ghaffari q_data[6][i]}, 355a515125bSLeila Ghaffari {q_data[7][i], 356a515125bSLeila Ghaffari q_data[8][i], 357a515125bSLeila Ghaffari q_data[9][i]} 358a515125bSLeila Ghaffari }; 359a515125bSLeila Ghaffari // *INDENT-ON* 360139613f2SLeila Ghaffari // dU/dx 361139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 362139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 363139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 364139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 365*493642f1SJames Wright for (CeedInt j=0; j<3; j++) { 366*493642f1SJames Wright for (CeedInt k=0; k<3; k++) { 367139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 368139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 369*493642f1SJames Wright for (CeedInt l=0; l<3; l++) { 370139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 371139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 372139613f2SLeila Ghaffari } 373139613f2SLeila Ghaffari } 374139613f2SLeila Ghaffari } 375139613f2SLeila Ghaffari // Pressure 376a515125bSLeila Ghaffari const CeedScalar 377a515125bSLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 378a515125bSLeila Ghaffari E_internal = E - E_kinetic, 379139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 380a515125bSLeila Ghaffari 381a515125bSLeila Ghaffari // The Physics 382a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 383*493642f1SJames Wright for (CeedInt j=0; j<5; j++) { 384139613f2SLeila Ghaffari v[j][i] = 0.; 385*493642f1SJames Wright for (CeedInt k=0; k<3; k++) 386139613f2SLeila Ghaffari dv[k][j][i] = 0.; 387a515125bSLeila Ghaffari } 388a515125bSLeila Ghaffari 389a515125bSLeila Ghaffari // -- Density 390a515125bSLeila Ghaffari // ---- u rho 391*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 392a515125bSLeila Ghaffari dv[j][0][i] += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 393a515125bSLeila Ghaffari rho*u[2]*dXdx[j][2]); 394a515125bSLeila Ghaffari // -- Momentum 395a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 396*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 397*493642f1SJames Wright for (CeedInt k=0; k<3; k++) 398139613f2SLeila Ghaffari dv[k][j+1][i] += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 399139613f2SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 400139613f2SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 401a515125bSLeila Ghaffari // -- Total Energy Density 402a515125bSLeila Ghaffari // ---- (E + P) u 403*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 404a515125bSLeila Ghaffari dv[j][4][i] += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 405a515125bSLeila Ghaffari u[2]*dXdx[j][2]); 406139613f2SLeila Ghaffari 407139613f2SLeila Ghaffari // --Stabilization terms 408139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 409139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 410d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 411139613f2SLeila Ghaffari 412139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 413139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 414*493642f1SJames Wright for (CeedInt j=0; j<3; j++) { 415139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 416139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 417*493642f1SJames Wright for (CeedInt k=0; k<3; k++) 418139613f2SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 419139613f2SLeila Ghaffari } 420139613f2SLeila Ghaffari 421139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 422139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 423*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 424*493642f1SJames Wright for (CeedInt k=0; k<5; k++) 425*493642f1SJames Wright for (CeedInt l=0; l<5; l++) 426139613f2SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 427139613f2SLeila Ghaffari 428d8a22b9eSJed Brown // Stabilization 429d8a22b9eSJed Brown // -- Tau elements 430d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 431d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 432d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 433139613f2SLeila Ghaffari 434d8a22b9eSJed Brown // -- Stabilization method: none or SU 435bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 436139613f2SLeila Ghaffari switch (context->stabilization) { 437139613f2SLeila Ghaffari case 0: // Galerkin 438139613f2SLeila Ghaffari break; 439139613f2SLeila Ghaffari case 1: // SU 440*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 441*493642f1SJames Wright for (CeedInt k=0; k<5; k++) 442*493642f1SJames Wright for (CeedInt l=0; l<5; l++) 443bb8a0c61SJames Wright stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 444139613f2SLeila Ghaffari 445*493642f1SJames Wright for (CeedInt j=0; j<5; j++) 446*493642f1SJames Wright for (CeedInt k=0; k<3; k++) 447139613f2SLeila Ghaffari dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] + 448139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 449139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 450139613f2SLeila Ghaffari break; 451139613f2SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 452139613f2SLeila Ghaffari break; 453139613f2SLeila Ghaffari } 454139613f2SLeila Ghaffari 455a515125bSLeila Ghaffari } // End Quadrature Point Loop 456a515125bSLeila Ghaffari 457a515125bSLeila Ghaffari // Return 458a515125bSLeila Ghaffari return 0; 459a515125bSLeila Ghaffari } 460a515125bSLeila Ghaffari 461a515125bSLeila Ghaffari // ***************************************************************************** 462a515125bSLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above) 463a515125bSLeila Ghaffari // with implicit time stepping method 464a515125bSLeila Ghaffari // 465a515125bSLeila Ghaffari // ***************************************************************************** 466a515125bSLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, 467a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 468a515125bSLeila Ghaffari // *INDENT-OFF* 469a515125bSLeila Ghaffari // Inputs 470a515125bSLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 471139613f2SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 472a515125bSLeila Ghaffari (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 473a515125bSLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 474a515125bSLeila Ghaffari // Outputs 475a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 476a515125bSLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 477a515125bSLeila Ghaffari 478139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 479d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 480a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 481a515125bSLeila Ghaffari 482a515125bSLeila Ghaffari CeedPragmaSIMD 483a515125bSLeila Ghaffari // Quadrature Point Loop 484a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 485a515125bSLeila Ghaffari // *INDENT-OFF* 486a515125bSLeila Ghaffari // Setup 487a515125bSLeila Ghaffari // -- Interp in 488a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 489a515125bSLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 490a515125bSLeila Ghaffari q[2][i] / rho, 491a515125bSLeila Ghaffari q[3][i] / rho 492a515125bSLeila Ghaffari }; 493a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 494139613f2SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 495139613f2SLeila Ghaffari dq[1][0][i], 496139613f2SLeila Ghaffari dq[2][0][i] 497139613f2SLeila Ghaffari }; 498139613f2SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 499139613f2SLeila Ghaffari dq[1][1][i], 500139613f2SLeila Ghaffari dq[2][1][i]}, 501139613f2SLeila Ghaffari {dq[0][2][i], 502139613f2SLeila Ghaffari dq[1][2][i], 503139613f2SLeila Ghaffari dq[2][2][i]}, 504139613f2SLeila Ghaffari {dq[0][3][i], 505139613f2SLeila Ghaffari dq[1][3][i], 506139613f2SLeila Ghaffari dq[2][3][i]} 507139613f2SLeila Ghaffari }; 508139613f2SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 509139613f2SLeila Ghaffari dq[1][4][i], 510139613f2SLeila Ghaffari dq[2][4][i] 511139613f2SLeila Ghaffari }; 512a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 513a515125bSLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 514a515125bSLeila Ghaffari // -- Interp-to-Grad q_data 515a515125bSLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 516a515125bSLeila Ghaffari // *INDENT-OFF* 517a515125bSLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 518a515125bSLeila Ghaffari q_data[2][i], 519a515125bSLeila Ghaffari q_data[3][i]}, 520a515125bSLeila Ghaffari {q_data[4][i], 521a515125bSLeila Ghaffari q_data[5][i], 522a515125bSLeila Ghaffari q_data[6][i]}, 523a515125bSLeila Ghaffari {q_data[7][i], 524a515125bSLeila Ghaffari q_data[8][i], 525a515125bSLeila Ghaffari q_data[9][i]} 526a515125bSLeila Ghaffari }; 527a515125bSLeila Ghaffari // *INDENT-ON* 528139613f2SLeila Ghaffari // dU/dx 529139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 530139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 531139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 532139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 533*493642f1SJames Wright for (CeedInt j=0; j<3; j++) { 534*493642f1SJames Wright for (CeedInt k=0; k<3; k++) { 535139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 536139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 537*493642f1SJames Wright for (CeedInt l=0; l<3; l++) { 538139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 539139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 540139613f2SLeila Ghaffari } 541139613f2SLeila Ghaffari } 542139613f2SLeila Ghaffari } 543a515125bSLeila Ghaffari const CeedScalar 544a515125bSLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 545a515125bSLeila Ghaffari E_internal = E - E_kinetic, 546139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 547a515125bSLeila Ghaffari 548a515125bSLeila Ghaffari // The Physics 549a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 550*493642f1SJames Wright for (CeedInt j=0; j<5; j++) { 551139613f2SLeila Ghaffari v[j][i] = 0.; 552*493642f1SJames Wright for (CeedInt k=0; k<3; k++) 553139613f2SLeila Ghaffari dv[k][j][i] = 0.; 554a515125bSLeila Ghaffari } 555a515125bSLeila Ghaffari //-----mass matrix 556*493642f1SJames Wright for (CeedInt j=0; j<5; j++) 557a515125bSLeila Ghaffari v[j][i] += wdetJ*q_dot[j][i]; 558a515125bSLeila Ghaffari 559a515125bSLeila Ghaffari // -- Density 560a515125bSLeila Ghaffari // ---- u rho 561*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 562a515125bSLeila Ghaffari dv[j][0][i] -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 563a515125bSLeila Ghaffari rho*u[2]*dXdx[j][2]); 564a515125bSLeila Ghaffari // -- Momentum 565a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 566*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 567*493642f1SJames Wright for (CeedInt k=0; k<3; k++) 568139613f2SLeila Ghaffari dv[k][j+1][i] -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 569139613f2SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 570139613f2SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 571a515125bSLeila Ghaffari // -- Total Energy Density 572a515125bSLeila Ghaffari // ---- (E + P) u 573*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 574a515125bSLeila Ghaffari dv[j][4][i] -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 575a515125bSLeila Ghaffari u[2]*dXdx[j][2]); 576139613f2SLeila Ghaffari 577139613f2SLeila Ghaffari // -- Stabilization terms 578139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 579139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 580d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 581139613f2SLeila Ghaffari 582139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 583139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 584*493642f1SJames Wright for (CeedInt j=0; j<3; j++) { 585139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 586139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 587*493642f1SJames Wright for (CeedInt k=0; k<3; k++) 588139613f2SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 589139613f2SLeila Ghaffari } 590139613f2SLeila Ghaffari 591139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 592139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 593*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 594*493642f1SJames Wright for (CeedInt k=0; k<5; k++) 595*493642f1SJames Wright for (CeedInt l=0; l<5; l++) 596139613f2SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 597139613f2SLeila Ghaffari 598139613f2SLeila Ghaffari // ---- Strong residual 599139613f2SLeila Ghaffari CeedScalar strong_res[5]; 600*493642f1SJames Wright for (CeedInt j=0; j<5; j++) 601139613f2SLeila Ghaffari strong_res[j] = q_dot[j][i] + strong_conv[j]; 602139613f2SLeila Ghaffari 603d8a22b9eSJed Brown // Stabilization 604d8a22b9eSJed Brown // -- Tau elements 605d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 606d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 607d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 608139613f2SLeila Ghaffari 609d8a22b9eSJed Brown // -- Stabilization method: none, SU, or SUPG 610bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 611139613f2SLeila Ghaffari switch (context->stabilization) { 612139613f2SLeila Ghaffari case 0: // Galerkin 613139613f2SLeila Ghaffari break; 614139613f2SLeila Ghaffari case 1: // SU 615*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 616*493642f1SJames Wright for (CeedInt k=0; k<5; k++) 617*493642f1SJames Wright for (CeedInt l=0; l<5; l++) 618bb8a0c61SJames Wright stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 619139613f2SLeila Ghaffari 620*493642f1SJames Wright for (CeedInt j=0; j<5; j++) 621*493642f1SJames Wright for (CeedInt k=0; k<3; k++) 622139613f2SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 623139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 624139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 625139613f2SLeila Ghaffari break; 626139613f2SLeila Ghaffari case 2: // SUPG 627*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 628*493642f1SJames Wright for (CeedInt k=0; k<5; k++) 629*493642f1SJames Wright for (CeedInt l=0; l<5; l++) 630bb8a0c61SJames Wright stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l]; 631139613f2SLeila Ghaffari 632*493642f1SJames Wright for (CeedInt j=0; j<5; j++) 633*493642f1SJames Wright for (CeedInt k=0; k<3; k++) 634139613f2SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 635139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 636139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 637139613f2SLeila Ghaffari break; 638139613f2SLeila Ghaffari } 639a515125bSLeila Ghaffari } // End Quadrature Point Loop 640a515125bSLeila Ghaffari 641a515125bSLeila Ghaffari // Return 642a515125bSLeila Ghaffari return 0; 643a515125bSLeila Ghaffari } 644a515125bSLeila Ghaffari // ***************************************************************************** 645002797a3SLeila Ghaffari // This QFunction sets the inflow boundary conditions for 646002797a3SLeila Ghaffari // the traveling vortex problem. 647a515125bSLeila Ghaffari // 648a515125bSLeila Ghaffari // Prescribed T_inlet and P_inlet are converted to conservative variables 649a515125bSLeila Ghaffari // and applied weakly. 650a515125bSLeila Ghaffari // 651a515125bSLeila Ghaffari // ***************************************************************************** 652002797a3SLeila Ghaffari CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, 653a515125bSLeila Ghaffari const CeedScalar *const *in, 654a515125bSLeila Ghaffari CeedScalar *const *out) { 655a515125bSLeila Ghaffari // *INDENT-OFF* 656a515125bSLeila Ghaffari // Inputs 657002797a3SLeila Ghaffari const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 658a515125bSLeila Ghaffari // Outputs 659a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 660a515125bSLeila Ghaffari // *INDENT-ON* 661a515125bSLeila Ghaffari EulerContext context = (EulerContext)ctx; 662a515125bSLeila Ghaffari const int euler_test = context->euler_test; 663a515125bSLeila Ghaffari const bool implicit = context->implicit; 664a515125bSLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 665a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 666a515125bSLeila Ghaffari const CeedScalar R = 1.; 667a515125bSLeila Ghaffari CeedScalar T_inlet; 668a515125bSLeila Ghaffari CeedScalar P_inlet; 669a515125bSLeila Ghaffari 670a515125bSLeila Ghaffari // For test cases 1 and 3 the background velocity is zero 671a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 3) 672a515125bSLeila Ghaffari for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.; 673a515125bSLeila Ghaffari 674a515125bSLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4 675a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 676a515125bSLeila Ghaffari else T_inlet = P_inlet = 1.; 677a515125bSLeila Ghaffari 678a515125bSLeila Ghaffari CeedPragmaSIMD 679a515125bSLeila Ghaffari // Quadrature Point Loop 680a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 681a515125bSLeila Ghaffari // Setup 682a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 683a515125bSLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 684a515125bSLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 685a515125bSLeila Ghaffari // We can effect this by swapping the sign on this weight 686a515125bSLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 687002797a3SLeila Ghaffari // ---- Normal vect 688a515125bSLeila Ghaffari const CeedScalar norm[3] = {q_data_sur[1][i], 689a515125bSLeila Ghaffari q_data_sur[2][i], 690a515125bSLeila Ghaffari q_data_sur[3][i] 691a515125bSLeila Ghaffari }; 692a515125bSLeila Ghaffari 693a515125bSLeila Ghaffari // face_normal = Normal vector of the face 694a515125bSLeila Ghaffari const CeedScalar face_normal = norm[0]*mean_velocity[0] + 695a515125bSLeila Ghaffari norm[1]*mean_velocity[1] + 696a515125bSLeila Ghaffari norm[2]*mean_velocity[2]; 697a515125bSLeila Ghaffari // The Physics 698a515125bSLeila Ghaffari // Zero v so all future terms can safely sum into it 699*493642f1SJames Wright for (CeedInt j=0; j<5; j++) v[j][i] = 0.; 700a515125bSLeila Ghaffari 701a515125bSLeila Ghaffari // Implementing in/outflow BCs 702002797a3SLeila Ghaffari if (face_normal > 0) { 703a515125bSLeila Ghaffari } else { // inflow 704a515125bSLeila Ghaffari const CeedScalar rho_inlet = P_inlet/(R*T_inlet); 705a515125bSLeila Ghaffari const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] + 706a515125bSLeila Ghaffari mean_velocity[1]*mean_velocity[1]) / 2.; 707a515125bSLeila Ghaffari // incoming total energy 708a515125bSLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 709a515125bSLeila Ghaffari 710a515125bSLeila Ghaffari // The Physics 711a515125bSLeila Ghaffari // -- Density 712a515125bSLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal; 713a515125bSLeila Ghaffari 714a515125bSLeila Ghaffari // -- Momentum 715*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 716a515125bSLeila Ghaffari v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] + 717a515125bSLeila Ghaffari norm[j] * P_inlet); 718a515125bSLeila Ghaffari 719a515125bSLeila Ghaffari // -- Total Energy Density 720a515125bSLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 721a515125bSLeila Ghaffari } 722a515125bSLeila Ghaffari 723a515125bSLeila Ghaffari } // End Quadrature Point Loop 724a515125bSLeila Ghaffari return 0; 725a515125bSLeila Ghaffari } 726a515125bSLeila Ghaffari 727a515125bSLeila Ghaffari // ***************************************************************************** 72868ef3d20SLeila Ghaffari // This QFunction sets the outflow boundary conditions for 72968ef3d20SLeila Ghaffari // the Euler solver. 73068ef3d20SLeila Ghaffari // 73168ef3d20SLeila Ghaffari // Outflow BCs: 73268ef3d20SLeila Ghaffari // The validity of the weak form of the governing equations is 73368ef3d20SLeila Ghaffari // extended to the outflow. 73468ef3d20SLeila Ghaffari // 73568ef3d20SLeila Ghaffari // ***************************************************************************** 73668ef3d20SLeila Ghaffari CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, 73768ef3d20SLeila Ghaffari const CeedScalar *const *in, 73868ef3d20SLeila Ghaffari CeedScalar *const *out) { 73968ef3d20SLeila Ghaffari // *INDENT-OFF* 74068ef3d20SLeila Ghaffari // Inputs 74168ef3d20SLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 74268ef3d20SLeila Ghaffari (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 74368ef3d20SLeila Ghaffari // Outputs 74468ef3d20SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 74568ef3d20SLeila Ghaffari // *INDENT-ON* 74668ef3d20SLeila Ghaffari EulerContext context = (EulerContext)ctx; 74768ef3d20SLeila Ghaffari const bool implicit = context->implicit; 74868ef3d20SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 74968ef3d20SLeila Ghaffari 75068ef3d20SLeila Ghaffari const CeedScalar gamma = 1.4; 75168ef3d20SLeila Ghaffari 75268ef3d20SLeila Ghaffari CeedPragmaSIMD 75368ef3d20SLeila Ghaffari // Quadrature Point Loop 75468ef3d20SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 75568ef3d20SLeila Ghaffari // Setup 75668ef3d20SLeila Ghaffari // -- Interp in 75768ef3d20SLeila Ghaffari const CeedScalar rho = q[0][i]; 75868ef3d20SLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 75968ef3d20SLeila Ghaffari q[2][i] / rho, 76068ef3d20SLeila Ghaffari q[3][i] / rho 76168ef3d20SLeila Ghaffari }; 76268ef3d20SLeila Ghaffari const CeedScalar E = q[4][i]; 76368ef3d20SLeila Ghaffari 76468ef3d20SLeila Ghaffari // -- Interp-to-Interp q_data 76568ef3d20SLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 76668ef3d20SLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 76768ef3d20SLeila Ghaffari // We can effect this by swapping the sign on this weight 76868ef3d20SLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 76968ef3d20SLeila Ghaffari // ---- Normal vectors 77068ef3d20SLeila Ghaffari const CeedScalar norm[3] = {q_data_sur[1][i], 77168ef3d20SLeila Ghaffari q_data_sur[2][i], 77268ef3d20SLeila Ghaffari q_data_sur[3][i] 77368ef3d20SLeila Ghaffari }; 77468ef3d20SLeila Ghaffari 77568ef3d20SLeila Ghaffari // face_normal = Normal vector of the face 77668ef3d20SLeila Ghaffari const CeedScalar face_normal = norm[0]*mean_velocity[0] + 77768ef3d20SLeila Ghaffari norm[1]*mean_velocity[1] + 77868ef3d20SLeila Ghaffari norm[2]*mean_velocity[2]; 77968ef3d20SLeila Ghaffari // The Physics 78068ef3d20SLeila Ghaffari // Zero v so all future terms can safely sum into it 781*493642f1SJames Wright for (CeedInt j=0; j<5; j++) v[j][i] = 0; 78268ef3d20SLeila Ghaffari 78368ef3d20SLeila Ghaffari // Implementing in/outflow BCs 78468ef3d20SLeila Ghaffari if (face_normal > 0) { // outflow 78568ef3d20SLeila Ghaffari const CeedScalar E_kinetic = (u[0]*u[0] + u[1]*u[1]) / 2.; 78668ef3d20SLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 78768ef3d20SLeila Ghaffari const CeedScalar u_normal = norm[0]*u[0] + norm[1]*u[1] + 78868ef3d20SLeila Ghaffari norm[2]*u[2]; // Normal velocity 78968ef3d20SLeila Ghaffari // The Physics 79068ef3d20SLeila Ghaffari // -- Density 79168ef3d20SLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal; 79268ef3d20SLeila Ghaffari 79368ef3d20SLeila Ghaffari // -- Momentum 794*493642f1SJames Wright for (CeedInt j=0; j<3; j++) 79568ef3d20SLeila Ghaffari v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + norm[j] * P); 79668ef3d20SLeila Ghaffari 79768ef3d20SLeila Ghaffari // -- Total Energy Density 79868ef3d20SLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P); 79968ef3d20SLeila Ghaffari } 80068ef3d20SLeila Ghaffari } // End Quadrature Point Loop 80168ef3d20SLeila Ghaffari return 0; 80268ef3d20SLeila Ghaffari } 80368ef3d20SLeila Ghaffari 80468ef3d20SLeila Ghaffari // ***************************************************************************** 805a515125bSLeila Ghaffari 806a515125bSLeila Ghaffari #endif // eulervortex_h 807