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