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