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