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 0; 174 } 175 176 // ***************************************************************************** 177 // Helper function for computing flux Jacobian 178 // ***************************************************************************** 179 CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 180 const CeedScalar gamma) { 181 CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2]; // Velocity square 182 for (CeedInt i = 0; i < 3; i++) { // Jacobian matrices for 3 directions 183 for (CeedInt j = 0; j < 3; j++) { // Rows of each Jacobian matrix 184 dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j]; 185 for (CeedInt k = 0; k < 3; k++) { // Columns of each Jacobian matrix 186 dF[i][0][k + 1] = ((i == k) ? 1. : 0.); 187 dF[i][j + 1][k + 1] = ((j == k) ? u[i] : 0.) + ((i == k) ? u[j] : 0.) - ((i == j) ? u[k] : 0.) * (gamma - 1.); 188 dF[i][4][k + 1] = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k]; 189 } 190 dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.); 191 } 192 dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho); 193 dF[i][4][4] = u[i] * gamma; 194 } 195 } 196 197 // ***************************************************************************** 198 // Helper function for computing Tau elements (stabilization constant) 199 // Model from: 200 // Stabilized Methods for Compressible Flows, Hughes et al 2010 201 // 202 // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 203 // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 204 // 205 // Where 206 // c_tau = stabilization constant (0.5 is reported as "optimal") 207 // h[i] = 2 length(dxdX[i]) 208 // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 209 // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 210 // rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i 211 // ***************************************************************************** 212 CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], const CeedScalar dXdx[3][3], const CeedScalar u[3], const CeedScalar sound_speed, 213 const CeedScalar c_tau) { 214 for (CeedInt i = 0; i < 3; i++) { 215 // length of element in direction i 216 CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]); 217 // fastest wave in direction i 218 CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 219 Tau_x[i] = c_tau * h / fastest_wave; 220 } 221 } 222 223 // ***************************************************************************** 224 // This QFunction sets the initial conditions for Euler traveling vortex 225 // ***************************************************************************** 226 CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 227 const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 228 CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 229 230 const EulerContext context = (EulerContext)ctx; 231 232 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 233 const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 234 CeedScalar q[5] = {0.}; 235 236 Exact_Euler(3, context->curr_time, x, 5, q, ctx); 237 for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 238 } 239 return 0; 240 } 241 242 // ***************************************************************************** 243 // This QFunction implements the following formulation of Euler equations with explicit time stepping method 244 // 245 // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density. 246 // 247 // State Variables: q = ( rho, U1, U2, U3, E ) 248 // rho - Mass Density 249 // Ui - Momentum Density, Ui = rho ui 250 // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 251 // 252 // Euler Equations: 253 // drho/dt + div( U ) = 0 254 // dU/dt + div( rho (u x u) + P I3 ) = 0 255 // dE/dt + div( (E + P) u ) = 0 256 // 257 // Equation of State: 258 // P = (gamma - 1) (E - rho (u u) / 2) 259 // 260 // Constants: 261 // cv , Specific heat, constant volume 262 // cp , Specific heat, constant pressure 263 // g , Gravity 264 // gamma = cp / cv, Specific heat ratio 265 // ***************************************************************************** 266 CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 267 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 268 const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 269 const CeedScalar(*q_data) = in[2]; 270 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 271 CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 272 273 EulerContext context = (EulerContext)ctx; 274 const CeedScalar c_tau = context->c_tau; 275 const CeedScalar gamma = 1.4; 276 277 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 278 // Setup 279 // -- Interp in 280 const CeedScalar rho = q[0][i]; 281 const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 282 const CeedScalar E = q[4][i]; 283 const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 284 const CeedScalar dU[3][3] = { 285 {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 286 {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 287 {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 288 }; 289 const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 290 CeedScalar wdetJ, dXdx[3][3]; 291 QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx); 292 // dU/dx 293 CeedScalar drhodx[3] = {0.}; 294 CeedScalar dEdx[3] = {0.}; 295 CeedScalar dUdx[3][3] = {{0.}}; 296 CeedScalar dXdxdXdxT[3][3] = {{0.}}; 297 for (CeedInt j = 0; j < 3; j++) { 298 for (CeedInt k = 0; k < 3; k++) { 299 drhodx[j] += drho[k] * dXdx[k][j]; 300 dEdx[j] += dE[k] * dXdx[k][j]; 301 for (CeedInt l = 0; l < 3; l++) { 302 dUdx[j][k] += dU[j][l] * dXdx[l][k]; 303 dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 304 } 305 } 306 } 307 // Pressure 308 const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic, 309 P = E_internal * (gamma - 1.); // P = pressure 310 311 // The Physics 312 // Zero v and dv so all future terms can safely sum into it 313 for (CeedInt j = 0; j < 5; j++) { 314 v[j][i] = 0.; 315 for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 316 } 317 318 // -- Density 319 // ---- u rho 320 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]); 321 // -- Momentum 322 // ---- rho (u x u) + P I3 323 for (CeedInt j = 0; j < 3; j++) { 324 for (CeedInt k = 0; k < 3; k++) { 325 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] + 326 (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 327 } 328 } 329 // -- Total Energy Density 330 // ---- (E + P) u 331 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]); 332 333 // --Stabilization terms 334 // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 335 CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 336 ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 337 338 // ---- dqdx collects drhodx, dUdx and dEdx in one vector 339 CeedScalar dqdx[5][3]; 340 for (CeedInt j = 0; j < 3; j++) { 341 dqdx[0][j] = drhodx[j]; 342 dqdx[4][j] = dEdx[j]; 343 for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 344 } 345 346 // ---- strong_conv = dF/dq * dq/dx (Strong convection) 347 CeedScalar strong_conv[5] = {0.}; 348 for (CeedInt j = 0; j < 3; j++) { 349 for (CeedInt k = 0; k < 5; k++) { 350 for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 351 } 352 } 353 354 // Stabilization 355 // -- Tau elements 356 const CeedScalar sound_speed = sqrt(gamma * P / rho); 357 CeedScalar Tau_x[3] = {0.}; 358 Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 359 360 // -- Stabilization method: none or SU 361 CeedScalar stab[5][3] = {{0.}}; 362 switch (context->stabilization) { 363 case 0: // Galerkin 364 break; 365 case 1: // SU 366 for (CeedInt j = 0; j < 3; j++) { 367 for (CeedInt k = 0; k < 5; k++) { 368 for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 369 } 370 } 371 372 for (CeedInt j = 0; j < 5; j++) { 373 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]); 374 } 375 break; 376 case 2: // SUPG is not implemented for explicit scheme 377 break; 378 } 379 } 380 return 0; 381 } 382 383 // ***************************************************************************** 384 // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method 385 // ***************************************************************************** 386 CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 387 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 388 const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 389 const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 390 const CeedScalar(*q_data) = in[3]; 391 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 392 CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 393 CeedScalar *jac_data = out[2]; 394 395 EulerContext context = (EulerContext)ctx; 396 const CeedScalar c_tau = context->c_tau; 397 const CeedScalar gamma = 1.4; 398 const CeedScalar zeros[14] = {0.}; 399 400 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 401 // Setup 402 // -- Interp in 403 const CeedScalar rho = q[0][i]; 404 const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 405 const CeedScalar E = q[4][i]; 406 const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 407 const CeedScalar dU[3][3] = { 408 {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 409 {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 410 {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 411 }; 412 const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 413 CeedScalar wdetJ, dXdx[3][3]; 414 QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx); 415 // dU/dx 416 CeedScalar drhodx[3] = {0.}; 417 CeedScalar dEdx[3] = {0.}; 418 CeedScalar dUdx[3][3] = {{0.}}; 419 CeedScalar dXdxdXdxT[3][3] = {{0.}}; 420 for (CeedInt j = 0; j < 3; j++) { 421 for (CeedInt k = 0; k < 3; k++) { 422 drhodx[j] += drho[k] * dXdx[k][j]; 423 dEdx[j] += dE[k] * dXdx[k][j]; 424 for (CeedInt l = 0; l < 3; l++) { 425 dUdx[j][k] += dU[j][l] * dXdx[l][k]; 426 dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 427 } 428 } 429 } 430 const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic, 431 P = E_internal * (gamma - 1.); // P = pressure 432 433 // The Physics 434 // Zero v and dv so all future terms can safely sum into it 435 for (CeedInt j = 0; j < 5; j++) { 436 v[j][i] = 0.; 437 for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 438 } 439 //-----mass matrix 440 for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i]; 441 442 // -- Density 443 // ---- u rho 444 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]); 445 // -- Momentum 446 // ---- rho (u x u) + P I3 447 for (CeedInt j = 0; j < 3; j++) { 448 for (CeedInt k = 0; k < 3; k++) { 449 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] + 450 (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 451 } 452 } 453 // -- Total Energy Density 454 // ---- (E + P) u 455 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]); 456 457 // -- Stabilization terms 458 // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 459 CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 460 ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 461 462 // ---- dqdx collects drhodx, dUdx and dEdx in one vector 463 CeedScalar dqdx[5][3]; 464 for (CeedInt j = 0; j < 3; j++) { 465 dqdx[0][j] = drhodx[j]; 466 dqdx[4][j] = dEdx[j]; 467 for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 468 } 469 470 // ---- strong_conv = dF/dq * dq/dx (Strong convection) 471 CeedScalar strong_conv[5] = {0.}; 472 for (CeedInt j = 0; j < 3; j++) { 473 for (CeedInt k = 0; k < 5; k++) { 474 for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 475 } 476 } 477 478 // ---- Strong residual 479 CeedScalar strong_res[5]; 480 for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j]; 481 482 // Stabilization 483 // -- Tau elements 484 const CeedScalar sound_speed = sqrt(gamma * P / rho); 485 CeedScalar Tau_x[3] = {0.}; 486 Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 487 488 // -- Stabilization method: none, SU, or SUPG 489 CeedScalar stab[5][3] = {{0.}}; 490 switch (context->stabilization) { 491 case 0: // Galerkin 492 break; 493 case 1: // SU 494 for (CeedInt j = 0; j < 3; j++) { 495 for (CeedInt k = 0; k < 5; k++) { 496 for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 497 } 498 } 499 500 for (CeedInt j = 0; j < 5; j++) { 501 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]); 502 } 503 break; 504 case 2: // SUPG 505 for (CeedInt j = 0; j < 3; j++) { 506 for (CeedInt k = 0; k < 5; k++) { 507 for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l]; 508 } 509 } 510 511 for (CeedInt j = 0; j < 5; j++) { 512 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]); 513 } 514 break; 515 } 516 StoredValuesPack(Q, i, 0, 14, zeros, jac_data); 517 } 518 return 0; 519 } 520 // ***************************************************************************** 521 // This QFunction sets the inflow boundary conditions for the traveling vortex problem. 522 // 523 // Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly. 524 // ***************************************************************************** 525 CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 526 const CeedScalar(*q_data_sur) = in[2]; 527 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 528 529 EulerContext context = (EulerContext)ctx; 530 const int euler_test = context->euler_test; 531 const bool is_implicit = context->implicit; 532 CeedScalar *mean_velocity = context->mean_velocity; 533 const CeedScalar cv = 2.5; 534 const CeedScalar R = 1.; 535 CeedScalar T_inlet; 536 CeedScalar P_inlet; 537 538 // For test cases 1 and 3 the background velocity is zero 539 if (euler_test == 1 || euler_test == 3) { 540 for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.; 541 } 542 543 // For test cases 1 and 2, T_inlet = T_inlet = 0.4 544 if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 545 else T_inlet = P_inlet = 1.; 546 547 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 548 CeedScalar wdetJb, norm[3]; 549 QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm); 550 wdetJb *= is_implicit ? -1. : 1.; 551 552 // face_normal = Normal vector of the face 553 const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2]; 554 // The Physics 555 // Zero v so all future terms can safely sum into it 556 for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 557 558 // Implementing in/outflow BCs 559 if (face_normal > 0) { 560 } else { // inflow 561 const CeedScalar rho_inlet = P_inlet / (R * T_inlet); 562 const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.; 563 // incoming total energy 564 const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 565 566 // The Physics 567 // -- Density 568 v[0][i] -= wdetJb * rho_inlet * face_normal; 569 570 // -- Momentum 571 for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho_inlet * face_normal * mean_velocity[j] + norm[j] * P_inlet); 572 573 // -- Total Energy Density 574 v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 575 } 576 } 577 return 0; 578 } 579 580 // ***************************************************************************** 581 // This QFunction sets the outflow boundary conditions for the Euler solver. 582 // 583 // Outflow BCs: 584 // The validity of the weak form of the governing equations is extended to the outflow. 585 // ***************************************************************************** 586 CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 587 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 588 const CeedScalar(*q_data_sur) = in[2]; 589 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 590 591 EulerContext context = (EulerContext)ctx; 592 const bool is_implicit = context->implicit; 593 CeedScalar *mean_velocity = context->mean_velocity; 594 595 const CeedScalar gamma = 1.4; 596 597 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 598 // Setup 599 // -- Interp in 600 const CeedScalar rho = q[0][i]; 601 const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 602 const CeedScalar E = q[4][i]; 603 604 CeedScalar wdetJb, norm[3]; 605 QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm); 606 wdetJb *= is_implicit ? -1. : 1.; 607 608 // face_normal = Normal vector of the face 609 const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2]; 610 // The Physics 611 // Zero v so all future terms can safely sum into it 612 for (CeedInt j = 0; j < 5; j++) v[j][i] = 0; 613 614 // Implementing in/outflow BCs 615 if (face_normal > 0) { // outflow 616 const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.; 617 const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 618 const CeedScalar u_normal = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2]; // Normal velocity 619 // The Physics 620 // -- Density 621 v[0][i] -= wdetJb * rho * u_normal; 622 623 // -- Momentum 624 for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 625 626 // -- Total Energy Density 627 v[4][i] -= wdetJb * u_normal * (E + P); 628 } 629 } 630 return 0; 631 } 632