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 1052b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void SpectrumConstants(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 1062b730f8bSJeremy L Thompson const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke, CeedScalar *keta, CeedScalar *kcut) { 107e159aeacSJames Wright *hmax = Max(Max(h[0], h[1]), h[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); 110175f00a6SJames Wright *kcut = M_PI / Min(Max(Max(h[1], h[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 121ba6664aeSJames Wright * @param[in] h 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 */ 1262b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void CalcSpectrum(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h[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 132175f00a6SJames Wright SpectrumConstants(wall_dist, eps, lt, h, 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 189175f00a6SJames Wright * @param[in] h 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], 197cbef7084SJames Wright const CeedScalar Ektot, const CeedScalar h[3], const CeedScalar wall_dist, const CeedScalar eps, 198cbef7084SJames Wright const CeedScalar lt, const CeedScalar nu, CeedScalar u[3], const StgShur14Context stg_ctx) { 199b277271eSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 200b277271eSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 201b277271eSJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 202b277271eSJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 203b277271eSJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 204b277271eSJames Wright CeedScalar hmax, ke, keta, kcut; 205175f00a6SJames Wright SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 206b277271eSJames Wright CeedScalar xdotd, vp[3] = {0.}; 207b277271eSJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 208b277271eSJames Wright 2092b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 210b277271eSJames Wright xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 211b277271eSJames Wright xdotd = 0.; 212b277271eSJames Wright for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 213b277271eSJames Wright const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 214b277271eSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 215b277271eSJames Wright const CeedScalar qn = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot); 216b277271eSJames Wright vp[0] += sqrt(qn) * sigma[0 * nmodes + n] * cos_kxdp; 217b277271eSJames Wright vp[1] += sqrt(qn) * sigma[1 * nmodes + n] * cos_kxdp; 218b277271eSJames Wright vp[2] += sqrt(qn) * sigma[2 * nmodes + n] * cos_kxdp; 219b277271eSJames Wright } 220b277271eSJames Wright for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 221b277271eSJames Wright 222b277271eSJames Wright u[0] = ubar[0] + cij[0] * vp[0]; 223b277271eSJames Wright u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 224b277271eSJames Wright u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 225b277271eSJames Wright } 226b277271eSJames Wright 22762e628f8SJames Wright // Create preprocessed input for the stg calculation 22862e628f8SJames Wright // 22962e628f8SJames Wright // stg_data[0] = 1 / Ektot (inverse of total spectrum energy) 230cbef7084SJames Wright CEED_QFUNCTION(StgShur14Preprocess)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 23166170c20SJames Wright const CeedScalar(*dXdx_q)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0]; 23246603fc5SJames Wright const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 233b277271eSJames Wright 234b277271eSJames Wright CeedScalar(*stg_data) = (CeedScalar(*))out[0]; 235b277271eSJames Wright 236b277271eSJames Wright CeedScalar ubar[3], cij[6], eps, lt; 237cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 238b277271eSJames Wright const CeedScalar dx = stg_ctx->dx; 239b277271eSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 240b277271eSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 241b277271eSJames Wright const CeedScalar P0 = stg_ctx->P0; 24246603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 243b277271eSJames Wright const CeedScalar rho = P0 / (Rd * theta0); 244b277271eSJames Wright const CeedScalar nu = mu / rho; 245b277271eSJames Wright 246b277271eSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 247b277271eSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 2485dc40723SJames Wright CeedScalar hmax, ke, keta, kcut; 249b277271eSJames Wright 2502b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 251175f00a6SJames Wright const CeedScalar wall_dist = x[1][i]; 252b277271eSJames Wright const CeedScalar dXdx[2][3] = { 25366170c20SJames Wright {dXdx_q[0][0][i], dXdx_q[0][1][i], dXdx_q[0][2][i]}, 25466170c20SJames Wright {dXdx_q[1][0][i], dXdx_q[1][1][i], dXdx_q[1][2][i]}, 255b277271eSJames Wright }; 256b277271eSJames Wright 257b277271eSJames Wright CeedScalar h[3]; 258b277271eSJames Wright h[0] = dx; 2592b730f8bSJeremy 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]); 260b277271eSJames Wright 261175f00a6SJames Wright InterpolateProfile(wall_dist, ubar, cij, &eps, <, stg_ctx); 262175f00a6SJames Wright SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 263b277271eSJames Wright 264b277271eSJames Wright // Calculate total TKE per spectrum 265d97dc904SJames Wright CeedScalar Ek_tot = 0; 2662b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 267b277271eSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 268d97dc904SJames Wright Ek_tot += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 269b277271eSJames Wright } 270d97dc904SJames Wright // avoid underflowed and poorly defined spectrum coefficients 271d97dc904SJames Wright stg_data[i] = Ek_tot != 0 ? 1 / Ek_tot : 0; 272b277271eSJames Wright } 273b277271eSJames Wright return 0; 274b277271eSJames Wright } 275b277271eSJames Wright 276b77c53c9SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition 277cbef7084SJames Wright CEED_QFUNCTION(ICsStg)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 27846603fc5SJames Wright const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 2791c299e57SJeremy L Thompson const CeedScalar(*J)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[1]; 280b77c53c9SJames Wright CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 281b77c53c9SJames Wright 282cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 283*a2d72b6fSJames Wright const NewtonianIdealGasContext gas = &stg_ctx->newtonian_ctx; 28489060322SJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 28589060322SJames Wright const CeedScalar dx = stg_ctx->dx; 28689060322SJames Wright const CeedScalar time = stg_ctx->time; 287b77c53c9SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 288b77c53c9SJames Wright const CeedScalar P0 = stg_ctx->P0; 289*a2d72b6fSJames Wright const CeedScalar rho = P0 / (GasConstant(gas) * theta0); 290*a2d72b6fSJames Wright const CeedScalar nu = gas->mu / rho; 291b77c53c9SJames Wright 2922b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 29389060322SJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 2941c299e57SJeremy L Thompson CeedScalar dXdx[3][3]; 2958756a6ccSJames Wright InvertMappingJacobian_3D(Q, i, J, dXdx, NULL); 29689060322SJames Wright CeedScalar h[3]; 29789060322SJames Wright h[0] = dx; 2982b730f8bSJeremy 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])); 29989060322SJames Wright 30089060322SJames Wright InterpolateProfile(x_i[1], ubar, cij, &eps, <, stg_ctx); 30189060322SJames Wright if (stg_ctx->use_fluctuating_IC) { 30246603fc5SJames Wright CalcSpectrum(x_i[1], eps, lt, h, nu, qn, stg_ctx); 303cbef7084SJames Wright StgShur14Calc(x_i, time, ubar, cij, qn, u, stg_ctx); 30489060322SJames Wright } else { 30589060322SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 30689060322SJames Wright } 307b77c53c9SJames Wright 308*a2d72b6fSJames Wright CeedScalar Y[5] = {P0, u[0], u[1], u[2], theta0}, q[5]; 309*a2d72b6fSJames Wright State s = StateFromY(gas, Y); 310*a2d72b6fSJames Wright StateToQ(gas, s, q, gas->state_var); 311*a2d72b6fSJames Wright for (CeedInt j = 0; j < 5; j++) { 312*a2d72b6fSJames Wright q0[j][i] = q[j]; 3137c4551aaSJames Wright } 314f0b01153SJames Wright } 315b77c53c9SJames Wright return 0; 316b77c53c9SJames Wright } 317b77c53c9SJames Wright 318ba6664aeSJames Wright /******************************************************************** 319ba6664aeSJames Wright * @brief QFunction to calculate the inflow boundary condition 320ba6664aeSJames Wright * 321ba6664aeSJames Wright * This will loop through quadrature points, calculate the wavemode amplitudes 322ba6664aeSJames Wright * at each location, then calculate the actual velocity. 323ba6664aeSJames Wright */ 324cbef7084SJames Wright CEED_QFUNCTION(StgShur14Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 32546603fc5SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 326f3e15844SJames Wright const CeedScalar(*q_data_sur) = in[2]; 32746603fc5SJames Wright const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 328ba6664aeSJames Wright 32946603fc5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 330f3e15844SJames Wright CeedScalar(*jac_data_sur) = out[1]; 331ba6664aeSJames Wright 332cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 333ba6664aeSJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 334ba6664aeSJames Wright const bool is_implicit = stg_ctx->is_implicit; 335ba6664aeSJames Wright const bool mean_only = stg_ctx->mean_only; 336ba6664aeSJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 337ba6664aeSJames Wright const CeedScalar dx = stg_ctx->dx; 338ba6664aeSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 339ba6664aeSJames Wright const CeedScalar time = stg_ctx->time; 340ba6664aeSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 341ba6664aeSJames Wright const CeedScalar P0 = stg_ctx->P0; 342ba6664aeSJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 34346603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 34446603fc5SJames Wright const CeedScalar gamma = HeatCapacityRatio(&stg_ctx->newtonian_ctx); 345ba6664aeSJames Wright 3462b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 347ba6664aeSJames Wright const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0); 348ba6664aeSJames Wright const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 349f3e15844SJames Wright CeedScalar wdetJb, dXdx[2][3], norm[3]; 350f3e15844SJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, dXdx, norm); 351f3e15844SJames Wright wdetJb *= is_implicit ? -1. : 1.; 352ba6664aeSJames Wright 353ba6664aeSJames Wright CeedScalar h[3]; 354ba6664aeSJames Wright h[0] = dx; 3552b730f8bSJeremy L Thompson for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 356ba6664aeSJames Wright 357ba6664aeSJames Wright InterpolateProfile(X[1][i], ubar, cij, &eps, <, stg_ctx); 358ba6664aeSJames Wright if (!mean_only) { 359ba6664aeSJames Wright CalcSpectrum(X[1][i], eps, lt, h, mu / rho, qn, stg_ctx); 360cbef7084SJames Wright StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx); 361ba6664aeSJames Wright } else { 362ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 363ba6664aeSJames Wright } 364ba6664aeSJames Wright 3654dbab5e5SJames Wright const CeedScalar E_kinetic = .5 * rho * Dot3(u, u); 366ba6664aeSJames Wright CeedScalar E_internal, P; 367ba6664aeSJames Wright if (prescribe_T) { 368ba6664aeSJames Wright // Temperature is being set weakly (theta0) and for constant cv this sets E_internal 369ba6664aeSJames Wright E_internal = rho * cv * theta0; 370ba6664aeSJames Wright // Find pressure using 371ba6664aeSJames Wright P = rho * Rd * theta0; // interior rho with exterior T 372ba6664aeSJames Wright } else { 373ba6664aeSJames Wright E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution 374ba6664aeSJames Wright P = E_internal * (gamma - 1.); 375ba6664aeSJames Wright } 376ba6664aeSJames Wright 377ba6664aeSJames Wright const CeedScalar E = E_internal + E_kinetic; 378ba6664aeSJames Wright 379ba6664aeSJames Wright // Velocity normal to the boundary 3804dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, u); 3814dbab5e5SJames Wright 382ba6664aeSJames Wright // The Physics 383ba6664aeSJames Wright // Zero v so all future terms can safely sum into it 384ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 385ba6664aeSJames Wright 386ba6664aeSJames Wright // The Physics 387ba6664aeSJames Wright // -- Density 388ba6664aeSJames Wright v[0][i] -= wdetJb * rho * u_normal; 389ba6664aeSJames Wright 390ba6664aeSJames Wright // -- Momentum 3912b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 392ba6664aeSJames Wright 393ba6664aeSJames Wright // -- Total Energy Density 394ba6664aeSJames Wright v[4][i] -= wdetJb * u_normal * (E + P); 3954dbab5e5SJames Wright 396f3e15844SJames Wright const CeedScalar U[] = {rho, u[0], u[1], u[2], E}, kmstress[6] = {0.}; 397f3e15844SJames Wright StoredValuesPack(Q, i, 0, 5, U, jac_data_sur); 398f3e15844SJames Wright StoredValuesPack(Q, i, 5, 6, kmstress, jac_data_sur); 399ba6664aeSJames Wright } 400ba6664aeSJames Wright return 0; 401ba6664aeSJames Wright } 402ba6664aeSJames Wright 403cbef7084SJames Wright CEED_QFUNCTION(StgShur14Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 40446603fc5SJames Wright const CeedScalar(*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 40546603fc5SJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 40646603fc5SJames Wright const CeedScalar(*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 4074dbab5e5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 40846603fc5SJames Wright 409cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 4104dbab5e5SJames Wright const bool implicit = stg_ctx->is_implicit; 4114dbab5e5SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 41246603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 41346603fc5SJames Wright const CeedScalar gamma = HeatCapacityRatio(&stg_ctx->newtonian_ctx); 4144dbab5e5SJames Wright 4154dbab5e5SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 4164dbab5e5SJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 4174dbab5e5SJames Wright 418f0b01153SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 4194dbab5e5SJames Wright // Setup 4204dbab5e5SJames Wright // -- Interp-to-Interp q_data 4214dbab5e5SJames Wright // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 4224dbab5e5SJames Wright // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 4234dbab5e5SJames Wright // We can effect this by swapping the sign on this weight 4244dbab5e5SJames Wright const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 4254dbab5e5SJames Wright 4264dbab5e5SJames Wright // Calculate inflow values 4274dbab5e5SJames Wright CeedScalar velocity[3]; 4284dbab5e5SJames Wright for (CeedInt j = 0; j < 3; j++) velocity[j] = jac_data_sur[5 + j][i]; 429f3e15844SJames Wright // TODO This is almost certainly a bug. Velocity isn't stored here, only 0s. 4304dbab5e5SJames Wright 4314dbab5e5SJames Wright // enabling user to choose between weak T and weak rho inflow 4324dbab5e5SJames Wright CeedScalar drho, dE, dP; 4334dbab5e5SJames Wright if (prescribe_T) { 4344dbab5e5SJames Wright // rho should be from the current solution 4354dbab5e5SJames Wright drho = dq[0][i]; 4364dbab5e5SJames Wright CeedScalar dE_internal = drho * cv * theta0; 4374dbab5e5SJames Wright CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity); 4384dbab5e5SJames Wright dE = dE_internal + dE_kinetic; 4394dbab5e5SJames Wright dP = drho * Rd * theta0; // interior rho with exterior T 4404dbab5e5SJames Wright } else { // rho specified, E_internal from solution 4414dbab5e5SJames Wright drho = 0; 4424dbab5e5SJames Wright dE = dq[4][i]; 4434dbab5e5SJames Wright dP = dE * (gamma - 1.); 4444dbab5e5SJames Wright } 4452b730f8bSJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 4464dbab5e5SJames Wright 4474dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, velocity); 4484dbab5e5SJames Wright 4494dbab5e5SJames Wright v[0][i] = -wdetJb * drho * u_normal; 4502b730f8bSJeremy L Thompson for (int j = 0; j < 3; j++) v[j + 1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP); 4514dbab5e5SJames Wright v[4][i] = -wdetJb * u_normal * (dE + dP); 452f0b01153SJames Wright } 4534dbab5e5SJames Wright return 0; 4544dbab5e5SJames Wright } 4554dbab5e5SJames Wright 4560a6353c2SJames Wright /******************************************************************** 4570a6353c2SJames Wright * @brief QFunction to calculate the strongly enforce inflow BC 4580a6353c2SJames Wright * 4590a6353c2SJames Wright * This QF is for the strong application of STG via libCEED (rather than 4600a6353c2SJames Wright * through the native PETSc `DMAddBoundary` -> `bcFunc` method. 4610a6353c2SJames Wright */ 462cbef7084SJames Wright CEED_QFUNCTION(StgShur14InflowStrongQF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 46366170c20SJames Wright const CeedScalar(*dXdx_q)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0]; 46446603fc5SJames Wright const CeedScalar(*coords)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 46546603fc5SJames Wright const CeedScalar(*scale) = (const CeedScalar(*))in[2]; 466f8839eb4SJames Wright const CeedScalar(*inv_Ektotal) = (const CeedScalar(*))in[3]; 4670a6353c2SJames Wright CeedScalar(*bcval)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 4680a6353c2SJames Wright 469cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 470*a2d72b6fSJames Wright const NewtonianIdealGasContext gas = &stg_ctx->newtonian_ctx; 47162e628f8SJames Wright CeedScalar u[3], ubar[3], cij[6], eps, lt; 4720a6353c2SJames Wright const bool mean_only = stg_ctx->mean_only; 4730a6353c2SJames Wright const CeedScalar dx = stg_ctx->dx; 4740a6353c2SJames Wright const CeedScalar time = stg_ctx->time; 4750a6353c2SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 4760a6353c2SJames Wright const CeedScalar P0 = stg_ctx->P0; 477*a2d72b6fSJames Wright const CeedScalar rho = P0 / (GasConstant(gas) * theta0); 478*a2d72b6fSJames Wright const CeedScalar nu = gas->mu / rho; 4790a6353c2SJames Wright 4802b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 4810a6353c2SJames Wright const CeedScalar x[] = {coords[0][i], coords[1][i], coords[2][i]}; 4820a6353c2SJames Wright const CeedScalar dXdx[2][3] = { 48366170c20SJames Wright {dXdx_q[0][0][i], dXdx_q[0][1][i], dXdx_q[0][2][i]}, 48466170c20SJames Wright {dXdx_q[1][0][i], dXdx_q[1][1][i], dXdx_q[1][2][i]}, 4850a6353c2SJames Wright }; 4860a6353c2SJames Wright 4870a6353c2SJames Wright CeedScalar h[3]; 4880a6353c2SJames Wright h[0] = dx; 4892b730f8bSJeremy L Thompson for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 4900a6353c2SJames Wright 4910a6353c2SJames Wright InterpolateProfile(coords[1][i], ubar, cij, &eps, <, stg_ctx); 4920a6353c2SJames Wright if (!mean_only) { 49362e628f8SJames Wright if (1) { 494cbef7084SJames Wright StgShur14Calc_PrecompEktot(x, time, ubar, cij, inv_Ektotal[i], h, x[1], eps, lt, nu, u, stg_ctx); 49562e628f8SJames Wright } else { // Original way 49662e628f8SJames Wright CeedScalar qn[STG_NMODES_MAX]; 49746603fc5SJames Wright CalcSpectrum(coords[1][i], eps, lt, h, nu, qn, stg_ctx); 498cbef7084SJames Wright StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx); 49962e628f8SJames Wright } 5000a6353c2SJames Wright } else { 5010a6353c2SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 5020a6353c2SJames Wright } 5030a6353c2SJames Wright 504*a2d72b6fSJames Wright CeedScalar Y[5] = {P0, u[0], u[1], u[2], theta0}, q[5]; 505*a2d72b6fSJames Wright State s = StateFromY(gas, Y); 506*a2d72b6fSJames Wright StateToQ(gas, s, q, gas->state_var); 507*a2d72b6fSJames Wright switch (gas->state_var) { 50897baf651SJames Wright case STATEVAR_CONSERVATIVE: 509*a2d72b6fSJames Wright q[4] = 0.; // Don't set energy 51097baf651SJames Wright break; 51197baf651SJames Wright case STATEVAR_PRIMITIVE: 512*a2d72b6fSJames Wright q[0] = 0; // Don't set pressure 51397baf651SJames Wright break; 514*a2d72b6fSJames Wright case STATEVAR_ENTROPY: 515*a2d72b6fSJames Wright q[0] = 0; // Don't set V_density 516*a2d72b6fSJames Wright break; 517*a2d72b6fSJames Wright } 518*a2d72b6fSJames Wright for (CeedInt j = 0; j < 5; j++) { 519*a2d72b6fSJames Wright bcval[j][i] = scale[i] * q[j]; 5200a6353c2SJames Wright } 5217c4551aaSJames Wright } 5220a6353c2SJames Wright return 0; 5230a6353c2SJames Wright } 524