1ba6664aeSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2ba6664aeSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3ba6664aeSJames Wright // 4ba6664aeSJames Wright // SPDX-License-Identifier: BSD-2-Clause 5ba6664aeSJames Wright // 6ba6664aeSJames Wright // This file is part of CEED: http://github.com/ceed 7ba6664aeSJames Wright 8ba6664aeSJames Wright /// @file 9ba6664aeSJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm 10ba6664aeSJames Wright /// presented in Shur et al. 2014 11ba6664aeSJames Wright // 12ba6664aeSJames Wright /// SetupSTG_Rand reads in the input files and fills in STGShur14Context. Then 13ba6664aeSJames Wright /// STGShur14_CalcQF is run over quadrature points. Before the program exits, 14ba6664aeSJames Wright /// TearDownSTG is run to free the memory of the allocated arrays. 15ba6664aeSJames Wright 16ba6664aeSJames Wright #ifndef stg_shur14_h 17ba6664aeSJames Wright #define stg_shur14_h 18ba6664aeSJames Wright 19ba6664aeSJames Wright #include <ceed.h> 20c9c2c079SJeremy L Thompson #include <math.h> 21ba6664aeSJames Wright #include <stdlib.h> 222b730f8bSJeremy L Thompson 23*46603fc5SJames Wright #include "newtonian_state.h" 24ba6664aeSJames Wright #include "stg_shur14_type.h" 2513fa47b2SJames Wright #include "utils.h" 26ba6664aeSJames Wright 27ba6664aeSJames Wright #define STG_NMODES_MAX 1024 28ba6664aeSJames Wright 29ba6664aeSJames Wright /* 30ba6664aeSJames Wright * @brief Interpolate quantities from input profile to given location 31ba6664aeSJames Wright * 32175f00a6SJames Wright * Assumed that prof_wd[i+1] > prof_wd[i] and prof_wd[0] = 0 33175f00a6SJames Wright * If wall_dist > prof_wd[-1], then the interpolation takes the values at prof_wd[-1] 34ba6664aeSJames Wright * 35175f00a6SJames Wright * @param[in] wall_dist Distance to the nearest wall 36175f00a6SJames Wright * @param[out] ubar Mean velocity at wall_dist 37175f00a6SJames Wright * @param[out] cij Cholesky decomposition at wall_dist 38175f00a6SJames Wright * @param[out] eps Turbulent dissipation at wall_dist 39175f00a6SJames Wright * @param[out] lt Turbulent length scale at wall_dist 40ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 41ba6664aeSJames Wright */ 422b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar wall_dist, CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt, 43ba6664aeSJames Wright const STGShur14Context stg_ctx) { 44ba6664aeSJames Wright const CeedInt nprofs = stg_ctx->nprofs; 45175f00a6SJames Wright const CeedScalar *prof_wd = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 46ba6664aeSJames Wright const CeedScalar *prof_eps = &stg_ctx->data[stg_ctx->offsets.eps]; 47ba6664aeSJames Wright const CeedScalar *prof_lt = &stg_ctx->data[stg_ctx->offsets.lt]; 48ba6664aeSJames Wright const CeedScalar *prof_ubar = &stg_ctx->data[stg_ctx->offsets.ubar]; 49ba6664aeSJames Wright const CeedScalar *prof_cij = &stg_ctx->data[stg_ctx->offsets.cij]; 50ba6664aeSJames Wright CeedInt idx = -1; 51ba6664aeSJames Wright 52ba6664aeSJames Wright for (CeedInt i = 0; i < nprofs; i++) { 53175f00a6SJames Wright if (wall_dist < prof_wd[i]) { 54ba6664aeSJames Wright idx = i; 55ba6664aeSJames Wright break; 56ba6664aeSJames Wright } 57ba6664aeSJames Wright } 58ba6664aeSJames Wright 59175f00a6SJames Wright if (idx > 0) { // y within the bounds of prof_wd 60175f00a6SJames Wright CeedScalar coeff = (wall_dist - prof_wd[idx - 1]) / (prof_wd[idx] - prof_wd[idx - 1]); 61175f00a6SJames Wright 62ba6664aeSJames Wright ubar[0] = prof_ubar[0 * nprofs + idx - 1] + coeff * (prof_ubar[0 * nprofs + idx] - prof_ubar[0 * nprofs + idx - 1]); 63ba6664aeSJames Wright ubar[1] = prof_ubar[1 * nprofs + idx - 1] + coeff * (prof_ubar[1 * nprofs + idx] - prof_ubar[1 * nprofs + idx - 1]); 64ba6664aeSJames Wright ubar[2] = prof_ubar[2 * nprofs + idx - 1] + coeff * (prof_ubar[2 * nprofs + idx] - prof_ubar[2 * nprofs + idx - 1]); 65ba6664aeSJames Wright cij[0] = prof_cij[0 * nprofs + idx - 1] + coeff * (prof_cij[0 * nprofs + idx] - prof_cij[0 * nprofs + idx - 1]); 66ba6664aeSJames Wright cij[1] = prof_cij[1 * nprofs + idx - 1] + coeff * (prof_cij[1 * nprofs + idx] - prof_cij[1 * nprofs + idx - 1]); 67ba6664aeSJames Wright cij[2] = prof_cij[2 * nprofs + idx - 1] + coeff * (prof_cij[2 * nprofs + idx] - prof_cij[2 * nprofs + idx - 1]); 68ba6664aeSJames Wright cij[3] = prof_cij[3 * nprofs + idx - 1] + coeff * (prof_cij[3 * nprofs + idx] - prof_cij[3 * nprofs + idx - 1]); 69ba6664aeSJames Wright cij[4] = prof_cij[4 * nprofs + idx - 1] + coeff * (prof_cij[4 * nprofs + idx] - prof_cij[4 * nprofs + idx - 1]); 70ba6664aeSJames Wright cij[5] = prof_cij[5 * nprofs + idx - 1] + coeff * (prof_cij[5 * nprofs + idx] - prof_cij[5 * nprofs + idx - 1]); 71ba6664aeSJames Wright *eps = prof_eps[idx - 1] + coeff * (prof_eps[idx] - prof_eps[idx - 1]); 72ba6664aeSJames Wright *lt = prof_lt[idx - 1] + coeff * (prof_lt[idx] - prof_lt[idx - 1]); 73175f00a6SJames Wright } else { // y outside bounds of prof_wd 74ba6664aeSJames Wright ubar[0] = prof_ubar[1 * nprofs - 1]; 75ba6664aeSJames Wright ubar[1] = prof_ubar[2 * nprofs - 1]; 76ba6664aeSJames Wright ubar[2] = prof_ubar[3 * nprofs - 1]; 77ba6664aeSJames Wright cij[0] = prof_cij[1 * nprofs - 1]; 78ba6664aeSJames Wright cij[1] = prof_cij[2 * nprofs - 1]; 79ba6664aeSJames Wright cij[2] = prof_cij[3 * nprofs - 1]; 80ba6664aeSJames Wright cij[3] = prof_cij[4 * nprofs - 1]; 81ba6664aeSJames Wright cij[4] = prof_cij[5 * nprofs - 1]; 82ba6664aeSJames Wright cij[5] = prof_cij[6 * nprofs - 1]; 83ba6664aeSJames Wright *eps = prof_eps[nprofs - 1]; 84ba6664aeSJames Wright *lt = prof_lt[nprofs - 1]; 85ba6664aeSJames Wright } 86ba6664aeSJames Wright } 87ba6664aeSJames Wright 88ba6664aeSJames Wright /* 89e159aeacSJames Wright * @brief Calculate spectrum coefficient, qn 90e159aeacSJames Wright * 91e159aeacSJames Wright * Calculates q_n at a given distance to the wall 92e159aeacSJames Wright * 93e159aeacSJames Wright * @param[in] kappa nth wavenumber 94e159aeacSJames Wright * @param[in] dkappa Difference between wavenumbers 95e159aeacSJames Wright * @param[in] keta Dissipation wavenumber 96e159aeacSJames Wright * @param[in] kcut Mesh-induced cutoff wavenumber 97e159aeacSJames Wright * @param[in] ke Energy-containing wavenumber 98e159aeacSJames Wright * @param[in] Ektot Total turbulent kinetic energy of spectrum 99e159aeacSJames Wright * @returns qn Spectrum coefficient 100e159aeacSJames Wright */ 1012b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER CeedScalar Calc_qn(const CeedScalar kappa, const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut, 10262e628f8SJames Wright const CeedScalar ke, const CeedScalar Ektot_inv) { 1032b730f8bSJeremy L Thompson const CeedScalar feta_x_fcut = exp(-Square(12 * kappa / keta) - Cube(4 * Max(kappa - 0.9 * kcut, 0) / kcut)); 1042b730f8bSJeremy L Thompson return pow(kappa / ke, 4.) * pow(1 + 2.4 * Square(kappa / ke), -17. / 6) * feta_x_fcut * dkappa * Ektot_inv; 105e159aeacSJames Wright } 106e159aeacSJames Wright 107e159aeacSJames Wright // Calculate hmax, ke, keta, and kcut 1082b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void SpectrumConstants(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 1092b730f8bSJeremy L Thompson const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke, CeedScalar *keta, CeedScalar *kcut) { 110e159aeacSJames Wright *hmax = Max(Max(h[0], h[1]), h[2]); 111175f00a6SJames Wright *ke = wall_dist == 0 ? 1e16 : 2 * M_PI / Min(2 * wall_dist, 3 * lt); 112e159aeacSJames Wright *keta = 2 * M_PI * pow(Cube(nu) / eps, -0.25); 113175f00a6SJames Wright *kcut = M_PI / Min(Max(Max(h[1], h[2]), 0.3 * (*hmax)) + 0.1 * wall_dist, *hmax); 114e159aeacSJames Wright } 115e159aeacSJames Wright 116e159aeacSJames Wright /* 117ba6664aeSJames Wright * @brief Calculate spectrum coefficients for STG 118ba6664aeSJames Wright * 119ba6664aeSJames Wright * Calculates q_n at a given distance to the wall 120ba6664aeSJames Wright * 121175f00a6SJames Wright * @param[in] wall_dist Distance to the nearest wall 122175f00a6SJames Wright * @param[in] eps Turbulent dissipation w/rt wall_dist 123175f00a6SJames Wright * @param[in] lt Turbulent length scale w/rt wall_dist 124ba6664aeSJames Wright * @param[in] h Element lengths in coordinate directions 125ba6664aeSJames Wright * @param[in] nu Dynamic Viscosity; 126ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 127ba6664aeSJames Wright * @param[out] qn Spectrum coefficients, [nmodes] 128ba6664aeSJames Wright */ 1292b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void CalcSpectrum(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 130ba6664aeSJames Wright const CeedScalar nu, CeedScalar qn[], const STGShur14Context stg_ctx) { 131ba6664aeSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 132ba6664aeSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 133e159aeacSJames Wright CeedScalar hmax, ke, keta, kcut, Ektot = 0.0; 1342b730f8bSJeremy L Thompson 135175f00a6SJames Wright SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 136ba6664aeSJames Wright 137ba6664aeSJames Wright for (CeedInt n = 0; n < nmodes; n++) { 138e159aeacSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 139e159aeacSJames Wright qn[n] = Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 140ba6664aeSJames Wright Ektot += qn[n]; 141ba6664aeSJames Wright } 142ba6664aeSJames Wright 143961c9c98SJames Wright if (Ektot == 0) return; 144ba6664aeSJames Wright for (CeedInt n = 0; n < nmodes; n++) qn[n] /= Ektot; 145ba6664aeSJames Wright } 146ba6664aeSJames Wright 147ba6664aeSJames Wright /****************************************************** 148ba6664aeSJames Wright * @brief Calculate u(x,t) for STG inflow condition 149ba6664aeSJames Wright * 150ba6664aeSJames Wright * @param[in] X Location to evaluate u(X,t) 151ba6664aeSJames Wright * @param[in] t Time to evaluate u(X,t) 152ba6664aeSJames Wright * @param[in] ubar Mean velocity at X 153ba6664aeSJames Wright * @param[in] cij Cholesky decomposition at X 154ba6664aeSJames Wright * @param[in] qn Wavemode amplitudes at X, [nmodes] 155ba6664aeSJames Wright * @param[out] u Velocity at X and t 156ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 157ba6664aeSJames Wright */ 1582b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void STGShur14_Calc(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 1592b730f8bSJeremy L Thompson const CeedScalar qn[], CeedScalar u[3], const STGShur14Context stg_ctx) { 160ba6664aeSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 161ba6664aeSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 162ba6664aeSJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 163ba6664aeSJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 164ba6664aeSJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 165ba6664aeSJames Wright CeedScalar xdotd, vp[3] = {0.}; 166ba6664aeSJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 167ba6664aeSJames Wright 1682b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 169ba6664aeSJames Wright xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 170ba6664aeSJames Wright xdotd = 0.; 171ba6664aeSJames Wright for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 172ba6664aeSJames Wright const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 173961c9c98SJames Wright vp[0] += sqrt(qn[n]) * sigma[0 * nmodes + n] * cos_kxdp; 174961c9c98SJames Wright vp[1] += sqrt(qn[n]) * sigma[1 * nmodes + n] * cos_kxdp; 175961c9c98SJames Wright vp[2] += sqrt(qn[n]) * sigma[2 * nmodes + n] * cos_kxdp; 176ba6664aeSJames Wright } 177961c9c98SJames Wright for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 178ba6664aeSJames Wright 179ba6664aeSJames Wright u[0] = ubar[0] + cij[0] * vp[0]; 180ba6664aeSJames Wright u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 181ba6664aeSJames Wright u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 182ba6664aeSJames Wright } 183ba6664aeSJames Wright 184b277271eSJames Wright /****************************************************** 185b277271eSJames Wright * @brief Calculate u(x,t) for STG inflow condition 186b277271eSJames Wright * 187b277271eSJames Wright * @param[in] X Location to evaluate u(X,t) 188b277271eSJames Wright * @param[in] t Time to evaluate u(X,t) 189b277271eSJames Wright * @param[in] ubar Mean velocity at X 190b277271eSJames Wright * @param[in] cij Cholesky decomposition at X 191175f00a6SJames Wright * @param[in] Ektot Total spectrum energy at this location 192175f00a6SJames Wright * @param[in] h Element size in 3 directions 193175f00a6SJames Wright * @param[in] wall_dist Distance to closest wall 194175f00a6SJames Wright * @param[in] eps Turbulent dissipation 195175f00a6SJames Wright * @param[in] lt Turbulent length scale 196b277271eSJames Wright * @param[out] u Velocity at X and t 197b277271eSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 198b277271eSJames Wright */ 1992b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void STGShur14_Calc_PrecompEktot(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 200175f00a6SJames Wright const CeedScalar Ektot, const CeedScalar h[3], const CeedScalar wall_dist, 201b277271eSJames Wright const CeedScalar eps, const CeedScalar lt, const CeedScalar nu, CeedScalar u[3], 202b277271eSJames Wright const STGShur14Context stg_ctx) { 203b277271eSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 204b277271eSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 205b277271eSJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 206b277271eSJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 207b277271eSJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 208b277271eSJames Wright CeedScalar hmax, ke, keta, kcut; 209175f00a6SJames Wright SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 210b277271eSJames Wright CeedScalar xdotd, vp[3] = {0.}; 211b277271eSJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 212b277271eSJames Wright 2132b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 214b277271eSJames Wright xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 215b277271eSJames Wright xdotd = 0.; 216b277271eSJames Wright for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 217b277271eSJames Wright const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 218b277271eSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 219b277271eSJames Wright const CeedScalar qn = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot); 220b277271eSJames Wright vp[0] += sqrt(qn) * sigma[0 * nmodes + n] * cos_kxdp; 221b277271eSJames Wright vp[1] += sqrt(qn) * sigma[1 * nmodes + n] * cos_kxdp; 222b277271eSJames Wright vp[2] += sqrt(qn) * sigma[2 * nmodes + n] * cos_kxdp; 223b277271eSJames Wright } 224b277271eSJames Wright for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 225b277271eSJames Wright 226b277271eSJames Wright u[0] = ubar[0] + cij[0] * vp[0]; 227b277271eSJames Wright u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 228b277271eSJames Wright u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 229b277271eSJames Wright } 230b277271eSJames Wright 23162e628f8SJames Wright // Create preprocessed input for the stg calculation 23262e628f8SJames Wright // 23362e628f8SJames Wright // stg_data[0] = 1 / Ektot (inverse of total spectrum energy) 2342b730f8bSJeremy L Thompson CEED_QFUNCTION(Preprocess_STGShur14)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 235*46603fc5SJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 236*46603fc5SJames Wright const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 237b277271eSJames Wright 238b277271eSJames Wright CeedScalar(*stg_data) = (CeedScalar(*))out[0]; 239b277271eSJames Wright 240b277271eSJames Wright CeedScalar ubar[3], cij[6], eps, lt; 241b277271eSJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 242b277271eSJames Wright const CeedScalar dx = stg_ctx->dx; 243b277271eSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 244b277271eSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 245b277271eSJames Wright const CeedScalar P0 = stg_ctx->P0; 246*46603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 247b277271eSJames Wright const CeedScalar rho = P0 / (Rd * theta0); 248b277271eSJames Wright const CeedScalar nu = mu / rho; 249b277271eSJames Wright 250b277271eSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 251b277271eSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 2525dc40723SJames Wright CeedScalar hmax, ke, keta, kcut; 253b277271eSJames Wright 2542b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 255175f00a6SJames Wright const CeedScalar wall_dist = x[1][i]; 256b277271eSJames Wright const CeedScalar dXdx[2][3] = { 257b277271eSJames Wright {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 258b277271eSJames Wright {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 259b277271eSJames Wright }; 260b277271eSJames Wright 261b277271eSJames Wright CeedScalar h[3]; 262b277271eSJames Wright h[0] = dx; 2632b730f8bSJeremy L Thompson for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(dXdx[0][j] * dXdx[0][j] + dXdx[1][j] * dXdx[1][j]); 264b277271eSJames Wright 265175f00a6SJames Wright InterpolateProfile(wall_dist, ubar, cij, &eps, <, stg_ctx); 266175f00a6SJames Wright SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 267b277271eSJames Wright 268b277271eSJames Wright // Calculate total TKE per spectrum 269d97dc904SJames Wright CeedScalar Ek_tot = 0; 2702b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 271b277271eSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 272d97dc904SJames Wright Ek_tot += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 273b277271eSJames Wright } 274d97dc904SJames Wright // avoid underflowed and poorly defined spectrum coefficients 275d97dc904SJames Wright stg_data[i] = Ek_tot != 0 ? 1 / Ek_tot : 0; 276b277271eSJames Wright } 277b277271eSJames Wright return 0; 278b277271eSJames Wright } 279b277271eSJames Wright 280b77c53c9SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition 2812b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsSTG)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 282b77c53c9SJames Wright // Inputs 283*46603fc5SJames Wright const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 284*46603fc5SJames Wright const CeedScalar(*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 285*46603fc5SJames Wright 286b77c53c9SJames Wright // Outputs 287b77c53c9SJames Wright CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 288b77c53c9SJames Wright 289b77c53c9SJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 29089060322SJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 29189060322SJames Wright const CeedScalar dx = stg_ctx->dx; 29289060322SJames Wright const CeedScalar time = stg_ctx->time; 293b77c53c9SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 294b77c53c9SJames Wright const CeedScalar P0 = stg_ctx->P0; 295b77c53c9SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 296*46603fc5SJames Wright const CeedScalar rho = P0 / (GasConstant(&stg_ctx->newtonian_ctx) * theta0); 297*46603fc5SJames Wright const CeedScalar nu = stg_ctx->newtonian_ctx.mu / rho; 298b77c53c9SJames Wright 2992b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 30089060322SJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 3012b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 3022b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 30389060322SJames Wright {q_data[4][i], q_data[5][i], q_data[6][i]}, 30489060322SJames Wright {q_data[7][i], q_data[8][i], q_data[9][i]} 30589060322SJames Wright }; 30689060322SJames Wright 30789060322SJames Wright CeedScalar h[3]; 30889060322SJames Wright h[0] = dx; 3092b730f8bSJeremy L Thompson for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]) + Square(dXdx[2][j])); 31089060322SJames Wright 31189060322SJames Wright InterpolateProfile(x_i[1], ubar, cij, &eps, <, stg_ctx); 31289060322SJames Wright if (stg_ctx->use_fluctuating_IC) { 313*46603fc5SJames Wright CalcSpectrum(x_i[1], eps, lt, h, nu, qn, stg_ctx); 31489060322SJames Wright STGShur14_Calc(x_i, time, ubar, cij, qn, u, stg_ctx); 31589060322SJames Wright } else { 31689060322SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 31789060322SJames Wright } 318b77c53c9SJames Wright 31997baf651SJames Wright switch (stg_ctx->newtonian_ctx.state_var) { 32097baf651SJames Wright case STATEVAR_CONSERVATIVE: 321b77c53c9SJames Wright q0[0][i] = rho; 322b77c53c9SJames Wright q0[1][i] = u[0] * rho; 323b77c53c9SJames Wright q0[2][i] = u[1] * rho; 324b77c53c9SJames Wright q0[3][i] = u[2] * rho; 325b77c53c9SJames Wright q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0); 32697baf651SJames Wright break; 32797baf651SJames Wright 32897baf651SJames Wright case STATEVAR_PRIMITIVE: 32997baf651SJames Wright q0[0][i] = P0; 33097baf651SJames Wright q0[1][i] = u[0]; 33197baf651SJames Wright q0[2][i] = u[1]; 33297baf651SJames Wright q0[3][i] = u[2]; 33397baf651SJames Wright q0[4][i] = theta0; 33497baf651SJames Wright break; 3357c4551aaSJames Wright } 336b77c53c9SJames Wright } // End of Quadrature Point Loop 337b77c53c9SJames Wright return 0; 338b77c53c9SJames Wright } 339b77c53c9SJames Wright 340ba6664aeSJames Wright /******************************************************************** 341ba6664aeSJames Wright * @brief QFunction to calculate the inflow boundary condition 342ba6664aeSJames Wright * 343ba6664aeSJames Wright * This will loop through quadrature points, calculate the wavemode amplitudes 344ba6664aeSJames Wright * at each location, then calculate the actual velocity. 345ba6664aeSJames Wright */ 3462b730f8bSJeremy L Thompson CEED_QFUNCTION(STGShur14_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 347*46603fc5SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 348*46603fc5SJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 349*46603fc5SJames Wright const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 350ba6664aeSJames Wright 351*46603fc5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 352*46603fc5SJames Wright CeedScalar(*jac_data_sur)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[1]; 353ba6664aeSJames Wright 354ba6664aeSJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 355ba6664aeSJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 356ba6664aeSJames Wright const bool is_implicit = stg_ctx->is_implicit; 357ba6664aeSJames Wright const bool mean_only = stg_ctx->mean_only; 358ba6664aeSJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 359ba6664aeSJames Wright const CeedScalar dx = stg_ctx->dx; 360ba6664aeSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 361ba6664aeSJames Wright const CeedScalar time = stg_ctx->time; 362ba6664aeSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 363ba6664aeSJames Wright const CeedScalar P0 = stg_ctx->P0; 364ba6664aeSJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 365*46603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 366*46603fc5SJames Wright const CeedScalar gamma = HeatCapacityRatio(&stg_ctx->newtonian_ctx); 367ba6664aeSJames Wright 3682b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 369ba6664aeSJames Wright const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0); 370ba6664aeSJames Wright const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 371ba6664aeSJames Wright const CeedScalar dXdx[2][3] = { 372ba6664aeSJames Wright {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 373ba6664aeSJames Wright {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 374ba6664aeSJames Wright }; 375ba6664aeSJames Wright 376ba6664aeSJames Wright CeedScalar h[3]; 377ba6664aeSJames Wright h[0] = dx; 3782b730f8bSJeremy L Thompson for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 379ba6664aeSJames Wright 380ba6664aeSJames Wright InterpolateProfile(X[1][i], ubar, cij, &eps, <, stg_ctx); 381ba6664aeSJames Wright if (!mean_only) { 382ba6664aeSJames Wright CalcSpectrum(X[1][i], eps, lt, h, mu / rho, qn, stg_ctx); 383ba6664aeSJames Wright STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 384ba6664aeSJames Wright } else { 385ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 386ba6664aeSJames Wright } 387ba6664aeSJames Wright 3884dbab5e5SJames Wright const CeedScalar E_kinetic = .5 * rho * Dot3(u, u); 389ba6664aeSJames Wright CeedScalar E_internal, P; 390ba6664aeSJames Wright if (prescribe_T) { 391ba6664aeSJames Wright // Temperature is being set weakly (theta0) and for constant cv this sets E_internal 392ba6664aeSJames Wright E_internal = rho * cv * theta0; 393ba6664aeSJames Wright // Find pressure using 394ba6664aeSJames Wright P = rho * Rd * theta0; // interior rho with exterior T 395ba6664aeSJames Wright } else { 396ba6664aeSJames Wright E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution 397ba6664aeSJames Wright P = E_internal * (gamma - 1.); 398ba6664aeSJames Wright } 399ba6664aeSJames Wright 400ba6664aeSJames Wright const CeedScalar wdetJb = (is_implicit ? -1. : 1.) * q_data_sur[0][i]; 401ba6664aeSJames Wright // ---- Normal vect 4022b730f8bSJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 403ba6664aeSJames Wright 404ba6664aeSJames Wright const CeedScalar E = E_internal + E_kinetic; 405ba6664aeSJames Wright 406ba6664aeSJames Wright // Velocity normal to the boundary 4074dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, u); 4084dbab5e5SJames Wright 409ba6664aeSJames Wright // The Physics 410ba6664aeSJames Wright // Zero v so all future terms can safely sum into it 411ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 412ba6664aeSJames Wright 413ba6664aeSJames Wright // The Physics 414ba6664aeSJames Wright // -- Density 415ba6664aeSJames Wright v[0][i] -= wdetJb * rho * u_normal; 416ba6664aeSJames Wright 417ba6664aeSJames Wright // -- Momentum 4182b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 419ba6664aeSJames Wright 420ba6664aeSJames Wright // -- Total Energy Density 421ba6664aeSJames Wright v[4][i] -= wdetJb * u_normal * (E + P); 4224dbab5e5SJames Wright 4234dbab5e5SJames Wright jac_data_sur[0][i] = rho; 4244dbab5e5SJames Wright jac_data_sur[1][i] = u[0]; 4254dbab5e5SJames Wright jac_data_sur[2][i] = u[1]; 4264dbab5e5SJames Wright jac_data_sur[3][i] = u[2]; 4274dbab5e5SJames Wright jac_data_sur[4][i] = E; 4284dbab5e5SJames Wright for (int j = 0; j < 6; j++) jac_data_sur[5 + j][i] = 0.; 429ba6664aeSJames Wright } 430ba6664aeSJames Wright return 0; 431ba6664aeSJames Wright } 432ba6664aeSJames Wright 4332b730f8bSJeremy L Thompson CEED_QFUNCTION(STGShur14_Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 4344dbab5e5SJames Wright // Inputs 435*46603fc5SJames Wright const CeedScalar(*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 436*46603fc5SJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 437*46603fc5SJames Wright const CeedScalar(*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 4384dbab5e5SJames Wright // Outputs 4394dbab5e5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 440*46603fc5SJames Wright 4414dbab5e5SJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 4424dbab5e5SJames Wright const bool implicit = stg_ctx->is_implicit; 4434dbab5e5SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 444*46603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 445*46603fc5SJames Wright const CeedScalar gamma = HeatCapacityRatio(&stg_ctx->newtonian_ctx); 4464dbab5e5SJames Wright 4474dbab5e5SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 4484dbab5e5SJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 4494dbab5e5SJames Wright 4504dbab5e5SJames Wright CeedPragmaSIMD 4514dbab5e5SJames Wright // Quadrature Point Loop 4524dbab5e5SJames Wright for (CeedInt i = 0; i < Q; i++) { 4534dbab5e5SJames Wright // Setup 4544dbab5e5SJames Wright // -- Interp-to-Interp q_data 4554dbab5e5SJames Wright // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 4564dbab5e5SJames Wright // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 4574dbab5e5SJames Wright // We can effect this by swapping the sign on this weight 4584dbab5e5SJames Wright const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 4594dbab5e5SJames Wright 4604dbab5e5SJames Wright // Calculate inflow values 4614dbab5e5SJames Wright CeedScalar velocity[3]; 4624dbab5e5SJames Wright for (CeedInt j = 0; j < 3; j++) velocity[j] = jac_data_sur[5 + j][i]; 4634dbab5e5SJames Wright 4644dbab5e5SJames Wright // enabling user to choose between weak T and weak rho inflow 4654dbab5e5SJames Wright CeedScalar drho, dE, dP; 4664dbab5e5SJames Wright if (prescribe_T) { 4674dbab5e5SJames Wright // rho should be from the current solution 4684dbab5e5SJames Wright drho = dq[0][i]; 4694dbab5e5SJames Wright CeedScalar dE_internal = drho * cv * theta0; 4704dbab5e5SJames Wright CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity); 4714dbab5e5SJames Wright dE = dE_internal + dE_kinetic; 4724dbab5e5SJames Wright dP = drho * Rd * theta0; // interior rho with exterior T 4734dbab5e5SJames Wright } else { // rho specified, E_internal from solution 4744dbab5e5SJames Wright drho = 0; 4754dbab5e5SJames Wright dE = dq[4][i]; 4764dbab5e5SJames Wright dP = dE * (gamma - 1.); 4774dbab5e5SJames Wright } 4782b730f8bSJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 4794dbab5e5SJames Wright 4804dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, velocity); 4814dbab5e5SJames Wright 4824dbab5e5SJames Wright v[0][i] = -wdetJb * drho * u_normal; 4832b730f8bSJeremy L Thompson for (int j = 0; j < 3; j++) v[j + 1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP); 4844dbab5e5SJames Wright v[4][i] = -wdetJb * u_normal * (dE + dP); 4854dbab5e5SJames Wright } // End Quadrature Point Loop 4864dbab5e5SJames Wright return 0; 4874dbab5e5SJames Wright } 4884dbab5e5SJames Wright 4890a6353c2SJames Wright /******************************************************************** 4900a6353c2SJames Wright * @brief QFunction to calculate the strongly enforce inflow BC 4910a6353c2SJames Wright * 4920a6353c2SJames Wright * This QF is for the strong application of STG via libCEED (rather than 4930a6353c2SJames Wright * through the native PETSc `DMAddBoundary` -> `bcFunc` method. 4940a6353c2SJames Wright */ 4952b730f8bSJeremy L Thompson CEED_QFUNCTION(STGShur14_Inflow_StrongQF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 496*46603fc5SJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 497*46603fc5SJames Wright const CeedScalar(*coords)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 498*46603fc5SJames Wright const CeedScalar(*scale) = (const CeedScalar(*))in[2]; 499*46603fc5SJames Wright const CeedScalar(*stg_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 5000a6353c2SJames Wright 5010a6353c2SJames Wright CeedScalar(*bcval)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 5020a6353c2SJames Wright 5030a6353c2SJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 50462e628f8SJames Wright CeedScalar u[3], ubar[3], cij[6], eps, lt; 5050a6353c2SJames Wright const bool mean_only = stg_ctx->mean_only; 5060a6353c2SJames Wright const CeedScalar dx = stg_ctx->dx; 5070a6353c2SJames Wright const CeedScalar time = stg_ctx->time; 5080a6353c2SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 5090a6353c2SJames Wright const CeedScalar P0 = stg_ctx->P0; 510*46603fc5SJames Wright const CeedScalar rho = P0 / (GasConstant(&stg_ctx->newtonian_ctx) * theta0); 511*46603fc5SJames Wright const CeedScalar nu = stg_ctx->newtonian_ctx.mu / rho; 5120a6353c2SJames Wright 5132b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 5140a6353c2SJames Wright const CeedScalar x[] = {coords[0][i], coords[1][i], coords[2][i]}; 5150a6353c2SJames Wright const CeedScalar dXdx[2][3] = { 5160a6353c2SJames Wright {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 5170a6353c2SJames Wright {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 5180a6353c2SJames Wright }; 5190a6353c2SJames Wright 5200a6353c2SJames Wright CeedScalar h[3]; 5210a6353c2SJames Wright h[0] = dx; 5222b730f8bSJeremy L Thompson for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 5230a6353c2SJames Wright 5240a6353c2SJames Wright InterpolateProfile(coords[1][i], ubar, cij, &eps, <, stg_ctx); 5250a6353c2SJames Wright if (!mean_only) { 52662e628f8SJames Wright if (1) { 527*46603fc5SJames Wright STGShur14_Calc_PrecompEktot(x, time, ubar, cij, stg_data[0][i], h, x[1], eps, lt, nu, u, stg_ctx); 52862e628f8SJames Wright } else { // Original way 52962e628f8SJames Wright CeedScalar qn[STG_NMODES_MAX]; 530*46603fc5SJames Wright CalcSpectrum(coords[1][i], eps, lt, h, nu, qn, stg_ctx); 53162e628f8SJames Wright STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 53262e628f8SJames Wright } 5330a6353c2SJames Wright } else { 5340a6353c2SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 5350a6353c2SJames Wright } 5360a6353c2SJames Wright 53797baf651SJames Wright switch (stg_ctx->newtonian_ctx.state_var) { 53897baf651SJames Wright case STATEVAR_CONSERVATIVE: 5390a6353c2SJames Wright bcval[0][i] = scale[i] * rho; 5400a6353c2SJames Wright bcval[1][i] = scale[i] * rho * u[0]; 5410a6353c2SJames Wright bcval[2][i] = scale[i] * rho * u[1]; 5420a6353c2SJames Wright bcval[3][i] = scale[i] * rho * u[2]; 543cf3d54ffSJames Wright bcval[4][i] = 0.; 54497baf651SJames Wright break; 54597baf651SJames Wright 54697baf651SJames Wright case STATEVAR_PRIMITIVE: 54797baf651SJames Wright bcval[0][i] = 0; 54897baf651SJames Wright bcval[1][i] = scale[i] * u[0]; 54997baf651SJames Wright bcval[2][i] = scale[i] * u[1]; 55097baf651SJames Wright bcval[3][i] = scale[i] * u[2]; 55197baf651SJames Wright bcval[4][i] = scale[i] * theta0; 55297baf651SJames Wright break; 5530a6353c2SJames Wright } 5547c4551aaSJames Wright } 5550a6353c2SJames Wright return 0; 5560a6353c2SJames Wright } 5570a6353c2SJames Wright 558ba6664aeSJames Wright #endif // stg_shur14_h 559