1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Euler traveling vortex initial condition and operator for Navier-Stokes 6 /// example using PETSc 7 8 // Model from: 9 // On the Order of Accuracy and Numerical Performance of Two Classes of Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011). 10 #include <ceed/types.h> 11 #ifndef CEED_RUNNING_JIT_PASS 12 #include <stdbool.h> 13 #endif 14 15 #include "utils.h" 16 17 typedef struct EulerContext_ *EulerContext; 18 struct EulerContext_ { 19 CeedScalar center[3]; 20 CeedScalar curr_time; 21 CeedScalar vortex_strength; 22 CeedScalar c_tau; 23 CeedScalar mean_velocity[3]; 24 bool implicit; 25 int euler_test; 26 int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 27 }; 28 29 // ***************************************************************************** 30 // This function sets the initial conditions 31 // 32 // Temperature: 33 // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 34 // Density: 35 // rho = (T/S_vortex)^(1 / (gamma - 1)) 36 // Pressure: 37 // P = rho * T 38 // Velocity: 39 // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 40 // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 41 // Velocity/Momentum Density: 42 // Ui = rho ui 43 // Total Energy: 44 // E = P / (gamma - 1) + rho (u u)/2 45 // 46 // Constants: 47 // cv , Specific heat, constant volume 48 // cp , Specific heat, constant pressure 49 // vortex_strength , Strength of vortex 50 // center , Location of bubble center 51 // gamma = cp / cv, Specific heat ratio 52 // 53 // ***************************************************************************** 54 55 // ***************************************************************************** 56 // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling 57 // vortex 58 // ***************************************************************************** 59 CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) { 60 // Context 61 const EulerContext context = (EulerContext)ctx; 62 const CeedScalar vortex_strength = context->vortex_strength; 63 const CeedScalar *center = context->center; // Center of the domain 64 const CeedScalar *mean_velocity = context->mean_velocity; 65 66 // Setup 67 const CeedScalar gamma = 1.4; 68 const CeedScalar cv = 2.5; 69 const CeedScalar R = 1.; 70 const CeedScalar x = X[0], y = X[1]; // Coordinates 71 // Vortex center 72 const CeedScalar xc = center[0] + mean_velocity[0] * time; 73 const CeedScalar yc = center[1] + mean_velocity[1] * time; 74 75 const CeedScalar x0 = x - xc; 76 const CeedScalar y0 = y - yc; 77 const CeedScalar r = sqrt(Square(x0) + Square(y0)); 78 const CeedScalar C = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI); 79 const CeedScalar delta_T = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI); 80 const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 81 const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI); 82 CeedScalar rho, P, T, E, u[3] = {0.}; 83 84 // Initial Conditions 85 switch (context->euler_test) { 86 case 0: // Traveling vortex 87 T = 1 + delta_T; 88 // P = rho * T 89 // P = S * rho^gamma 90 // Solve for rho, then substitute for P 91 rho = pow(T / S_vortex, 1 / (gamma - 1.)); 92 P = rho * T; 93 u[0] = mean_velocity[0] - C * y0; 94 u[1] = mean_velocity[1] + C * x0; 95 96 // Assign exact solution 97 q[0] = rho; 98 q[1] = rho * u[0]; 99 q[2] = rho * u[1]; 100 q[3] = rho * u[2]; 101 q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.; 102 break; 103 case 1: // Constant zero velocity, density constant, total energy constant 104 rho = 1.; 105 E = 2.; 106 107 // Assign exact solution 108 q[0] = rho; 109 q[1] = rho * u[0]; 110 q[2] = rho * u[1]; 111 q[3] = rho * u[2]; 112 q[4] = E; 113 break; 114 case 2: // Constant nonzero velocity, density constant, total energy constant 115 rho = 1.; 116 E = 2.; 117 u[0] = mean_velocity[0]; 118 u[1] = mean_velocity[1]; 119 120 // Assign exact solution 121 q[0] = rho; 122 q[1] = rho * u[0]; 123 q[2] = rho * u[1]; 124 q[3] = rho * u[2]; 125 q[4] = E; 126 break; 127 case 3: // Velocity zero, pressure constant (so density and internal energy will be non-constant), but the velocity should stay zero and the 128 // bubble won't diffuse 129 // (for Euler, where there is no thermal conductivity) 130 P = 1.; 131 T = 1. - S_bubble * exp(1. - r * r); 132 rho = P / (R * T); 133 134 // Assign exact solution 135 q[0] = rho; 136 q[1] = rho * u[0]; 137 q[2] = rho * u[1]; 138 q[3] = rho * u[2]; 139 q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 140 break; 141 case 4: // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant), 142 // It should be transported across the domain, but velocity stays constant 143 P = 1.; 144 T = 1. - S_bubble * exp(1. - r * r); 145 rho = P / (R * T); 146 u[0] = mean_velocity[0]; 147 u[1] = mean_velocity[1]; 148 149 // Assign exact solution 150 q[0] = rho; 151 q[1] = rho * u[0]; 152 q[2] = rho * u[1]; 153 q[3] = rho * u[2]; 154 q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 155 break; 156 case 5: // non-smooth thermal bubble - cylinder 157 P = 1.; 158 T = 1. - (r < 1. ? S_bubble : 0.); 159 rho = P / (R * T); 160 u[0] = mean_velocity[0]; 161 u[1] = mean_velocity[1]; 162 163 // Assign exact solution 164 q[0] = rho; 165 q[1] = rho * u[0]; 166 q[2] = rho * u[1]; 167 q[3] = rho * u[2]; 168 q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 169 break; 170 } 171 return 0; 172 } 173 174 // ***************************************************************************** 175 // Helper function for computing flux Jacobian 176 // ***************************************************************************** 177 CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 178 const CeedScalar gamma) { 179 CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2]; // Velocity square 180 for (CeedInt i = 0; i < 3; i++) { // Jacobian matrices for 3 directions 181 for (CeedInt j = 0; j < 3; j++) { // Rows of each Jacobian matrix 182 dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j]; 183 for (CeedInt k = 0; k < 3; k++) { // Columns of each Jacobian matrix 184 dF[i][0][k + 1] = ((i == k) ? 1. : 0.); 185 dF[i][j + 1][k + 1] = ((j == k) ? u[i] : 0.) + ((i == k) ? u[j] : 0.) - ((i == j) ? u[k] : 0.) * (gamma - 1.); 186 dF[i][4][k + 1] = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k]; 187 } 188 dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.); 189 } 190 dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho); 191 dF[i][4][4] = u[i] * gamma; 192 } 193 } 194 195 // ***************************************************************************** 196 // Helper function for computing Tau elements (stabilization constant) 197 // Model from: 198 // Stabilized Methods for Compressible Flows, Hughes et al 2010 199 // 200 // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 201 // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 202 // 203 // Where 204 // c_tau = stabilization constant (0.5 is reported as "optimal") 205 // h[i] = 2 length(dxdX[i]) 206 // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 207 // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 208 // rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i 209 // ***************************************************************************** 210 CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], const CeedScalar dXdx[3][3], const CeedScalar u[3], const CeedScalar sound_speed, 211 const CeedScalar c_tau) { 212 for (CeedInt i = 0; i < 3; i++) { 213 // length of element in direction i 214 CeedScalar h = 2 / sqrt(Square(dXdx[0][i]) + Square(dXdx[1][i]) + Square(dXdx[2][i])); 215 // fastest wave in direction i 216 CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 217 Tau_x[i] = c_tau * h / fastest_wave; 218 } 219 } 220 221 // ***************************************************************************** 222 // This QFunction sets the initial conditions for Euler traveling vortex 223 // ***************************************************************************** 224 CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 225 const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 226 CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 227 228 const EulerContext context = (EulerContext)ctx; 229 230 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 231 const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 232 CeedScalar q[5] = {0.}; 233 234 Exact_Euler(3, context->curr_time, x, 5, q, ctx); 235 for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 236 } 237 return 0; 238 } 239 240 // ***************************************************************************** 241 // This QFunction implements the following formulation of Euler equations with explicit time stepping method 242 // 243 // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density. 244 // 245 // State Variables: q = ( rho, U1, U2, U3, E ) 246 // rho - Mass Density 247 // Ui - Momentum Density, Ui = rho ui 248 // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 249 // 250 // Euler Equations: 251 // drho/dt + div( U ) = 0 252 // dU/dt + div( rho (u x u) + P I3 ) = 0 253 // dE/dt + div( (E + P) u ) = 0 254 // 255 // Equation of State: 256 // P = (gamma - 1) (E - rho (u u) / 2) 257 // 258 // Constants: 259 // cv , Specific heat, constant volume 260 // cp , Specific heat, constant pressure 261 // g , Gravity 262 // gamma = cp / cv, Specific heat ratio 263 // ***************************************************************************** 264 CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 265 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 266 const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 267 const CeedScalar(*q_data) = in[2]; 268 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 269 CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 270 271 EulerContext context = (EulerContext)ctx; 272 const CeedScalar c_tau = context->c_tau; 273 const CeedScalar gamma = 1.4; 274 275 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 276 // Setup 277 // -- Interp in 278 const CeedScalar rho = q[0][i]; 279 const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 280 const CeedScalar E = q[4][i]; 281 const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 282 const CeedScalar dU[3][3] = { 283 {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 284 {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 285 {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 286 }; 287 const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 288 CeedScalar wdetJ, dXdx[3][3]; 289 QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx); 290 // dU/dx 291 CeedScalar drhodx[3] = {0.}; 292 CeedScalar dEdx[3] = {0.}; 293 CeedScalar dUdx[3][3] = {{0.}}; 294 CeedScalar dXdxdXdxT[3][3] = {{0.}}; 295 for (CeedInt j = 0; j < 3; j++) { 296 for (CeedInt k = 0; k < 3; k++) { 297 drhodx[j] += drho[k] * dXdx[k][j]; 298 dEdx[j] += dE[k] * dXdx[k][j]; 299 for (CeedInt l = 0; l < 3; l++) { 300 dUdx[j][k] += dU[j][l] * dXdx[l][k]; 301 dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 302 } 303 } 304 } 305 // Pressure 306 const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic, 307 P = E_internal * (gamma - 1.); // P = pressure 308 309 // The Physics 310 // Zero v and dv so all future terms can safely sum into it 311 for (CeedInt j = 0; j < 5; j++) { 312 v[j][i] = 0.; 313 for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 314 } 315 316 // -- Density 317 // ---- u rho 318 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]); 319 // -- Momentum 320 // ---- rho (u x u) + P I3 321 for (CeedInt j = 0; j < 3; j++) { 322 for (CeedInt k = 0; k < 3; k++) { 323 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] + 324 (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 325 } 326 } 327 // -- Total Energy Density 328 // ---- (E + P) u 329 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]); 330 331 // --Stabilization terms 332 // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 333 CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 334 ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 335 336 // ---- dqdx collects drhodx, dUdx and dEdx in one vector 337 CeedScalar dqdx[5][3]; 338 for (CeedInt j = 0; j < 3; j++) { 339 dqdx[0][j] = drhodx[j]; 340 dqdx[4][j] = dEdx[j]; 341 for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 342 } 343 344 // ---- strong_conv = dF/dq * dq/dx (Strong convection) 345 CeedScalar strong_conv[5] = {0.}; 346 for (CeedInt j = 0; j < 3; j++) { 347 for (CeedInt k = 0; k < 5; k++) { 348 for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 349 } 350 } 351 352 // Stabilization 353 // -- Tau elements 354 const CeedScalar sound_speed = sqrt(gamma * P / rho); 355 CeedScalar Tau_x[3] = {0.}; 356 Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 357 358 // -- Stabilization method: none or SU 359 CeedScalar stab[5][3] = {{0.}}; 360 switch (context->stabilization) { 361 case 0: // Galerkin 362 break; 363 case 1: // SU 364 for (CeedInt j = 0; j < 3; j++) { 365 for (CeedInt k = 0; k < 5; k++) { 366 for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 367 } 368 } 369 370 for (CeedInt j = 0; j < 5; j++) { 371 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]); 372 } 373 break; 374 case 2: // SUPG is not implemented for explicit scheme 375 break; 376 } 377 } 378 return 0; 379 } 380 381 // ***************************************************************************** 382 // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method 383 // ***************************************************************************** 384 CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 385 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 386 const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 387 const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 388 const CeedScalar(*q_data) = in[3]; 389 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 390 CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 391 392 EulerContext context = (EulerContext)ctx; 393 const CeedScalar c_tau = context->c_tau; 394 const CeedScalar gamma = 1.4; 395 396 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 397 // Setup 398 // -- Interp in 399 const CeedScalar rho = q[0][i]; 400 const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 401 const CeedScalar E = q[4][i]; 402 const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 403 const CeedScalar dU[3][3] = { 404 {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 405 {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 406 {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 407 }; 408 const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 409 CeedScalar wdetJ, dXdx[3][3]; 410 QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx); 411 // dU/dx 412 CeedScalar drhodx[3] = {0.}; 413 CeedScalar dEdx[3] = {0.}; 414 CeedScalar dUdx[3][3] = {{0.}}; 415 CeedScalar dXdxdXdxT[3][3] = {{0.}}; 416 for (CeedInt j = 0; j < 3; j++) { 417 for (CeedInt k = 0; k < 3; k++) { 418 drhodx[j] += drho[k] * dXdx[k][j]; 419 dEdx[j] += dE[k] * dXdx[k][j]; 420 for (CeedInt l = 0; l < 3; l++) { 421 dUdx[j][k] += dU[j][l] * dXdx[l][k]; 422 dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 423 } 424 } 425 } 426 const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic, 427 P = E_internal * (gamma - 1.); // P = pressure 428 429 // The Physics 430 // Zero v and dv so all future terms can safely sum into it 431 for (CeedInt j = 0; j < 5; j++) { 432 v[j][i] = 0.; 433 for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 434 } 435 //-----mass matrix 436 for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i]; 437 438 // -- Density 439 // ---- u rho 440 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]); 441 // -- Momentum 442 // ---- rho (u x u) + P I3 443 for (CeedInt j = 0; j < 3; j++) { 444 for (CeedInt k = 0; k < 3; k++) { 445 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] + 446 (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 447 } 448 } 449 // -- Total Energy Density 450 // ---- (E + P) u 451 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]); 452 453 // -- Stabilization terms 454 // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 455 CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 456 ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 457 458 // ---- dqdx collects drhodx, dUdx and dEdx in one vector 459 CeedScalar dqdx[5][3]; 460 for (CeedInt j = 0; j < 3; j++) { 461 dqdx[0][j] = drhodx[j]; 462 dqdx[4][j] = dEdx[j]; 463 for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 464 } 465 466 // ---- strong_conv = dF/dq * dq/dx (Strong convection) 467 CeedScalar strong_conv[5] = {0.}; 468 for (CeedInt j = 0; j < 3; j++) { 469 for (CeedInt k = 0; k < 5; k++) { 470 for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 471 } 472 } 473 474 // ---- Strong residual 475 CeedScalar strong_res[5]; 476 for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j]; 477 478 // Stabilization 479 // -- Tau elements 480 const CeedScalar sound_speed = sqrt(gamma * P / rho); 481 CeedScalar Tau_x[3] = {0.}; 482 Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 483 484 // -- Stabilization method: none, SU, or SUPG 485 CeedScalar stab[5][3] = {{0.}}; 486 switch (context->stabilization) { 487 case 0: // Galerkin 488 break; 489 case 1: // SU 490 for (CeedInt j = 0; j < 3; j++) { 491 for (CeedInt k = 0; k < 5; k++) { 492 for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 493 } 494 } 495 496 for (CeedInt j = 0; j < 5; j++) { 497 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]); 498 } 499 break; 500 case 2: // SUPG 501 for (CeedInt j = 0; j < 3; j++) { 502 for (CeedInt k = 0; k < 5; k++) { 503 for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l]; 504 } 505 } 506 507 for (CeedInt j = 0; j < 5; j++) { 508 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]); 509 } 510 break; 511 } 512 } 513 return 0; 514 } 515 // ***************************************************************************** 516 // This QFunction sets the inflow boundary conditions for the traveling vortex problem. 517 // 518 // Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly. 519 // ***************************************************************************** 520 CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 521 const CeedScalar(*q_data_sur) = in[2]; 522 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 523 524 EulerContext context = (EulerContext)ctx; 525 const int euler_test = context->euler_test; 526 const bool is_implicit = context->implicit; 527 CeedScalar *mean_velocity = context->mean_velocity; 528 const CeedScalar cv = 2.5; 529 const CeedScalar R = 1.; 530 CeedScalar T_inlet; 531 CeedScalar P_inlet; 532 533 // For test cases 1 and 3 the background velocity is zero 534 if (euler_test == 1 || euler_test == 3) { 535 for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.; 536 } 537 538 // For test cases 1 and 2, T_inlet = T_inlet = 0.4 539 if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 540 else T_inlet = P_inlet = 1.; 541 542 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 543 CeedScalar wdetJb, normal[3]; 544 QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, normal); 545 wdetJb *= is_implicit ? -1. : 1.; 546 547 // face_normal = Normal vector of the face 548 const CeedScalar face_normal = Dot3(normal, mean_velocity); 549 // The Physics 550 // Zero v so all future terms can safely sum into it 551 for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 552 553 // Implementing in/outflow BCs 554 if (face_normal > 0) { 555 } else { // inflow 556 const CeedScalar rho_inlet = P_inlet / (R * T_inlet); 557 const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.; 558 // incoming total energy 559 const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 560 561 // The Physics 562 // -- Density 563 v[0][i] -= wdetJb * rho_inlet * face_normal; 564 565 // -- Momentum 566 for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho_inlet * face_normal * mean_velocity[j] + normal[j] * P_inlet); 567 568 // -- Total Energy Density 569 v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 570 } 571 } 572 return 0; 573 } 574 575 // ***************************************************************************** 576 // This QFunction sets the outflow boundary conditions for the Euler solver. 577 // 578 // Outflow BCs: 579 // The validity of the weak form of the governing equations is extended to the outflow. 580 // ***************************************************************************** 581 CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 582 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 583 const CeedScalar(*q_data_sur) = in[2]; 584 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 585 586 EulerContext context = (EulerContext)ctx; 587 const bool is_implicit = context->implicit; 588 CeedScalar *mean_velocity = context->mean_velocity; 589 590 const CeedScalar gamma = 1.4; 591 592 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 593 // Setup 594 // -- Interp in 595 const CeedScalar rho = q[0][i]; 596 const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 597 const CeedScalar E = q[4][i]; 598 599 CeedScalar wdetJb, normal[3]; 600 QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, normal); 601 wdetJb *= is_implicit ? -1. : 1.; 602 603 // face_normal = Normal vector of the face 604 const CeedScalar face_normal = Dot3(normal, mean_velocity); 605 // The Physics 606 // Zero v so all future terms can safely sum into it 607 for (CeedInt j = 0; j < 5; j++) v[j][i] = 0; 608 609 // Implementing in/outflow BCs 610 if (face_normal > 0) { // outflow 611 const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.; 612 const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 613 const CeedScalar u_normal = Dot3(normal, u); // Normal velocity 614 // The Physics 615 // -- Density 616 v[0][i] -= wdetJb * rho * u_normal; 617 618 // -- Momentum 619 for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + normal[j] * P); 620 621 // -- Total Energy Density 622 v[4][i] -= wdetJb * u_normal * (E + P); 623 } 624 } 625 return 0; 626 } 627