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> 22*2b730f8bSJeremy L Thompson 23ba6664aeSJames Wright #include "stg_shur14_type.h" 2413fa47b2SJames Wright #include "utils.h" 25ba6664aeSJames Wright 26ba6664aeSJames Wright #define STG_NMODES_MAX 1024 27ba6664aeSJames Wright 28ba6664aeSJames Wright /* 29ba6664aeSJames Wright * @brief Interpolate quantities from input profile to given location 30ba6664aeSJames Wright * 31175f00a6SJames Wright * Assumed that prof_wd[i+1] > prof_wd[i] and prof_wd[0] = 0 32175f00a6SJames Wright * If wall_dist > prof_wd[-1], then the interpolation takes the values at prof_wd[-1] 33ba6664aeSJames Wright * 34175f00a6SJames Wright * @param[in] wall_dist Distance to the nearest wall 35175f00a6SJames Wright * @param[out] ubar Mean velocity at wall_dist 36175f00a6SJames Wright * @param[out] cij Cholesky decomposition at wall_dist 37175f00a6SJames Wright * @param[out] eps Turbulent dissipation at wall_dist 38175f00a6SJames Wright * @param[out] lt Turbulent length scale at wall_dist 39ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 40ba6664aeSJames Wright */ 41*2b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar wall_dist, CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt, 42ba6664aeSJames Wright const STGShur14Context stg_ctx) { 43ba6664aeSJames Wright const CeedInt nprofs = stg_ctx->nprofs; 44175f00a6SJames Wright const CeedScalar *prof_wd = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 45ba6664aeSJames Wright const CeedScalar *prof_eps = &stg_ctx->data[stg_ctx->offsets.eps]; 46ba6664aeSJames Wright const CeedScalar *prof_lt = &stg_ctx->data[stg_ctx->offsets.lt]; 47ba6664aeSJames Wright const CeedScalar *prof_ubar = &stg_ctx->data[stg_ctx->offsets.ubar]; 48ba6664aeSJames Wright const CeedScalar *prof_cij = &stg_ctx->data[stg_ctx->offsets.cij]; 49ba6664aeSJames Wright CeedInt idx = -1; 50ba6664aeSJames Wright 51ba6664aeSJames Wright for (CeedInt i = 0; i < nprofs; i++) { 52175f00a6SJames Wright if (wall_dist < prof_wd[i]) { 53ba6664aeSJames Wright idx = i; 54ba6664aeSJames Wright break; 55ba6664aeSJames Wright } 56ba6664aeSJames Wright } 57ba6664aeSJames Wright 58175f00a6SJames Wright if (idx > 0) { // y within the bounds of prof_wd 59ba6664aeSJames Wright //*INDENT-OFF* 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]); 73ba6664aeSJames Wright //*INDENT-ON* 74175f00a6SJames Wright } else { // y outside bounds of prof_wd 75ba6664aeSJames Wright ubar[0] = prof_ubar[1 * nprofs - 1]; 76ba6664aeSJames Wright ubar[1] = prof_ubar[2 * nprofs - 1]; 77ba6664aeSJames Wright ubar[2] = prof_ubar[3 * nprofs - 1]; 78ba6664aeSJames Wright cij[0] = prof_cij[1 * nprofs - 1]; 79ba6664aeSJames Wright cij[1] = prof_cij[2 * nprofs - 1]; 80ba6664aeSJames Wright cij[2] = prof_cij[3 * nprofs - 1]; 81ba6664aeSJames Wright cij[3] = prof_cij[4 * nprofs - 1]; 82ba6664aeSJames Wright cij[4] = prof_cij[5 * nprofs - 1]; 83ba6664aeSJames Wright cij[5] = prof_cij[6 * nprofs - 1]; 84ba6664aeSJames Wright *eps = prof_eps[nprofs - 1]; 85ba6664aeSJames Wright *lt = prof_lt[nprofs - 1]; 86ba6664aeSJames Wright } 87ba6664aeSJames Wright } 88ba6664aeSJames Wright 89ba6664aeSJames Wright /* 90e159aeacSJames Wright * @brief Calculate spectrum coefficient, qn 91e159aeacSJames Wright * 92e159aeacSJames Wright * Calculates q_n at a given distance to the wall 93e159aeacSJames Wright * 94e159aeacSJames Wright * @param[in] kappa nth wavenumber 95e159aeacSJames Wright * @param[in] dkappa Difference between wavenumbers 96e159aeacSJames Wright * @param[in] keta Dissipation wavenumber 97e159aeacSJames Wright * @param[in] kcut Mesh-induced cutoff wavenumber 98e159aeacSJames Wright * @param[in] ke Energy-containing wavenumber 99e159aeacSJames Wright * @param[in] Ektot Total turbulent kinetic energy of spectrum 100e159aeacSJames Wright * @returns qn Spectrum coefficient 101e159aeacSJames Wright */ 102*2b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER CeedScalar Calc_qn(const CeedScalar kappa, const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut, 10362e628f8SJames Wright const CeedScalar ke, const CeedScalar Ektot_inv) { 104*2b730f8bSJeremy L Thompson const CeedScalar feta_x_fcut = exp(-Square(12 * kappa / keta) - Cube(4 * Max(kappa - 0.9 * kcut, 0) / kcut)); 105*2b730f8bSJeremy L Thompson return pow(kappa / ke, 4.) * pow(1 + 2.4 * Square(kappa / ke), -17. / 6) * feta_x_fcut * dkappa * Ektot_inv; 106e159aeacSJames Wright } 107e159aeacSJames Wright 108e159aeacSJames Wright // Calculate hmax, ke, keta, and kcut 109*2b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void SpectrumConstants(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 110*2b730f8bSJeremy L Thompson const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke, CeedScalar *keta, CeedScalar *kcut) { 111e159aeacSJames Wright *hmax = Max(Max(h[0], h[1]), h[2]); 112175f00a6SJames Wright *ke = wall_dist == 0 ? 1e16 : 2 * M_PI / Min(2 * wall_dist, 3 * lt); 113e159aeacSJames Wright *keta = 2 * M_PI * pow(Cube(nu) / eps, -0.25); 114175f00a6SJames Wright *kcut = M_PI / Min(Max(Max(h[1], h[2]), 0.3 * (*hmax)) + 0.1 * wall_dist, *hmax); 115e159aeacSJames Wright } 116e159aeacSJames Wright 117e159aeacSJames Wright /* 118ba6664aeSJames Wright * @brief Calculate spectrum coefficients for STG 119ba6664aeSJames Wright * 120ba6664aeSJames Wright * Calculates q_n at a given distance to the wall 121ba6664aeSJames Wright * 122175f00a6SJames Wright * @param[in] wall_dist Distance to the nearest wall 123175f00a6SJames Wright * @param[in] eps Turbulent dissipation w/rt wall_dist 124175f00a6SJames Wright * @param[in] lt Turbulent length scale w/rt wall_dist 125ba6664aeSJames Wright * @param[in] h Element lengths in coordinate directions 126ba6664aeSJames Wright * @param[in] nu Dynamic Viscosity; 127ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 128ba6664aeSJames Wright * @param[out] qn Spectrum coefficients, [nmodes] 129ba6664aeSJames Wright */ 130*2b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void CalcSpectrum(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 131ba6664aeSJames Wright const CeedScalar nu, CeedScalar qn[], const STGShur14Context stg_ctx) { 132ba6664aeSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 133ba6664aeSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 134e159aeacSJames Wright CeedScalar hmax, ke, keta, kcut, Ektot = 0.0; 135*2b730f8bSJeremy L Thompson 136175f00a6SJames Wright SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 137ba6664aeSJames Wright 138ba6664aeSJames Wright for (CeedInt n = 0; n < nmodes; n++) { 139e159aeacSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 140e159aeacSJames Wright qn[n] = Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 141ba6664aeSJames Wright Ektot += qn[n]; 142ba6664aeSJames Wright } 143ba6664aeSJames Wright 144961c9c98SJames Wright if (Ektot == 0) return; 145ba6664aeSJames Wright for (CeedInt n = 0; n < nmodes; n++) qn[n] /= Ektot; 146ba6664aeSJames Wright } 147ba6664aeSJames Wright 148ba6664aeSJames Wright /****************************************************** 149ba6664aeSJames Wright * @brief Calculate u(x,t) for STG inflow condition 150ba6664aeSJames Wright * 151ba6664aeSJames Wright * @param[in] X Location to evaluate u(X,t) 152ba6664aeSJames Wright * @param[in] t Time to evaluate u(X,t) 153ba6664aeSJames Wright * @param[in] ubar Mean velocity at X 154ba6664aeSJames Wright * @param[in] cij Cholesky decomposition at X 155ba6664aeSJames Wright * @param[in] qn Wavemode amplitudes at X, [nmodes] 156ba6664aeSJames Wright * @param[out] u Velocity at X and t 157ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 158ba6664aeSJames Wright */ 159*2b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void STGShur14_Calc(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 160*2b730f8bSJeremy L Thompson const CeedScalar qn[], CeedScalar u[3], const STGShur14Context stg_ctx) { 161ba6664aeSJames Wright //*INDENT-OFF* 162ba6664aeSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 163ba6664aeSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 164ba6664aeSJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 165ba6664aeSJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 166ba6664aeSJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 167ba6664aeSJames Wright //*INDENT-ON* 168ba6664aeSJames Wright CeedScalar xdotd, vp[3] = {0.}; 169ba6664aeSJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 170ba6664aeSJames Wright 171*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 172ba6664aeSJames Wright xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 173ba6664aeSJames Wright xdotd = 0.; 174ba6664aeSJames Wright for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 175ba6664aeSJames Wright const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 176961c9c98SJames Wright vp[0] += sqrt(qn[n]) * sigma[0 * nmodes + n] * cos_kxdp; 177961c9c98SJames Wright vp[1] += sqrt(qn[n]) * sigma[1 * nmodes + n] * cos_kxdp; 178961c9c98SJames Wright vp[2] += sqrt(qn[n]) * sigma[2 * nmodes + n] * cos_kxdp; 179ba6664aeSJames Wright } 180961c9c98SJames Wright for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 181ba6664aeSJames Wright 182ba6664aeSJames Wright u[0] = ubar[0] + cij[0] * vp[0]; 183ba6664aeSJames Wright u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 184ba6664aeSJames Wright u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 185ba6664aeSJames Wright } 186ba6664aeSJames Wright 187b277271eSJames Wright /****************************************************** 188b277271eSJames Wright * @brief Calculate u(x,t) for STG inflow condition 189b277271eSJames Wright * 190b277271eSJames Wright * @param[in] X Location to evaluate u(X,t) 191b277271eSJames Wright * @param[in] t Time to evaluate u(X,t) 192b277271eSJames Wright * @param[in] ubar Mean velocity at X 193b277271eSJames Wright * @param[in] cij Cholesky decomposition at X 194175f00a6SJames Wright * @param[in] Ektot Total spectrum energy at this location 195175f00a6SJames Wright * @param[in] h Element size in 3 directions 196175f00a6SJames Wright * @param[in] wall_dist Distance to closest wall 197175f00a6SJames Wright * @param[in] eps Turbulent dissipation 198175f00a6SJames Wright * @param[in] lt Turbulent length scale 199b277271eSJames Wright * @param[out] u Velocity at X and t 200b277271eSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 201b277271eSJames Wright */ 202*2b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void STGShur14_Calc_PrecompEktot(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 203175f00a6SJames Wright const CeedScalar Ektot, const CeedScalar h[3], const CeedScalar wall_dist, 204b277271eSJames Wright const CeedScalar eps, const CeedScalar lt, const CeedScalar nu, CeedScalar u[3], 205b277271eSJames Wright const STGShur14Context stg_ctx) { 206b277271eSJames Wright //*INDENT-OFF* 207b277271eSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 208b277271eSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 209b277271eSJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 210b277271eSJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 211b277271eSJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 212b277271eSJames Wright //*INDENT-ON* 213b277271eSJames Wright CeedScalar hmax, ke, keta, kcut; 214175f00a6SJames Wright SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 215b277271eSJames Wright CeedScalar xdotd, vp[3] = {0.}; 216b277271eSJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 217b277271eSJames Wright 218*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 219b277271eSJames Wright xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 220b277271eSJames Wright xdotd = 0.; 221b277271eSJames Wright for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 222b277271eSJames Wright const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 223b277271eSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 224b277271eSJames Wright const CeedScalar qn = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot); 225b277271eSJames Wright vp[0] += sqrt(qn) * sigma[0 * nmodes + n] * cos_kxdp; 226b277271eSJames Wright vp[1] += sqrt(qn) * sigma[1 * nmodes + n] * cos_kxdp; 227b277271eSJames Wright vp[2] += sqrt(qn) * sigma[2 * nmodes + n] * cos_kxdp; 228b277271eSJames Wright } 229b277271eSJames Wright for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 230b277271eSJames Wright 231b277271eSJames Wright u[0] = ubar[0] + cij[0] * vp[0]; 232b277271eSJames Wright u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 233b277271eSJames Wright u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 234b277271eSJames Wright } 235b277271eSJames Wright 23662e628f8SJames Wright // Create preprocessed input for the stg calculation 23762e628f8SJames Wright // 23862e628f8SJames Wright // stg_data[0] = 1 / Ektot (inverse of total spectrum energy) 239*2b730f8bSJeremy L Thompson CEED_QFUNCTION(Preprocess_STGShur14)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 240b277271eSJames Wright //*INDENT-OFF* 241*2b730f8bSJeremy L Thompson const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 242b277271eSJames Wright 243b277271eSJames Wright CeedScalar(*stg_data) = (CeedScalar(*))out[0]; 244b277271eSJames Wright 245b277271eSJames Wright //*INDENT-ON* 246b277271eSJames Wright CeedScalar ubar[3], cij[6], eps, lt; 247b277271eSJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 248b277271eSJames Wright const CeedScalar dx = stg_ctx->dx; 249b277271eSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 250b277271eSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 251b277271eSJames Wright const CeedScalar P0 = stg_ctx->P0; 252b277271eSJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 253b277271eSJames Wright const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 254b277271eSJames Wright const CeedScalar Rd = cp - cv; 255b277271eSJames Wright const CeedScalar rho = P0 / (Rd * theta0); 256b277271eSJames Wright const CeedScalar nu = mu / rho; 257b277271eSJames Wright 258b277271eSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 259b277271eSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 2605dc40723SJames Wright CeedScalar hmax, ke, keta, kcut; 261b277271eSJames Wright 262*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 263175f00a6SJames Wright const CeedScalar wall_dist = x[1][i]; 264b277271eSJames Wright const CeedScalar dXdx[2][3] = { 265b277271eSJames Wright {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 266b277271eSJames Wright {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 267b277271eSJames Wright }; 268b277271eSJames Wright 269b277271eSJames Wright CeedScalar h[3]; 270b277271eSJames Wright h[0] = dx; 271*2b730f8bSJeremy 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]); 272b277271eSJames Wright 273175f00a6SJames Wright InterpolateProfile(wall_dist, ubar, cij, &eps, <, stg_ctx); 274175f00a6SJames Wright SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 275b277271eSJames Wright 276b277271eSJames Wright // Calculate total TKE per spectrum 277d97dc904SJames Wright CeedScalar Ek_tot = 0; 278*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 279b277271eSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 280d97dc904SJames Wright Ek_tot += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 281b277271eSJames Wright } 282d97dc904SJames Wright // avoid underflowed and poorly defined spectrum coefficients 283d97dc904SJames Wright stg_data[i] = Ek_tot != 0 ? 1 / Ek_tot : 0; 284b277271eSJames Wright } 285b277271eSJames Wright return 0; 286b277271eSJames Wright } 287b277271eSJames Wright 288b77c53c9SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition 289*2b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsSTG)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 290b77c53c9SJames Wright // Inputs 291*2b730f8bSJeremy L Thompson const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 292b77c53c9SJames Wright // Outputs 293b77c53c9SJames Wright CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 294b77c53c9SJames Wright 295b77c53c9SJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 29689060322SJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 29789060322SJames Wright const CeedScalar dx = stg_ctx->dx; 29889060322SJames Wright const CeedScalar time = stg_ctx->time; 299b77c53c9SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 300b77c53c9SJames Wright const CeedScalar P0 = stg_ctx->P0; 30189060322SJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 302b77c53c9SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 303b77c53c9SJames Wright const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 304b77c53c9SJames Wright const CeedScalar Rd = cp - cv; 305b77c53c9SJames Wright const CeedScalar rho = P0 / (Rd * theta0); 306b77c53c9SJames Wright 307*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 30889060322SJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 30989060322SJames Wright // *INDENT-OFF* 310*2b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 311*2b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 31289060322SJames Wright {q_data[4][i], q_data[5][i], q_data[6][i]}, 31389060322SJames Wright {q_data[7][i], q_data[8][i], q_data[9][i]} 31489060322SJames Wright }; 31589060322SJames Wright // *INDENT-ON* 31689060322SJames Wright 31789060322SJames Wright CeedScalar h[3]; 31889060322SJames Wright h[0] = dx; 319*2b730f8bSJeremy 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])); 32089060322SJames Wright 32189060322SJames Wright InterpolateProfile(x_i[1], ubar, cij, &eps, <, stg_ctx); 32289060322SJames Wright if (stg_ctx->use_fluctuating_IC) { 32389060322SJames Wright CalcSpectrum(x_i[1], eps, lt, h, mu / rho, qn, stg_ctx); 32489060322SJames Wright STGShur14_Calc(x_i, time, ubar, cij, qn, u, stg_ctx); 32589060322SJames Wright } else { 32689060322SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 32789060322SJames Wright } 328b77c53c9SJames Wright 32997baf651SJames Wright switch (stg_ctx->newtonian_ctx.state_var) { 33097baf651SJames Wright case STATEVAR_CONSERVATIVE: 331b77c53c9SJames Wright q0[0][i] = rho; 332b77c53c9SJames Wright q0[1][i] = u[0] * rho; 333b77c53c9SJames Wright q0[2][i] = u[1] * rho; 334b77c53c9SJames Wright q0[3][i] = u[2] * rho; 335b77c53c9SJames Wright q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0); 33697baf651SJames Wright break; 33797baf651SJames Wright 33897baf651SJames Wright case STATEVAR_PRIMITIVE: 33997baf651SJames Wright q0[0][i] = P0; 34097baf651SJames Wright q0[1][i] = u[0]; 34197baf651SJames Wright q0[2][i] = u[1]; 34297baf651SJames Wright q0[3][i] = u[2]; 34397baf651SJames Wright q0[4][i] = theta0; 34497baf651SJames Wright break; 3457c4551aaSJames Wright } 346b77c53c9SJames Wright } // End of Quadrature Point Loop 347b77c53c9SJames Wright return 0; 348b77c53c9SJames Wright } 349b77c53c9SJames Wright 350ba6664aeSJames Wright /******************************************************************** 351ba6664aeSJames Wright * @brief QFunction to calculate the inflow boundary condition 352ba6664aeSJames Wright * 353ba6664aeSJames Wright * This will loop through quadrature points, calculate the wavemode amplitudes 354ba6664aeSJames Wright * at each location, then calculate the actual velocity. 355ba6664aeSJames Wright */ 356*2b730f8bSJeremy L Thompson CEED_QFUNCTION(STGShur14_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 357ba6664aeSJames Wright //*INDENT-OFF* 358*2b730f8bSJeremy L Thompson const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 359e8b03feeSJames Wright (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 360ba6664aeSJames Wright 361*2b730f8bSJeremy L Thompson CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], (*jac_data_sur)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[1]; 362ba6664aeSJames Wright 363ba6664aeSJames Wright //*INDENT-ON* 364ba6664aeSJames Wright 365ba6664aeSJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 366ba6664aeSJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 367ba6664aeSJames Wright const bool is_implicit = stg_ctx->is_implicit; 368ba6664aeSJames Wright const bool mean_only = stg_ctx->mean_only; 369ba6664aeSJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 370ba6664aeSJames Wright const CeedScalar dx = stg_ctx->dx; 371ba6664aeSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 372ba6664aeSJames Wright const CeedScalar time = stg_ctx->time; 373ba6664aeSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 374ba6664aeSJames Wright const CeedScalar P0 = stg_ctx->P0; 375ba6664aeSJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 376ba6664aeSJames Wright const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 377ba6664aeSJames Wright const CeedScalar Rd = cp - cv; 378ba6664aeSJames Wright const CeedScalar gamma = cp / cv; 379ba6664aeSJames Wright 380*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 381ba6664aeSJames Wright const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0); 382ba6664aeSJames Wright const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 383ba6664aeSJames Wright const CeedScalar dXdx[2][3] = { 384ba6664aeSJames Wright {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 385ba6664aeSJames Wright {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 386ba6664aeSJames Wright }; 387ba6664aeSJames Wright 388ba6664aeSJames Wright CeedScalar h[3]; 389ba6664aeSJames Wright h[0] = dx; 390*2b730f8bSJeremy L Thompson for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 391ba6664aeSJames Wright 392ba6664aeSJames Wright InterpolateProfile(X[1][i], ubar, cij, &eps, <, stg_ctx); 393ba6664aeSJames Wright if (!mean_only) { 394ba6664aeSJames Wright CalcSpectrum(X[1][i], eps, lt, h, mu / rho, qn, stg_ctx); 395ba6664aeSJames Wright STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 396ba6664aeSJames Wright } else { 397ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 398ba6664aeSJames Wright } 399ba6664aeSJames Wright 4004dbab5e5SJames Wright const CeedScalar E_kinetic = .5 * rho * Dot3(u, u); 401ba6664aeSJames Wright CeedScalar E_internal, P; 402ba6664aeSJames Wright if (prescribe_T) { 403ba6664aeSJames Wright // Temperature is being set weakly (theta0) and for constant cv this sets E_internal 404ba6664aeSJames Wright E_internal = rho * cv * theta0; 405ba6664aeSJames Wright // Find pressure using 406ba6664aeSJames Wright P = rho * Rd * theta0; // interior rho with exterior T 407ba6664aeSJames Wright } else { 408ba6664aeSJames Wright E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution 409ba6664aeSJames Wright P = E_internal * (gamma - 1.); 410ba6664aeSJames Wright } 411ba6664aeSJames Wright 412ba6664aeSJames Wright const CeedScalar wdetJb = (is_implicit ? -1. : 1.) * q_data_sur[0][i]; 413ba6664aeSJames Wright // ---- Normal vect 414*2b730f8bSJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 415ba6664aeSJames Wright 416ba6664aeSJames Wright const CeedScalar E = E_internal + E_kinetic; 417ba6664aeSJames Wright 418ba6664aeSJames Wright // Velocity normal to the boundary 4194dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, u); 4204dbab5e5SJames Wright 421ba6664aeSJames Wright // The Physics 422ba6664aeSJames Wright // Zero v so all future terms can safely sum into it 423ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 424ba6664aeSJames Wright 425ba6664aeSJames Wright // The Physics 426ba6664aeSJames Wright // -- Density 427ba6664aeSJames Wright v[0][i] -= wdetJb * rho * u_normal; 428ba6664aeSJames Wright 429ba6664aeSJames Wright // -- Momentum 430*2b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 431ba6664aeSJames Wright 432ba6664aeSJames Wright // -- Total Energy Density 433ba6664aeSJames Wright v[4][i] -= wdetJb * u_normal * (E + P); 4344dbab5e5SJames Wright 4354dbab5e5SJames Wright jac_data_sur[0][i] = rho; 4364dbab5e5SJames Wright jac_data_sur[1][i] = u[0]; 4374dbab5e5SJames Wright jac_data_sur[2][i] = u[1]; 4384dbab5e5SJames Wright jac_data_sur[3][i] = u[2]; 4394dbab5e5SJames Wright jac_data_sur[4][i] = E; 4404dbab5e5SJames Wright for (int j = 0; j < 6; j++) jac_data_sur[5 + j][i] = 0.; 441ba6664aeSJames Wright } 442ba6664aeSJames Wright return 0; 443ba6664aeSJames Wright } 444ba6664aeSJames Wright 445*2b730f8bSJeremy L Thompson CEED_QFUNCTION(STGShur14_Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 4464dbab5e5SJames Wright // *INDENT-OFF* 4474dbab5e5SJames Wright // Inputs 448*2b730f8bSJeremy L Thompson const CeedScalar(*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 4494dbab5e5SJames Wright (*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 4504dbab5e5SJames Wright // Outputs 4514dbab5e5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 4524dbab5e5SJames Wright // *INDENT-ON* 4534dbab5e5SJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 4544dbab5e5SJames Wright const bool implicit = stg_ctx->is_implicit; 4554dbab5e5SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 4564dbab5e5SJames Wright const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 4574dbab5e5SJames Wright const CeedScalar Rd = cp - cv; 4584dbab5e5SJames Wright const CeedScalar gamma = cp / cv; 4594dbab5e5SJames Wright 4604dbab5e5SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 4614dbab5e5SJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 4624dbab5e5SJames Wright 4634dbab5e5SJames Wright CeedPragmaSIMD 4644dbab5e5SJames Wright // Quadrature Point Loop 4654dbab5e5SJames Wright for (CeedInt i = 0; i < Q; i++) { 4664dbab5e5SJames Wright // Setup 4674dbab5e5SJames Wright // -- Interp-to-Interp q_data 4684dbab5e5SJames Wright // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 4694dbab5e5SJames Wright // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 4704dbab5e5SJames Wright // We can effect this by swapping the sign on this weight 4714dbab5e5SJames Wright const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 4724dbab5e5SJames Wright 4734dbab5e5SJames Wright // Calculate inflow values 4744dbab5e5SJames Wright CeedScalar velocity[3]; 4754dbab5e5SJames Wright for (CeedInt j = 0; j < 3; j++) velocity[j] = jac_data_sur[5 + j][i]; 4764dbab5e5SJames Wright 4774dbab5e5SJames Wright // enabling user to choose between weak T and weak rho inflow 4784dbab5e5SJames Wright CeedScalar drho, dE, dP; 4794dbab5e5SJames Wright if (prescribe_T) { 4804dbab5e5SJames Wright // rho should be from the current solution 4814dbab5e5SJames Wright drho = dq[0][i]; 4824dbab5e5SJames Wright CeedScalar dE_internal = drho * cv * theta0; 4834dbab5e5SJames Wright CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity); 4844dbab5e5SJames Wright dE = dE_internal + dE_kinetic; 4854dbab5e5SJames Wright dP = drho * Rd * theta0; // interior rho with exterior T 4864dbab5e5SJames Wright } else { // rho specified, E_internal from solution 4874dbab5e5SJames Wright drho = 0; 4884dbab5e5SJames Wright dE = dq[4][i]; 4894dbab5e5SJames Wright dP = dE * (gamma - 1.); 4904dbab5e5SJames Wright } 491*2b730f8bSJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 4924dbab5e5SJames Wright 4934dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, velocity); 4944dbab5e5SJames Wright 4954dbab5e5SJames Wright v[0][i] = -wdetJb * drho * u_normal; 496*2b730f8bSJeremy L Thompson for (int j = 0; j < 3; j++) v[j + 1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP); 4974dbab5e5SJames Wright v[4][i] = -wdetJb * u_normal * (dE + dP); 4984dbab5e5SJames Wright } // End Quadrature Point Loop 4994dbab5e5SJames Wright return 0; 5004dbab5e5SJames Wright } 5014dbab5e5SJames Wright 5020a6353c2SJames Wright /******************************************************************** 5030a6353c2SJames Wright * @brief QFunction to calculate the strongly enforce inflow BC 5040a6353c2SJames Wright * 5050a6353c2SJames Wright * This QF is for the strong application of STG via libCEED (rather than 5060a6353c2SJames Wright * through the native PETSc `DMAddBoundary` -> `bcFunc` method. 5070a6353c2SJames Wright */ 508*2b730f8bSJeremy L Thompson CEED_QFUNCTION(STGShur14_Inflow_StrongQF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 5090a6353c2SJames Wright //*INDENT-OFF* 510*2b730f8bSJeremy L Thompson const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*coords)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1], 511*2b730f8bSJeremy L Thompson (*scale) = (const CeedScalar(*))in[2], (*stg_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 5120a6353c2SJames Wright 5130a6353c2SJames Wright CeedScalar(*bcval)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 5140a6353c2SJames Wright //*INDENT-ON* 5150a6353c2SJames Wright 5160a6353c2SJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 51762e628f8SJames Wright CeedScalar u[3], ubar[3], cij[6], eps, lt; 5180a6353c2SJames Wright const bool mean_only = stg_ctx->mean_only; 5190a6353c2SJames Wright const CeedScalar dx = stg_ctx->dx; 5200a6353c2SJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 5210a6353c2SJames Wright const CeedScalar time = stg_ctx->time; 5220a6353c2SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 5230a6353c2SJames Wright const CeedScalar P0 = stg_ctx->P0; 5240a6353c2SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 5250a6353c2SJames Wright const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 5260a6353c2SJames Wright const CeedScalar Rd = cp - cv; 5270a6353c2SJames Wright const CeedScalar rho = P0 / (Rd * theta0); 5280a6353c2SJames Wright 529*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 5300a6353c2SJames Wright const CeedScalar x[] = {coords[0][i], coords[1][i], coords[2][i]}; 5310a6353c2SJames Wright const CeedScalar dXdx[2][3] = { 5320a6353c2SJames Wright {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 5330a6353c2SJames Wright {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 5340a6353c2SJames Wright }; 5350a6353c2SJames Wright 5360a6353c2SJames Wright CeedScalar h[3]; 5370a6353c2SJames Wright h[0] = dx; 538*2b730f8bSJeremy L Thompson for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 5390a6353c2SJames Wright 5400a6353c2SJames Wright InterpolateProfile(coords[1][i], ubar, cij, &eps, <, stg_ctx); 5410a6353c2SJames Wright if (!mean_only) { 54262e628f8SJames Wright if (1) { 543*2b730f8bSJeremy L Thompson STGShur14_Calc_PrecompEktot(x, time, ubar, cij, stg_data[0][i], h, x[1], eps, lt, mu / rho, u, stg_ctx); 54462e628f8SJames Wright } else { // Original way 54562e628f8SJames Wright CeedScalar qn[STG_NMODES_MAX]; 54662e628f8SJames Wright CalcSpectrum(coords[1][i], eps, lt, h, mu / rho, qn, stg_ctx); 54762e628f8SJames Wright STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 54862e628f8SJames Wright } 5490a6353c2SJames Wright } else { 5500a6353c2SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 5510a6353c2SJames Wright } 5520a6353c2SJames Wright 55397baf651SJames Wright switch (stg_ctx->newtonian_ctx.state_var) { 55497baf651SJames Wright case STATEVAR_CONSERVATIVE: 5550a6353c2SJames Wright bcval[0][i] = scale[i] * rho; 5560a6353c2SJames Wright bcval[1][i] = scale[i] * rho * u[0]; 5570a6353c2SJames Wright bcval[2][i] = scale[i] * rho * u[1]; 5580a6353c2SJames Wright bcval[3][i] = scale[i] * rho * u[2]; 559cf3d54ffSJames Wright bcval[4][i] = 0.; 56097baf651SJames Wright break; 56197baf651SJames Wright 56297baf651SJames Wright case STATEVAR_PRIMITIVE: 56397baf651SJames Wright bcval[0][i] = 0; 56497baf651SJames Wright bcval[1][i] = scale[i] * u[0]; 56597baf651SJames Wright bcval[2][i] = scale[i] * u[1]; 56697baf651SJames Wright bcval[3][i] = scale[i] * u[2]; 56797baf651SJames Wright bcval[4][i] = scale[i] * theta0; 56897baf651SJames Wright break; 5690a6353c2SJames Wright } 5707c4551aaSJames Wright } 5710a6353c2SJames Wright return 0; 5720a6353c2SJames Wright } 5730a6353c2SJames Wright 574ba6664aeSJames Wright #endif // stg_shur14_h 575