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