15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 // 12ea61e9acSJeremy L Thompson /// SetupSTG_Rand reads in the input files and fills in STGShur14Context. 13ea61e9acSJeremy L Thompson /// Then STGShur14_CalcQF is run over quadrature points. 14ea61e9acSJeremy L Thompson /// Before the program exits, TearDownSTG is run to free the memory of the allocated arrays. 15ba6664aeSJames Wright #include <ceed.h> 16c9c2c079SJeremy L Thompson #include <math.h> 17ba6664aeSJames Wright #include <stdlib.h> 182b730f8bSJeremy L Thompson 1946603fc5SJames Wright #include "newtonian_state.h" 208756a6ccSJames Wright #include "setupgeo_helpers.h" 21ba6664aeSJames Wright #include "stg_shur14_type.h" 2213fa47b2SJames Wright #include "utils.h" 23ba6664aeSJames Wright 24ba6664aeSJames Wright #define STG_NMODES_MAX 1024 25ba6664aeSJames Wright 26ba6664aeSJames Wright /* 27ba6664aeSJames Wright * @brief Interpolate quantities from input profile to given location 28ba6664aeSJames Wright * 29175f00a6SJames Wright * Assumed that prof_wd[i+1] > prof_wd[i] and prof_wd[0] = 0 30175f00a6SJames Wright * If wall_dist > prof_wd[-1], then the interpolation takes the values at prof_wd[-1] 31ba6664aeSJames Wright * 32175f00a6SJames Wright * @param[in] wall_dist Distance to the nearest wall 33175f00a6SJames Wright * @param[out] ubar Mean velocity at wall_dist 34175f00a6SJames Wright * @param[out] cij Cholesky decomposition at wall_dist 35175f00a6SJames Wright * @param[out] eps Turbulent dissipation at wall_dist 36175f00a6SJames Wright * @param[out] lt Turbulent length scale at wall_dist 37ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 38ba6664aeSJames Wright */ 392b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar wall_dist, CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt, 40cbef7084SJames Wright const StgShur14Context stg_ctx) { 41ba6664aeSJames Wright const CeedInt nprofs = stg_ctx->nprofs; 42175f00a6SJames Wright const CeedScalar *prof_wd = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 43ba6664aeSJames Wright const CeedScalar *prof_eps = &stg_ctx->data[stg_ctx->offsets.eps]; 44ba6664aeSJames Wright const CeedScalar *prof_lt = &stg_ctx->data[stg_ctx->offsets.lt]; 45ba6664aeSJames Wright const CeedScalar *prof_ubar = &stg_ctx->data[stg_ctx->offsets.ubar]; 46ba6664aeSJames Wright const CeedScalar *prof_cij = &stg_ctx->data[stg_ctx->offsets.cij]; 47ba6664aeSJames Wright CeedInt idx = -1; 48ba6664aeSJames Wright 49ba6664aeSJames Wright for (CeedInt i = 0; i < nprofs; i++) { 50175f00a6SJames Wright if (wall_dist < prof_wd[i]) { 51ba6664aeSJames Wright idx = i; 52ba6664aeSJames Wright break; 53ba6664aeSJames Wright } 54ba6664aeSJames Wright } 55ba6664aeSJames Wright 56175f00a6SJames Wright if (idx > 0) { // y within the bounds of prof_wd 57175f00a6SJames Wright CeedScalar coeff = (wall_dist - prof_wd[idx - 1]) / (prof_wd[idx] - prof_wd[idx - 1]); 58175f00a6SJames Wright 59ba6664aeSJames Wright ubar[0] = prof_ubar[0 * nprofs + idx - 1] + coeff * (prof_ubar[0 * nprofs + idx] - prof_ubar[0 * nprofs + idx - 1]); 60ba6664aeSJames Wright ubar[1] = prof_ubar[1 * nprofs + idx - 1] + coeff * (prof_ubar[1 * nprofs + idx] - prof_ubar[1 * nprofs + idx - 1]); 61ba6664aeSJames Wright ubar[2] = prof_ubar[2 * nprofs + idx - 1] + coeff * (prof_ubar[2 * nprofs + idx] - prof_ubar[2 * nprofs + idx - 1]); 62ba6664aeSJames Wright cij[0] = prof_cij[0 * nprofs + idx - 1] + coeff * (prof_cij[0 * nprofs + idx] - prof_cij[0 * nprofs + idx - 1]); 63ba6664aeSJames Wright cij[1] = prof_cij[1 * nprofs + idx - 1] + coeff * (prof_cij[1 * nprofs + idx] - prof_cij[1 * nprofs + idx - 1]); 64ba6664aeSJames Wright cij[2] = prof_cij[2 * nprofs + idx - 1] + coeff * (prof_cij[2 * nprofs + idx] - prof_cij[2 * nprofs + idx - 1]); 65ba6664aeSJames Wright cij[3] = prof_cij[3 * nprofs + idx - 1] + coeff * (prof_cij[3 * nprofs + idx] - prof_cij[3 * nprofs + idx - 1]); 66ba6664aeSJames Wright cij[4] = prof_cij[4 * nprofs + idx - 1] + coeff * (prof_cij[4 * nprofs + idx] - prof_cij[4 * nprofs + idx - 1]); 67ba6664aeSJames Wright cij[5] = prof_cij[5 * nprofs + idx - 1] + coeff * (prof_cij[5 * nprofs + idx] - prof_cij[5 * nprofs + idx - 1]); 68ba6664aeSJames Wright *eps = prof_eps[idx - 1] + coeff * (prof_eps[idx] - prof_eps[idx - 1]); 69ba6664aeSJames Wright *lt = prof_lt[idx - 1] + coeff * (prof_lt[idx] - prof_lt[idx - 1]); 70175f00a6SJames Wright } else { // y outside bounds of prof_wd 71ba6664aeSJames Wright ubar[0] = prof_ubar[1 * nprofs - 1]; 72ba6664aeSJames Wright ubar[1] = prof_ubar[2 * nprofs - 1]; 73ba6664aeSJames Wright ubar[2] = prof_ubar[3 * nprofs - 1]; 74ba6664aeSJames Wright cij[0] = prof_cij[1 * nprofs - 1]; 75ba6664aeSJames Wright cij[1] = prof_cij[2 * nprofs - 1]; 76ba6664aeSJames Wright cij[2] = prof_cij[3 * nprofs - 1]; 77ba6664aeSJames Wright cij[3] = prof_cij[4 * nprofs - 1]; 78ba6664aeSJames Wright cij[4] = prof_cij[5 * nprofs - 1]; 79ba6664aeSJames Wright cij[5] = prof_cij[6 * nprofs - 1]; 80ba6664aeSJames Wright *eps = prof_eps[nprofs - 1]; 81ba6664aeSJames Wright *lt = prof_lt[nprofs - 1]; 82ba6664aeSJames Wright } 83ba6664aeSJames Wright } 84ba6664aeSJames Wright 85ba6664aeSJames Wright /* 86e159aeacSJames Wright * @brief Calculate spectrum coefficient, qn 87e159aeacSJames Wright * 88e159aeacSJames Wright * Calculates q_n at a given distance to the wall 89e159aeacSJames Wright * 90e159aeacSJames Wright * @param[in] kappa nth wavenumber 91e159aeacSJames Wright * @param[in] dkappa Difference between wavenumbers 92e159aeacSJames Wright * @param[in] keta Dissipation wavenumber 93e159aeacSJames Wright * @param[in] kcut Mesh-induced cutoff wavenumber 94e159aeacSJames Wright * @param[in] ke Energy-containing wavenumber 95f8839eb4SJames Wright * @param[in] Ektot_inv Inverse of total turbulent kinetic energy of spectrum 96e159aeacSJames Wright * @returns qn Spectrum coefficient 97e159aeacSJames Wright */ 982b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER CeedScalar Calc_qn(const CeedScalar kappa, const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut, 9962e628f8SJames Wright const CeedScalar ke, const CeedScalar Ektot_inv) { 1002b730f8bSJeremy L Thompson const CeedScalar feta_x_fcut = exp(-Square(12 * kappa / keta) - Cube(4 * Max(kappa - 0.9 * kcut, 0) / kcut)); 1012b730f8bSJeremy L Thompson return pow(kappa / ke, 4.) * pow(1 + 2.4 * Square(kappa / ke), -17. / 6) * feta_x_fcut * dkappa * Ektot_inv; 102e159aeacSJames Wright } 103e159aeacSJames Wright 104e159aeacSJames Wright // Calculate hmax, ke, keta, and kcut 105*831dbe9eSJames Wright CEED_QFUNCTION_HELPER void SpectrumConstants(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar hNodSep[3], 1062b730f8bSJeremy L Thompson const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke, CeedScalar *keta, CeedScalar *kcut) { 107*831dbe9eSJames Wright *hmax = Max(Max(hNodSep[0], hNodSep[1]), hNodSep[2]); 108175f00a6SJames Wright *ke = wall_dist == 0 ? 1e16 : 2 * M_PI / Min(2 * wall_dist, 3 * lt); 109e159aeacSJames Wright *keta = 2 * M_PI * pow(Cube(nu) / eps, -0.25); 110*831dbe9eSJames Wright *kcut = M_PI / Min(Max(Max(hNodSep[1], hNodSep[2]), 0.3 * (*hmax)) + 0.1 * wall_dist, *hmax); 111e159aeacSJames Wright } 112e159aeacSJames Wright 113e159aeacSJames Wright /* 114ba6664aeSJames Wright * @brief Calculate spectrum coefficients for STG 115ba6664aeSJames Wright * 116ba6664aeSJames Wright * Calculates q_n at a given distance to the wall 117ba6664aeSJames Wright * 118175f00a6SJames Wright * @param[in] wall_dist Distance to the nearest wall 119175f00a6SJames Wright * @param[in] eps Turbulent dissipation w/rt wall_dist 120175f00a6SJames Wright * @param[in] lt Turbulent length scale w/rt wall_dist 121*831dbe9eSJames Wright * @param[in] h_node_sep Element lengths in coordinate directions 122ba6664aeSJames Wright * @param[in] nu Dynamic Viscosity; 123ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 124ba6664aeSJames Wright * @param[out] qn Spectrum coefficients, [nmodes] 125ba6664aeSJames Wright */ 126*831dbe9eSJames Wright CEED_QFUNCTION_HELPER void CalcSpectrum(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h_node_sep[3], 127cbef7084SJames Wright const CeedScalar nu, CeedScalar qn[], const StgShur14Context stg_ctx) { 128ba6664aeSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 129ba6664aeSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 130e159aeacSJames Wright CeedScalar hmax, ke, keta, kcut, Ektot = 0.0; 1312b730f8bSJeremy L Thompson 132*831dbe9eSJames Wright SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut); 133ba6664aeSJames Wright 134ba6664aeSJames Wright for (CeedInt n = 0; n < nmodes; n++) { 135e159aeacSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 136e159aeacSJames Wright qn[n] = Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 137ba6664aeSJames Wright Ektot += qn[n]; 138ba6664aeSJames Wright } 139ba6664aeSJames Wright 140961c9c98SJames Wright if (Ektot == 0) return; 141ba6664aeSJames Wright for (CeedInt n = 0; n < nmodes; n++) qn[n] /= Ektot; 142ba6664aeSJames Wright } 143ba6664aeSJames Wright 144ba6664aeSJames Wright /****************************************************** 145ba6664aeSJames Wright * @brief Calculate u(x,t) for STG inflow condition 146ba6664aeSJames Wright * 147ba6664aeSJames Wright * @param[in] X Location to evaluate u(X,t) 148ba6664aeSJames Wright * @param[in] t Time to evaluate u(X,t) 149ba6664aeSJames Wright * @param[in] ubar Mean velocity at X 150ba6664aeSJames Wright * @param[in] cij Cholesky decomposition at X 151ba6664aeSJames Wright * @param[in] qn Wavemode amplitudes at X, [nmodes] 152ba6664aeSJames Wright * @param[out] u Velocity at X and t 153ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 154ba6664aeSJames Wright */ 155cbef7084SJames Wright CEED_QFUNCTION_HELPER void StgShur14Calc(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 156cbef7084SJames Wright const CeedScalar qn[], CeedScalar u[3], const StgShur14Context stg_ctx) { 157ba6664aeSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 158ba6664aeSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 159ba6664aeSJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 160ba6664aeSJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 161ba6664aeSJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 162ba6664aeSJames Wright CeedScalar xdotd, vp[3] = {0.}; 163ba6664aeSJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 164ba6664aeSJames Wright 1652b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 166ba6664aeSJames Wright xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 167ba6664aeSJames Wright xdotd = 0.; 168ba6664aeSJames Wright for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 169ba6664aeSJames Wright const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 170961c9c98SJames Wright vp[0] += sqrt(qn[n]) * sigma[0 * nmodes + n] * cos_kxdp; 171961c9c98SJames Wright vp[1] += sqrt(qn[n]) * sigma[1 * nmodes + n] * cos_kxdp; 172961c9c98SJames Wright vp[2] += sqrt(qn[n]) * sigma[2 * nmodes + n] * cos_kxdp; 173ba6664aeSJames Wright } 174961c9c98SJames Wright for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 175ba6664aeSJames Wright 176ba6664aeSJames Wright u[0] = ubar[0] + cij[0] * vp[0]; 177ba6664aeSJames Wright u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 178ba6664aeSJames Wright u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 179ba6664aeSJames Wright } 180ba6664aeSJames Wright 181b277271eSJames Wright /****************************************************** 182b277271eSJames Wright * @brief Calculate u(x,t) for STG inflow condition 183b277271eSJames Wright * 184b277271eSJames Wright * @param[in] X Location to evaluate u(X,t) 185b277271eSJames Wright * @param[in] t Time to evaluate u(X,t) 186b277271eSJames Wright * @param[in] ubar Mean velocity at X 187b277271eSJames Wright * @param[in] cij Cholesky decomposition at X 188175f00a6SJames Wright * @param[in] Ektot Total spectrum energy at this location 189*831dbe9eSJames Wright * @param[in] h_node_sep Element size in 3 directions 190175f00a6SJames Wright * @param[in] wall_dist Distance to closest wall 191175f00a6SJames Wright * @param[in] eps Turbulent dissipation 192175f00a6SJames Wright * @param[in] lt Turbulent length scale 193b277271eSJames Wright * @param[out] u Velocity at X and t 194b277271eSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 195b277271eSJames Wright */ 196cbef7084SJames Wright CEED_QFUNCTION_HELPER void StgShur14Calc_PrecompEktot(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 197*831dbe9eSJames Wright const CeedScalar Ektot, const CeedScalar h_node_sep[3], const CeedScalar wall_dist, 198*831dbe9eSJames Wright const CeedScalar eps, const CeedScalar lt, const CeedScalar nu, CeedScalar u[3], 199*831dbe9eSJames Wright const StgShur14Context stg_ctx) { 200b277271eSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 201b277271eSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 202b277271eSJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 203b277271eSJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 204b277271eSJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 205b277271eSJames Wright CeedScalar hmax, ke, keta, kcut; 206*831dbe9eSJames Wright SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut); 207b277271eSJames Wright CeedScalar xdotd, vp[3] = {0.}; 208b277271eSJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 209b277271eSJames Wright 2102b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 211b277271eSJames Wright xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 212b277271eSJames Wright xdotd = 0.; 213b277271eSJames Wright for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 214b277271eSJames Wright const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 215b277271eSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 216b277271eSJames Wright const CeedScalar qn = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot); 217b277271eSJames Wright vp[0] += sqrt(qn) * sigma[0 * nmodes + n] * cos_kxdp; 218b277271eSJames Wright vp[1] += sqrt(qn) * sigma[1 * nmodes + n] * cos_kxdp; 219b277271eSJames Wright vp[2] += sqrt(qn) * sigma[2 * nmodes + n] * cos_kxdp; 220b277271eSJames Wright } 221b277271eSJames Wright for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 222b277271eSJames Wright 223b277271eSJames Wright u[0] = ubar[0] + cij[0] * vp[0]; 224b277271eSJames Wright u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 225b277271eSJames Wright u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 226b277271eSJames Wright } 227b277271eSJames Wright 22862e628f8SJames Wright // Create preprocessed input for the stg calculation 22962e628f8SJames Wright // 23062e628f8SJames Wright // stg_data[0] = 1 / Ektot (inverse of total spectrum energy) 231cbef7084SJames Wright CEED_QFUNCTION(StgShur14Preprocess)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 23266170c20SJames Wright const CeedScalar(*dXdx_q)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0]; 23346603fc5SJames Wright const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 234b277271eSJames Wright 235b277271eSJames Wright CeedScalar(*stg_data) = (CeedScalar(*))out[0]; 236b277271eSJames Wright 237b277271eSJames Wright CeedScalar ubar[3], cij[6], eps, lt; 238cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 239b277271eSJames Wright const CeedScalar dx = stg_ctx->dx; 240b277271eSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 241b277271eSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 242b277271eSJames Wright const CeedScalar P0 = stg_ctx->P0; 24346603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 244b277271eSJames Wright const CeedScalar rho = P0 / (Rd * theta0); 245b277271eSJames Wright const CeedScalar nu = mu / rho; 246b277271eSJames Wright 247b277271eSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 248b277271eSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 2495dc40723SJames Wright CeedScalar hmax, ke, keta, kcut; 250b277271eSJames Wright 2512b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 252175f00a6SJames Wright const CeedScalar wall_dist = x[1][i]; 253b277271eSJames Wright const CeedScalar dXdx[2][3] = { 25466170c20SJames Wright {dXdx_q[0][0][i], dXdx_q[0][1][i], dXdx_q[0][2][i]}, 25566170c20SJames Wright {dXdx_q[1][0][i], dXdx_q[1][1][i], dXdx_q[1][2][i]}, 256b277271eSJames Wright }; 257b277271eSJames Wright 258*831dbe9eSJames Wright CeedScalar h_node_sep[3]; 259*831dbe9eSJames Wright h_node_sep[0] = dx; 260*831dbe9eSJames Wright for (CeedInt j = 1; j < 3; j++) h_node_sep[j] = 2 / sqrt(dXdx[0][j] * dXdx[0][j] + dXdx[1][j] * dXdx[1][j]); 261*831dbe9eSJames Wright ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3); 262b277271eSJames Wright 263175f00a6SJames Wright InterpolateProfile(wall_dist, ubar, cij, &eps, <, stg_ctx); 264*831dbe9eSJames Wright SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut); 265b277271eSJames Wright 266b277271eSJames Wright // Calculate total TKE per spectrum 267d97dc904SJames Wright CeedScalar Ek_tot = 0; 2682b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 269b277271eSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 270d97dc904SJames Wright Ek_tot += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 271b277271eSJames Wright } 272d97dc904SJames Wright // avoid underflowed and poorly defined spectrum coefficients 273d97dc904SJames Wright stg_data[i] = Ek_tot != 0 ? 1 / Ek_tot : 0; 274b277271eSJames Wright } 275b277271eSJames Wright return 0; 276b277271eSJames Wright } 277b277271eSJames Wright 278b77c53c9SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition 279cbef7084SJames Wright CEED_QFUNCTION(ICsStg)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 28046603fc5SJames Wright const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 2811c299e57SJeremy L Thompson const CeedScalar(*J)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[1]; 282b77c53c9SJames Wright CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 283b77c53c9SJames Wright 284cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 28589060322SJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 28689060322SJames Wright const CeedScalar dx = stg_ctx->dx; 28789060322SJames Wright const CeedScalar time = stg_ctx->time; 288b77c53c9SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 289b77c53c9SJames Wright const CeedScalar P0 = stg_ctx->P0; 290b77c53c9SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 29146603fc5SJames Wright const CeedScalar rho = P0 / (GasConstant(&stg_ctx->newtonian_ctx) * theta0); 29246603fc5SJames Wright const CeedScalar nu = stg_ctx->newtonian_ctx.mu / rho; 293b77c53c9SJames Wright 2942b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 29589060322SJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 2961c299e57SJeremy L Thompson CeedScalar dXdx[3][3]; 2978756a6ccSJames Wright InvertMappingJacobian_3D(Q, i, J, dXdx, NULL); 298*831dbe9eSJames Wright CeedScalar h_node_sep[3]; 299*831dbe9eSJames Wright h_node_sep[0] = dx; 300*831dbe9eSJames Wright for (CeedInt j = 1; j < 3; j++) h_node_sep[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]) + Square(dXdx[2][j])); 301*831dbe9eSJames Wright ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3); 30289060322SJames Wright 30389060322SJames Wright InterpolateProfile(x_i[1], ubar, cij, &eps, <, stg_ctx); 30489060322SJames Wright if (stg_ctx->use_fluctuating_IC) { 305*831dbe9eSJames Wright CalcSpectrum(x_i[1], eps, lt, h_node_sep, nu, qn, stg_ctx); 306cbef7084SJames Wright StgShur14Calc(x_i, time, ubar, cij, qn, u, stg_ctx); 30789060322SJames Wright } else { 30889060322SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 30989060322SJames Wright } 310b77c53c9SJames Wright 31197baf651SJames Wright switch (stg_ctx->newtonian_ctx.state_var) { 31297baf651SJames Wright case STATEVAR_CONSERVATIVE: 313b77c53c9SJames Wright q0[0][i] = rho; 314b77c53c9SJames Wright q0[1][i] = u[0] * rho; 315b77c53c9SJames Wright q0[2][i] = u[1] * rho; 316b77c53c9SJames Wright q0[3][i] = u[2] * rho; 317b77c53c9SJames Wright q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0); 31897baf651SJames Wright break; 31997baf651SJames Wright 32097baf651SJames Wright case STATEVAR_PRIMITIVE: 32197baf651SJames Wright q0[0][i] = P0; 32297baf651SJames Wright q0[1][i] = u[0]; 32397baf651SJames Wright q0[2][i] = u[1]; 32497baf651SJames Wright q0[3][i] = u[2]; 32597baf651SJames Wright q0[4][i] = theta0; 32697baf651SJames Wright break; 3277c4551aaSJames Wright } 328f0b01153SJames Wright } 329b77c53c9SJames Wright return 0; 330b77c53c9SJames Wright } 331b77c53c9SJames Wright 332ba6664aeSJames Wright /******************************************************************** 333ba6664aeSJames Wright * @brief QFunction to calculate the inflow boundary condition 334ba6664aeSJames Wright * 335ba6664aeSJames Wright * This will loop through quadrature points, calculate the wavemode amplitudes 336ba6664aeSJames Wright * at each location, then calculate the actual velocity. 337ba6664aeSJames Wright */ 338cbef7084SJames Wright CEED_QFUNCTION(StgShur14Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 33946603fc5SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 340f3e15844SJames Wright const CeedScalar(*q_data_sur) = in[2]; 34146603fc5SJames Wright const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 342ba6664aeSJames Wright 34346603fc5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 344f3e15844SJames Wright CeedScalar(*jac_data_sur) = out[1]; 345ba6664aeSJames Wright 346cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 347ba6664aeSJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 348ba6664aeSJames Wright const bool is_implicit = stg_ctx->is_implicit; 349ba6664aeSJames Wright const bool mean_only = stg_ctx->mean_only; 350ba6664aeSJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 351ba6664aeSJames Wright const CeedScalar dx = stg_ctx->dx; 352ba6664aeSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 353ba6664aeSJames Wright const CeedScalar time = stg_ctx->time; 354ba6664aeSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 355ba6664aeSJames Wright const CeedScalar P0 = stg_ctx->P0; 356ba6664aeSJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 35746603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 35846603fc5SJames Wright const CeedScalar gamma = HeatCapacityRatio(&stg_ctx->newtonian_ctx); 359ba6664aeSJames Wright 3602b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 361ba6664aeSJames Wright const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0); 362ba6664aeSJames Wright const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 363f3e15844SJames Wright CeedScalar wdetJb, dXdx[2][3], norm[3]; 364f3e15844SJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, dXdx, norm); 365f3e15844SJames Wright wdetJb *= is_implicit ? -1. : 1.; 366ba6664aeSJames Wright 367*831dbe9eSJames Wright CeedScalar h_node_sep[3]; 368*831dbe9eSJames Wright h_node_sep[0] = dx; 369*831dbe9eSJames Wright for (CeedInt j = 1; j < 3; j++) h_node_sep[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 370*831dbe9eSJames Wright ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3); 371ba6664aeSJames Wright 372ba6664aeSJames Wright InterpolateProfile(X[1][i], ubar, cij, &eps, <, stg_ctx); 373ba6664aeSJames Wright if (!mean_only) { 374*831dbe9eSJames Wright CalcSpectrum(X[1][i], eps, lt, h_node_sep, mu / rho, qn, stg_ctx); 375cbef7084SJames Wright StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx); 376ba6664aeSJames Wright } else { 377ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 378ba6664aeSJames Wright } 379ba6664aeSJames Wright 3804dbab5e5SJames Wright const CeedScalar E_kinetic = .5 * rho * Dot3(u, u); 381ba6664aeSJames Wright CeedScalar E_internal, P; 382ba6664aeSJames Wright if (prescribe_T) { 383ba6664aeSJames Wright // Temperature is being set weakly (theta0) and for constant cv this sets E_internal 384ba6664aeSJames Wright E_internal = rho * cv * theta0; 385ba6664aeSJames Wright // Find pressure using 386ba6664aeSJames Wright P = rho * Rd * theta0; // interior rho with exterior T 387ba6664aeSJames Wright } else { 388ba6664aeSJames Wright E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution 389ba6664aeSJames Wright P = E_internal * (gamma - 1.); 390ba6664aeSJames Wright } 391ba6664aeSJames Wright 392ba6664aeSJames Wright const CeedScalar E = E_internal + E_kinetic; 393ba6664aeSJames Wright 394ba6664aeSJames Wright // Velocity normal to the boundary 3954dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, u); 3964dbab5e5SJames Wright 397ba6664aeSJames Wright // The Physics 398ba6664aeSJames Wright // Zero v so all future terms can safely sum into it 399ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 400ba6664aeSJames Wright 401ba6664aeSJames Wright // The Physics 402ba6664aeSJames Wright // -- Density 403ba6664aeSJames Wright v[0][i] -= wdetJb * rho * u_normal; 404ba6664aeSJames Wright 405ba6664aeSJames Wright // -- Momentum 4062b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 407ba6664aeSJames Wright 408ba6664aeSJames Wright // -- Total Energy Density 409ba6664aeSJames Wright v[4][i] -= wdetJb * u_normal * (E + P); 4104dbab5e5SJames Wright 411f3e15844SJames Wright const CeedScalar U[] = {rho, u[0], u[1], u[2], E}, kmstress[6] = {0.}; 412f3e15844SJames Wright StoredValuesPack(Q, i, 0, 5, U, jac_data_sur); 413f3e15844SJames Wright StoredValuesPack(Q, i, 5, 6, kmstress, jac_data_sur); 414ba6664aeSJames Wright } 415ba6664aeSJames Wright return 0; 416ba6664aeSJames Wright } 417ba6664aeSJames Wright 418cbef7084SJames Wright CEED_QFUNCTION(StgShur14Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 41946603fc5SJames Wright const CeedScalar(*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 42046603fc5SJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 42146603fc5SJames Wright const CeedScalar(*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 4224dbab5e5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 42346603fc5SJames Wright 424cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 4254dbab5e5SJames Wright const bool implicit = stg_ctx->is_implicit; 4264dbab5e5SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 42746603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 42846603fc5SJames Wright const CeedScalar gamma = HeatCapacityRatio(&stg_ctx->newtonian_ctx); 4294dbab5e5SJames Wright 4304dbab5e5SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 4314dbab5e5SJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 4324dbab5e5SJames Wright 433f0b01153SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 4344dbab5e5SJames Wright // Setup 4354dbab5e5SJames Wright // -- Interp-to-Interp q_data 4364dbab5e5SJames Wright // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 4374dbab5e5SJames Wright // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 4384dbab5e5SJames Wright // We can effect this by swapping the sign on this weight 4394dbab5e5SJames Wright const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 4404dbab5e5SJames Wright 4414dbab5e5SJames Wright // Calculate inflow values 4424dbab5e5SJames Wright CeedScalar velocity[3]; 4434dbab5e5SJames Wright for (CeedInt j = 0; j < 3; j++) velocity[j] = jac_data_sur[5 + j][i]; 444f3e15844SJames Wright // TODO This is almost certainly a bug. Velocity isn't stored here, only 0s. 4454dbab5e5SJames Wright 4464dbab5e5SJames Wright // enabling user to choose between weak T and weak rho inflow 4474dbab5e5SJames Wright CeedScalar drho, dE, dP; 4484dbab5e5SJames Wright if (prescribe_T) { 4494dbab5e5SJames Wright // rho should be from the current solution 4504dbab5e5SJames Wright drho = dq[0][i]; 4514dbab5e5SJames Wright CeedScalar dE_internal = drho * cv * theta0; 4524dbab5e5SJames Wright CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity); 4534dbab5e5SJames Wright dE = dE_internal + dE_kinetic; 4544dbab5e5SJames Wright dP = drho * Rd * theta0; // interior rho with exterior T 4554dbab5e5SJames Wright } else { // rho specified, E_internal from solution 4564dbab5e5SJames Wright drho = 0; 4574dbab5e5SJames Wright dE = dq[4][i]; 4584dbab5e5SJames Wright dP = dE * (gamma - 1.); 4594dbab5e5SJames Wright } 4602b730f8bSJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 4614dbab5e5SJames Wright 4624dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, velocity); 4634dbab5e5SJames Wright 4644dbab5e5SJames Wright v[0][i] = -wdetJb * drho * u_normal; 4652b730f8bSJeremy L Thompson for (int j = 0; j < 3; j++) v[j + 1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP); 4664dbab5e5SJames Wright v[4][i] = -wdetJb * u_normal * (dE + dP); 467f0b01153SJames Wright } 4684dbab5e5SJames Wright return 0; 4694dbab5e5SJames Wright } 4704dbab5e5SJames Wright 4710a6353c2SJames Wright /******************************************************************** 4720a6353c2SJames Wright * @brief QFunction to calculate the strongly enforce inflow BC 4730a6353c2SJames Wright * 4740a6353c2SJames Wright * This QF is for the strong application of STG via libCEED (rather than 4750a6353c2SJames Wright * through the native PETSc `DMAddBoundary` -> `bcFunc` method. 4760a6353c2SJames Wright */ 477cbef7084SJames Wright CEED_QFUNCTION(StgShur14InflowStrongQF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 47866170c20SJames Wright const CeedScalar(*dXdx_q)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0]; 47946603fc5SJames Wright const CeedScalar(*coords)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 48046603fc5SJames Wright const CeedScalar(*scale) = (const CeedScalar(*))in[2]; 481f8839eb4SJames Wright const CeedScalar(*inv_Ektotal) = (const CeedScalar(*))in[3]; 4820a6353c2SJames Wright CeedScalar(*bcval)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 4830a6353c2SJames Wright 484cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 48562e628f8SJames Wright CeedScalar u[3], ubar[3], cij[6], eps, lt; 4860a6353c2SJames Wright const bool mean_only = stg_ctx->mean_only; 4870a6353c2SJames Wright const CeedScalar dx = stg_ctx->dx; 4880a6353c2SJames Wright const CeedScalar time = stg_ctx->time; 4890a6353c2SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 4900a6353c2SJames Wright const CeedScalar P0 = stg_ctx->P0; 49146603fc5SJames Wright const CeedScalar rho = P0 / (GasConstant(&stg_ctx->newtonian_ctx) * theta0); 49246603fc5SJames Wright const CeedScalar nu = stg_ctx->newtonian_ctx.mu / rho; 4930a6353c2SJames Wright 4942b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 4950a6353c2SJames Wright const CeedScalar x[] = {coords[0][i], coords[1][i], coords[2][i]}; 4960a6353c2SJames Wright const CeedScalar dXdx[2][3] = { 49766170c20SJames Wright {dXdx_q[0][0][i], dXdx_q[0][1][i], dXdx_q[0][2][i]}, 49866170c20SJames Wright {dXdx_q[1][0][i], dXdx_q[1][1][i], dXdx_q[1][2][i]}, 4990a6353c2SJames Wright }; 5000a6353c2SJames Wright 501*831dbe9eSJames Wright CeedScalar h_node_sep[3]; 502*831dbe9eSJames Wright h_node_sep[0] = dx; 503*831dbe9eSJames Wright for (CeedInt j = 1; j < 3; j++) h_node_sep[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 504*831dbe9eSJames Wright ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3); 5050a6353c2SJames Wright 5060a6353c2SJames Wright InterpolateProfile(coords[1][i], ubar, cij, &eps, <, stg_ctx); 5070a6353c2SJames Wright if (!mean_only) { 50862e628f8SJames Wright if (1) { 509*831dbe9eSJames Wright StgShur14Calc_PrecompEktot(x, time, ubar, cij, inv_Ektotal[i], h_node_sep, x[1], eps, lt, nu, u, stg_ctx); 51062e628f8SJames Wright } else { // Original way 51162e628f8SJames Wright CeedScalar qn[STG_NMODES_MAX]; 512*831dbe9eSJames Wright CalcSpectrum(coords[1][i], eps, lt, h_node_sep, nu, qn, stg_ctx); 513cbef7084SJames Wright StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx); 51462e628f8SJames Wright } 5150a6353c2SJames Wright } else { 5160a6353c2SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 5170a6353c2SJames Wright } 5180a6353c2SJames Wright 51997baf651SJames Wright switch (stg_ctx->newtonian_ctx.state_var) { 52097baf651SJames Wright case STATEVAR_CONSERVATIVE: 5210a6353c2SJames Wright bcval[0][i] = scale[i] * rho; 5220a6353c2SJames Wright bcval[1][i] = scale[i] * rho * u[0]; 5230a6353c2SJames Wright bcval[2][i] = scale[i] * rho * u[1]; 5240a6353c2SJames Wright bcval[3][i] = scale[i] * rho * u[2]; 525cf3d54ffSJames Wright bcval[4][i] = 0.; 52697baf651SJames Wright break; 52797baf651SJames Wright 52897baf651SJames Wright case STATEVAR_PRIMITIVE: 52997baf651SJames Wright bcval[0][i] = 0; 53097baf651SJames Wright bcval[1][i] = scale[i] * u[0]; 53197baf651SJames Wright bcval[2][i] = scale[i] * u[1]; 53297baf651SJames Wright bcval[3][i] = scale[i] * u[2]; 53397baf651SJames Wright bcval[4][i] = scale[i] * theta0; 53497baf651SJames Wright break; 5350a6353c2SJames Wright } 5367c4551aaSJames Wright } 5370a6353c2SJames Wright return 0; 5380a6353c2SJames Wright } 539