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 /// Implementation of the Synthetic Turbulence Generation (STG) algorithm 10 /// presented in Shur et al. 2014 11 // 12 /// SetupSTG_Rand reads in the input files and fills in STGShur14Context. Then 13 /// STGShur14_CalcQF is run over quadrature points. Before the program exits, 14 /// TearDownSTG is run to free the memory of the allocated arrays. 15 16 #ifndef stg_shur14_h 17 #define stg_shur14_h 18 19 #include <ceed.h> 20 #include <math.h> 21 #include <stdlib.h> 22 23 #include "stg_shur14_type.h" 24 #include "utils.h" 25 26 #define STG_NMODES_MAX 1024 27 28 /* 29 * @brief Interpolate quantities from input profile to given location 30 * 31 * Assumed that prof_wd[i+1] > prof_wd[i] and prof_wd[0] = 0 32 * If wall_dist > prof_wd[-1], then the interpolation takes the values at prof_wd[-1] 33 * 34 * @param[in] wall_dist Distance to the nearest wall 35 * @param[out] ubar Mean velocity at wall_dist 36 * @param[out] cij Cholesky decomposition at wall_dist 37 * @param[out] eps Turbulent dissipation at wall_dist 38 * @param[out] lt Turbulent length scale at wall_dist 39 * @param[in] stg_ctx STGShur14Context for the problem 40 */ 41 CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar wall_dist, CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt, 42 const STGShur14Context stg_ctx) { 43 const CeedInt nprofs = stg_ctx->nprofs; 44 const CeedScalar *prof_wd = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 45 const CeedScalar *prof_eps = &stg_ctx->data[stg_ctx->offsets.eps]; 46 const CeedScalar *prof_lt = &stg_ctx->data[stg_ctx->offsets.lt]; 47 const CeedScalar *prof_ubar = &stg_ctx->data[stg_ctx->offsets.ubar]; 48 const CeedScalar *prof_cij = &stg_ctx->data[stg_ctx->offsets.cij]; 49 CeedInt idx = -1; 50 51 for (CeedInt i = 0; i < nprofs; i++) { 52 if (wall_dist < prof_wd[i]) { 53 idx = i; 54 break; 55 } 56 } 57 58 if (idx > 0) { // y within the bounds of prof_wd 59 CeedScalar coeff = (wall_dist - prof_wd[idx - 1]) / (prof_wd[idx] - prof_wd[idx - 1]); 60 61 ubar[0] = prof_ubar[0 * nprofs + idx - 1] + coeff * (prof_ubar[0 * nprofs + idx] - prof_ubar[0 * nprofs + idx - 1]); 62 ubar[1] = prof_ubar[1 * nprofs + idx - 1] + coeff * (prof_ubar[1 * nprofs + idx] - prof_ubar[1 * nprofs + idx - 1]); 63 ubar[2] = prof_ubar[2 * nprofs + idx - 1] + coeff * (prof_ubar[2 * nprofs + idx] - prof_ubar[2 * nprofs + idx - 1]); 64 cij[0] = prof_cij[0 * nprofs + idx - 1] + coeff * (prof_cij[0 * nprofs + idx] - prof_cij[0 * nprofs + idx - 1]); 65 cij[1] = prof_cij[1 * nprofs + idx - 1] + coeff * (prof_cij[1 * nprofs + idx] - prof_cij[1 * nprofs + idx - 1]); 66 cij[2] = prof_cij[2 * nprofs + idx - 1] + coeff * (prof_cij[2 * nprofs + idx] - prof_cij[2 * nprofs + idx - 1]); 67 cij[3] = prof_cij[3 * nprofs + idx - 1] + coeff * (prof_cij[3 * nprofs + idx] - prof_cij[3 * nprofs + idx - 1]); 68 cij[4] = prof_cij[4 * nprofs + idx - 1] + coeff * (prof_cij[4 * nprofs + idx] - prof_cij[4 * nprofs + idx - 1]); 69 cij[5] = prof_cij[5 * nprofs + idx - 1] + coeff * (prof_cij[5 * nprofs + idx] - prof_cij[5 * nprofs + idx - 1]); 70 *eps = prof_eps[idx - 1] + coeff * (prof_eps[idx] - prof_eps[idx - 1]); 71 *lt = prof_lt[idx - 1] + coeff * (prof_lt[idx] - prof_lt[idx - 1]); 72 } else { // y outside bounds of prof_wd 73 ubar[0] = prof_ubar[1 * nprofs - 1]; 74 ubar[1] = prof_ubar[2 * nprofs - 1]; 75 ubar[2] = prof_ubar[3 * nprofs - 1]; 76 cij[0] = prof_cij[1 * nprofs - 1]; 77 cij[1] = prof_cij[2 * nprofs - 1]; 78 cij[2] = prof_cij[3 * nprofs - 1]; 79 cij[3] = prof_cij[4 * nprofs - 1]; 80 cij[4] = prof_cij[5 * nprofs - 1]; 81 cij[5] = prof_cij[6 * nprofs - 1]; 82 *eps = prof_eps[nprofs - 1]; 83 *lt = prof_lt[nprofs - 1]; 84 } 85 } 86 87 /* 88 * @brief Calculate spectrum coefficient, qn 89 * 90 * Calculates q_n at a given distance to the wall 91 * 92 * @param[in] kappa nth wavenumber 93 * @param[in] dkappa Difference between wavenumbers 94 * @param[in] keta Dissipation wavenumber 95 * @param[in] kcut Mesh-induced cutoff wavenumber 96 * @param[in] ke Energy-containing wavenumber 97 * @param[in] Ektot Total turbulent kinetic energy of spectrum 98 * @returns qn Spectrum coefficient 99 */ 100 CEED_QFUNCTION_HELPER CeedScalar Calc_qn(const CeedScalar kappa, const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut, 101 const CeedScalar ke, const CeedScalar Ektot_inv) { 102 const CeedScalar feta_x_fcut = exp(-Square(12 * kappa / keta) - Cube(4 * Max(kappa - 0.9 * kcut, 0) / kcut)); 103 return pow(kappa / ke, 4.) * pow(1 + 2.4 * Square(kappa / ke), -17. / 6) * feta_x_fcut * dkappa * Ektot_inv; 104 } 105 106 // Calculate hmax, ke, keta, and kcut 107 CEED_QFUNCTION_HELPER void SpectrumConstants(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 108 const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke, CeedScalar *keta, CeedScalar *kcut) { 109 *hmax = Max(Max(h[0], h[1]), h[2]); 110 *ke = wall_dist == 0 ? 1e16 : 2 * M_PI / Min(2 * wall_dist, 3 * lt); 111 *keta = 2 * M_PI * pow(Cube(nu) / eps, -0.25); 112 *kcut = M_PI / Min(Max(Max(h[1], h[2]), 0.3 * (*hmax)) + 0.1 * wall_dist, *hmax); 113 } 114 115 /* 116 * @brief Calculate spectrum coefficients for STG 117 * 118 * Calculates q_n at a given distance to the wall 119 * 120 * @param[in] wall_dist Distance to the nearest wall 121 * @param[in] eps Turbulent dissipation w/rt wall_dist 122 * @param[in] lt Turbulent length scale w/rt wall_dist 123 * @param[in] h Element lengths in coordinate directions 124 * @param[in] nu Dynamic Viscosity; 125 * @param[in] stg_ctx STGShur14Context for the problem 126 * @param[out] qn Spectrum coefficients, [nmodes] 127 */ 128 CEED_QFUNCTION_HELPER void CalcSpectrum(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 129 const CeedScalar nu, CeedScalar qn[], const STGShur14Context stg_ctx) { 130 const CeedInt nmodes = stg_ctx->nmodes; 131 const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 132 CeedScalar hmax, ke, keta, kcut, Ektot = 0.0; 133 134 SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 135 136 for (CeedInt n = 0; n < nmodes; n++) { 137 const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 138 qn[n] = Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 139 Ektot += qn[n]; 140 } 141 142 if (Ektot == 0) return; 143 for (CeedInt n = 0; n < nmodes; n++) qn[n] /= Ektot; 144 } 145 146 /****************************************************** 147 * @brief Calculate u(x,t) for STG inflow condition 148 * 149 * @param[in] X Location to evaluate u(X,t) 150 * @param[in] t Time to evaluate u(X,t) 151 * @param[in] ubar Mean velocity at X 152 * @param[in] cij Cholesky decomposition at X 153 * @param[in] qn Wavemode amplitudes at X, [nmodes] 154 * @param[out] u Velocity at X and t 155 * @param[in] stg_ctx STGShur14Context for the problem 156 */ 157 CEED_QFUNCTION_HELPER void STGShur14_Calc(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 158 const CeedScalar qn[], CeedScalar u[3], const STGShur14Context stg_ctx) { 159 const CeedInt nmodes = stg_ctx->nmodes; 160 const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 161 const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 162 const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 163 const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 164 CeedScalar xdotd, vp[3] = {0.}; 165 CeedScalar xhat[] = {0., X[1], X[2]}; 166 167 CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 168 xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 169 xdotd = 0.; 170 for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 171 const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 172 vp[0] += sqrt(qn[n]) * sigma[0 * nmodes + n] * cos_kxdp; 173 vp[1] += sqrt(qn[n]) * sigma[1 * nmodes + n] * cos_kxdp; 174 vp[2] += sqrt(qn[n]) * sigma[2 * nmodes + n] * cos_kxdp; 175 } 176 for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 177 178 u[0] = ubar[0] + cij[0] * vp[0]; 179 u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 180 u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 181 } 182 183 /****************************************************** 184 * @brief Calculate u(x,t) for STG inflow condition 185 * 186 * @param[in] X Location to evaluate u(X,t) 187 * @param[in] t Time to evaluate u(X,t) 188 * @param[in] ubar Mean velocity at X 189 * @param[in] cij Cholesky decomposition at X 190 * @param[in] Ektot Total spectrum energy at this location 191 * @param[in] h Element size in 3 directions 192 * @param[in] wall_dist Distance to closest wall 193 * @param[in] eps Turbulent dissipation 194 * @param[in] lt Turbulent length scale 195 * @param[out] u Velocity at X and t 196 * @param[in] stg_ctx STGShur14Context for the problem 197 */ 198 CEED_QFUNCTION_HELPER void STGShur14_Calc_PrecompEktot(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 199 const CeedScalar Ektot, const CeedScalar h[3], const CeedScalar wall_dist, 200 const CeedScalar eps, const CeedScalar lt, const CeedScalar nu, CeedScalar u[3], 201 const STGShur14Context stg_ctx) { 202 const CeedInt nmodes = stg_ctx->nmodes; 203 const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 204 const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 205 const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 206 const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 207 CeedScalar hmax, ke, keta, kcut; 208 SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 209 CeedScalar xdotd, vp[3] = {0.}; 210 CeedScalar xhat[] = {0., X[1], X[2]}; 211 212 CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 213 xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 214 xdotd = 0.; 215 for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 216 const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 217 const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 218 const CeedScalar qn = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot); 219 vp[0] += sqrt(qn) * sigma[0 * nmodes + n] * cos_kxdp; 220 vp[1] += sqrt(qn) * sigma[1 * nmodes + n] * cos_kxdp; 221 vp[2] += sqrt(qn) * sigma[2 * nmodes + n] * cos_kxdp; 222 } 223 for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 224 225 u[0] = ubar[0] + cij[0] * vp[0]; 226 u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 227 u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 228 } 229 230 // Create preprocessed input for the stg calculation 231 // 232 // stg_data[0] = 1 / Ektot (inverse of total spectrum energy) 233 CEED_QFUNCTION(Preprocess_STGShur14)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 234 const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 235 236 CeedScalar(*stg_data) = (CeedScalar(*))out[0]; 237 238 CeedScalar ubar[3], cij[6], eps, lt; 239 const STGShur14Context stg_ctx = (STGShur14Context)ctx; 240 const CeedScalar dx = stg_ctx->dx; 241 const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 242 const CeedScalar theta0 = stg_ctx->theta0; 243 const CeedScalar P0 = stg_ctx->P0; 244 const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 245 const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 246 const CeedScalar Rd = cp - cv; 247 const CeedScalar rho = P0 / (Rd * theta0); 248 const CeedScalar nu = mu / rho; 249 250 const CeedInt nmodes = stg_ctx->nmodes; 251 const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 252 CeedScalar hmax, ke, keta, kcut; 253 254 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 255 const CeedScalar wall_dist = x[1][i]; 256 const CeedScalar dXdx[2][3] = { 257 {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 258 {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 259 }; 260 261 CeedScalar h[3]; 262 h[0] = dx; 263 for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(dXdx[0][j] * dXdx[0][j] + dXdx[1][j] * dXdx[1][j]); 264 265 InterpolateProfile(wall_dist, ubar, cij, &eps, <, stg_ctx); 266 SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 267 268 // Calculate total TKE per spectrum 269 CeedScalar Ek_tot = 0; 270 CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 271 const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 272 Ek_tot += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 273 } 274 // avoid underflowed and poorly defined spectrum coefficients 275 stg_data[i] = Ek_tot != 0 ? 1 / Ek_tot : 0; 276 } 277 return 0; 278 } 279 280 // Extrude the STGInflow profile through out the domain for an initial condition 281 CEED_QFUNCTION(ICsSTG)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 282 // Inputs 283 const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 284 // Outputs 285 CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 286 287 const STGShur14Context stg_ctx = (STGShur14Context)ctx; 288 CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 289 const CeedScalar dx = stg_ctx->dx; 290 const CeedScalar time = stg_ctx->time; 291 const CeedScalar theta0 = stg_ctx->theta0; 292 const CeedScalar P0 = stg_ctx->P0; 293 const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 294 const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 295 const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 296 const CeedScalar Rd = cp - cv; 297 const CeedScalar rho = P0 / (Rd * theta0); 298 299 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 300 const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 301 const CeedScalar dXdx[3][3] = { 302 {q_data[1][i], q_data[2][i], q_data[3][i]}, 303 {q_data[4][i], q_data[5][i], q_data[6][i]}, 304 {q_data[7][i], q_data[8][i], q_data[9][i]} 305 }; 306 307 CeedScalar h[3]; 308 h[0] = dx; 309 for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]) + Square(dXdx[2][j])); 310 311 InterpolateProfile(x_i[1], ubar, cij, &eps, <, stg_ctx); 312 if (stg_ctx->use_fluctuating_IC) { 313 CalcSpectrum(x_i[1], eps, lt, h, mu / rho, qn, stg_ctx); 314 STGShur14_Calc(x_i, time, ubar, cij, qn, u, stg_ctx); 315 } else { 316 for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 317 } 318 319 switch (stg_ctx->newtonian_ctx.state_var) { 320 case STATEVAR_CONSERVATIVE: 321 q0[0][i] = rho; 322 q0[1][i] = u[0] * rho; 323 q0[2][i] = u[1] * rho; 324 q0[3][i] = u[2] * rho; 325 q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0); 326 break; 327 328 case STATEVAR_PRIMITIVE: 329 q0[0][i] = P0; 330 q0[1][i] = u[0]; 331 q0[2][i] = u[1]; 332 q0[3][i] = u[2]; 333 q0[4][i] = theta0; 334 break; 335 } 336 } // End of Quadrature Point Loop 337 return 0; 338 } 339 340 /******************************************************************** 341 * @brief QFunction to calculate the inflow boundary condition 342 * 343 * This will loop through quadrature points, calculate the wavemode amplitudes 344 * at each location, then calculate the actual velocity. 345 */ 346 CEED_QFUNCTION(STGShur14_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 347 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 348 (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 349 350 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], (*jac_data_sur)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[1]; 351 352 const STGShur14Context stg_ctx = (STGShur14Context)ctx; 353 CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 354 const bool is_implicit = stg_ctx->is_implicit; 355 const bool mean_only = stg_ctx->mean_only; 356 const bool prescribe_T = stg_ctx->prescribe_T; 357 const CeedScalar dx = stg_ctx->dx; 358 const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 359 const CeedScalar time = stg_ctx->time; 360 const CeedScalar theta0 = stg_ctx->theta0; 361 const CeedScalar P0 = stg_ctx->P0; 362 const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 363 const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 364 const CeedScalar Rd = cp - cv; 365 const CeedScalar gamma = cp / cv; 366 367 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 368 const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0); 369 const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 370 const CeedScalar dXdx[2][3] = { 371 {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 372 {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 373 }; 374 375 CeedScalar h[3]; 376 h[0] = dx; 377 for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 378 379 InterpolateProfile(X[1][i], ubar, cij, &eps, <, stg_ctx); 380 if (!mean_only) { 381 CalcSpectrum(X[1][i], eps, lt, h, mu / rho, qn, stg_ctx); 382 STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 383 } else { 384 for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 385 } 386 387 const CeedScalar E_kinetic = .5 * rho * Dot3(u, u); 388 CeedScalar E_internal, P; 389 if (prescribe_T) { 390 // Temperature is being set weakly (theta0) and for constant cv this sets E_internal 391 E_internal = rho * cv * theta0; 392 // Find pressure using 393 P = rho * Rd * theta0; // interior rho with exterior T 394 } else { 395 E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution 396 P = E_internal * (gamma - 1.); 397 } 398 399 const CeedScalar wdetJb = (is_implicit ? -1. : 1.) * q_data_sur[0][i]; 400 // ---- Normal vect 401 const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 402 403 const CeedScalar E = E_internal + E_kinetic; 404 405 // Velocity normal to the boundary 406 const CeedScalar u_normal = Dot3(norm, u); 407 408 // The Physics 409 // Zero v so all future terms can safely sum into it 410 for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 411 412 // The Physics 413 // -- Density 414 v[0][i] -= wdetJb * rho * u_normal; 415 416 // -- Momentum 417 for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 418 419 // -- Total Energy Density 420 v[4][i] -= wdetJb * u_normal * (E + P); 421 422 jac_data_sur[0][i] = rho; 423 jac_data_sur[1][i] = u[0]; 424 jac_data_sur[2][i] = u[1]; 425 jac_data_sur[3][i] = u[2]; 426 jac_data_sur[4][i] = E; 427 for (int j = 0; j < 6; j++) jac_data_sur[5 + j][i] = 0.; 428 } 429 return 0; 430 } 431 432 CEED_QFUNCTION(STGShur14_Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 433 // Inputs 434 const CeedScalar(*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 435 (*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 436 // Outputs 437 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 438 const STGShur14Context stg_ctx = (STGShur14Context)ctx; 439 const bool implicit = stg_ctx->is_implicit; 440 const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 441 const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 442 const CeedScalar Rd = cp - cv; 443 const CeedScalar gamma = cp / cv; 444 445 const CeedScalar theta0 = stg_ctx->theta0; 446 const bool prescribe_T = stg_ctx->prescribe_T; 447 448 CeedPragmaSIMD 449 // Quadrature Point Loop 450 for (CeedInt i = 0; i < Q; i++) { 451 // Setup 452 // -- Interp-to-Interp q_data 453 // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 454 // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 455 // We can effect this by swapping the sign on this weight 456 const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 457 458 // Calculate inflow values 459 CeedScalar velocity[3]; 460 for (CeedInt j = 0; j < 3; j++) velocity[j] = jac_data_sur[5 + j][i]; 461 462 // enabling user to choose between weak T and weak rho inflow 463 CeedScalar drho, dE, dP; 464 if (prescribe_T) { 465 // rho should be from the current solution 466 drho = dq[0][i]; 467 CeedScalar dE_internal = drho * cv * theta0; 468 CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity); 469 dE = dE_internal + dE_kinetic; 470 dP = drho * Rd * theta0; // interior rho with exterior T 471 } else { // rho specified, E_internal from solution 472 drho = 0; 473 dE = dq[4][i]; 474 dP = dE * (gamma - 1.); 475 } 476 const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 477 478 const CeedScalar u_normal = Dot3(norm, velocity); 479 480 v[0][i] = -wdetJb * drho * u_normal; 481 for (int j = 0; j < 3; j++) v[j + 1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP); 482 v[4][i] = -wdetJb * u_normal * (dE + dP); 483 } // End Quadrature Point Loop 484 return 0; 485 } 486 487 /******************************************************************** 488 * @brief QFunction to calculate the strongly enforce inflow BC 489 * 490 * This QF is for the strong application of STG via libCEED (rather than 491 * through the native PETSc `DMAddBoundary` -> `bcFunc` method. 492 */ 493 CEED_QFUNCTION(STGShur14_Inflow_StrongQF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 494 const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*coords)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1], 495 (*scale) = (const CeedScalar(*))in[2], (*stg_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 496 497 CeedScalar(*bcval)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 498 499 const STGShur14Context stg_ctx = (STGShur14Context)ctx; 500 CeedScalar u[3], ubar[3], cij[6], eps, lt; 501 const bool mean_only = stg_ctx->mean_only; 502 const CeedScalar dx = stg_ctx->dx; 503 const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 504 const CeedScalar time = stg_ctx->time; 505 const CeedScalar theta0 = stg_ctx->theta0; 506 const CeedScalar P0 = stg_ctx->P0; 507 const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 508 const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 509 const CeedScalar Rd = cp - cv; 510 const CeedScalar rho = P0 / (Rd * theta0); 511 512 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 513 const CeedScalar x[] = {coords[0][i], coords[1][i], coords[2][i]}; 514 const CeedScalar dXdx[2][3] = { 515 {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 516 {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 517 }; 518 519 CeedScalar h[3]; 520 h[0] = dx; 521 for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 522 523 InterpolateProfile(coords[1][i], ubar, cij, &eps, <, stg_ctx); 524 if (!mean_only) { 525 if (1) { 526 STGShur14_Calc_PrecompEktot(x, time, ubar, cij, stg_data[0][i], h, x[1], eps, lt, mu / rho, u, stg_ctx); 527 } else { // Original way 528 CeedScalar qn[STG_NMODES_MAX]; 529 CalcSpectrum(coords[1][i], eps, lt, h, mu / rho, qn, stg_ctx); 530 STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 531 } 532 } else { 533 for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 534 } 535 536 switch (stg_ctx->newtonian_ctx.state_var) { 537 case STATEVAR_CONSERVATIVE: 538 bcval[0][i] = scale[i] * rho; 539 bcval[1][i] = scale[i] * rho * u[0]; 540 bcval[2][i] = scale[i] * rho * u[1]; 541 bcval[3][i] = scale[i] * rho * u[2]; 542 bcval[4][i] = 0.; 543 break; 544 545 case STATEVAR_PRIMITIVE: 546 bcval[0][i] = 0; 547 bcval[1][i] = scale[i] * u[0]; 548 bcval[2][i] = scale[i] * u[1]; 549 bcval[3][i] = scale[i] * u[2]; 550 bcval[4][i] = scale[i] * theta0; 551 break; 552 } 553 } 554 return 0; 555 } 556 557 #endif // stg_shur14_h 558