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. 15*c0b5abf0SJeremy L Thompson #include <ceed/types.h> 16*c0b5abf0SJeremy L Thompson #ifndef CEED_RUNNING_JIT_PASS 17c9c2c079SJeremy L Thompson #include <math.h> 18ba6664aeSJames Wright #include <stdlib.h> 19*c0b5abf0SJeremy L Thompson #endif 202b730f8bSJeremy L Thompson 2146603fc5SJames Wright #include "newtonian_state.h" 228756a6ccSJames Wright #include "setupgeo_helpers.h" 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 */ 412b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar wall_dist, CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt, 42cbef7084SJames 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 59175f00a6SJames Wright CeedScalar coeff = (wall_dist - prof_wd[idx - 1]) / (prof_wd[idx] - prof_wd[idx - 1]); 60175f00a6SJames Wright 61ba6664aeSJames Wright ubar[0] = prof_ubar[0 * nprofs + idx - 1] + coeff * (prof_ubar[0 * nprofs + idx] - prof_ubar[0 * nprofs + idx - 1]); 62ba6664aeSJames Wright ubar[1] = prof_ubar[1 * nprofs + idx - 1] + coeff * (prof_ubar[1 * nprofs + idx] - prof_ubar[1 * nprofs + idx - 1]); 63ba6664aeSJames Wright ubar[2] = prof_ubar[2 * nprofs + idx - 1] + coeff * (prof_ubar[2 * nprofs + idx] - prof_ubar[2 * nprofs + idx - 1]); 64ba6664aeSJames Wright cij[0] = prof_cij[0 * nprofs + idx - 1] + coeff * (prof_cij[0 * nprofs + idx] - prof_cij[0 * nprofs + idx - 1]); 65ba6664aeSJames Wright cij[1] = prof_cij[1 * nprofs + idx - 1] + coeff * (prof_cij[1 * nprofs + idx] - prof_cij[1 * nprofs + idx - 1]); 66ba6664aeSJames Wright cij[2] = prof_cij[2 * nprofs + idx - 1] + coeff * (prof_cij[2 * nprofs + idx] - prof_cij[2 * nprofs + idx - 1]); 67ba6664aeSJames Wright cij[3] = prof_cij[3 * nprofs + idx - 1] + coeff * (prof_cij[3 * nprofs + idx] - prof_cij[3 * nprofs + idx - 1]); 68ba6664aeSJames Wright cij[4] = prof_cij[4 * nprofs + idx - 1] + coeff * (prof_cij[4 * nprofs + idx] - prof_cij[4 * nprofs + idx - 1]); 69ba6664aeSJames Wright cij[5] = prof_cij[5 * nprofs + idx - 1] + coeff * (prof_cij[5 * nprofs + idx] - prof_cij[5 * nprofs + idx - 1]); 70ba6664aeSJames Wright *eps = prof_eps[idx - 1] + coeff * (prof_eps[idx] - prof_eps[idx - 1]); 71ba6664aeSJames Wright *lt = prof_lt[idx - 1] + coeff * (prof_lt[idx] - prof_lt[idx - 1]); 72175f00a6SJames Wright } else { // y outside bounds of prof_wd 73ba6664aeSJames Wright ubar[0] = prof_ubar[1 * nprofs - 1]; 74ba6664aeSJames Wright ubar[1] = prof_ubar[2 * nprofs - 1]; 75ba6664aeSJames Wright ubar[2] = prof_ubar[3 * nprofs - 1]; 76ba6664aeSJames Wright cij[0] = prof_cij[1 * nprofs - 1]; 77ba6664aeSJames Wright cij[1] = prof_cij[2 * nprofs - 1]; 78ba6664aeSJames Wright cij[2] = prof_cij[3 * nprofs - 1]; 79ba6664aeSJames Wright cij[3] = prof_cij[4 * nprofs - 1]; 80ba6664aeSJames Wright cij[4] = prof_cij[5 * nprofs - 1]; 81ba6664aeSJames Wright cij[5] = prof_cij[6 * nprofs - 1]; 82ba6664aeSJames Wright *eps = prof_eps[nprofs - 1]; 83ba6664aeSJames Wright *lt = prof_lt[nprofs - 1]; 84ba6664aeSJames Wright } 85ba6664aeSJames Wright } 86ba6664aeSJames Wright 87ba6664aeSJames Wright /* 88e159aeacSJames Wright * @brief Calculate spectrum coefficient, qn 89e159aeacSJames Wright * 90e159aeacSJames Wright * Calculates q_n at a given distance to the wall 91e159aeacSJames Wright * 92e159aeacSJames Wright * @param[in] kappa nth wavenumber 93e159aeacSJames Wright * @param[in] dkappa Difference between wavenumbers 94e159aeacSJames Wright * @param[in] keta Dissipation wavenumber 95e159aeacSJames Wright * @param[in] kcut Mesh-induced cutoff wavenumber 96e159aeacSJames Wright * @param[in] ke Energy-containing wavenumber 97f8839eb4SJames Wright * @param[in] Ektot_inv Inverse of total turbulent kinetic energy of spectrum 98e159aeacSJames Wright * @returns qn Spectrum coefficient 99e159aeacSJames Wright */ 1002b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER CeedScalar Calc_qn(const CeedScalar kappa, const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut, 10162e628f8SJames Wright const CeedScalar ke, const CeedScalar Ektot_inv) { 1022b730f8bSJeremy L Thompson const CeedScalar feta_x_fcut = exp(-Square(12 * kappa / keta) - Cube(4 * Max(kappa - 0.9 * kcut, 0) / kcut)); 1032b730f8bSJeremy L Thompson return pow(kappa / ke, 4.) * pow(1 + 2.4 * Square(kappa / ke), -17. / 6) * feta_x_fcut * dkappa * Ektot_inv; 104e159aeacSJames Wright } 105e159aeacSJames Wright 106e159aeacSJames Wright // Calculate hmax, ke, keta, and kcut 107831dbe9eSJames Wright CEED_QFUNCTION_HELPER void SpectrumConstants(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar hNodSep[3], 1082b730f8bSJeremy L Thompson const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke, CeedScalar *keta, CeedScalar *kcut) { 109831dbe9eSJames Wright *hmax = Max(Max(hNodSep[0], hNodSep[1]), hNodSep[2]); 110175f00a6SJames Wright *ke = wall_dist == 0 ? 1e16 : 2 * M_PI / Min(2 * wall_dist, 3 * lt); 111e159aeacSJames Wright *keta = 2 * M_PI * pow(Cube(nu) / eps, -0.25); 112831dbe9eSJames Wright *kcut = M_PI / Min(Max(Max(hNodSep[1], hNodSep[2]), 0.3 * (*hmax)) + 0.1 * wall_dist, *hmax); 113e159aeacSJames Wright } 114e159aeacSJames Wright 115e159aeacSJames Wright /* 116ba6664aeSJames Wright * @brief Calculate spectrum coefficients for STG 117ba6664aeSJames Wright * 118ba6664aeSJames Wright * Calculates q_n at a given distance to the wall 119ba6664aeSJames Wright * 120175f00a6SJames Wright * @param[in] wall_dist Distance to the nearest wall 121175f00a6SJames Wright * @param[in] eps Turbulent dissipation w/rt wall_dist 122175f00a6SJames Wright * @param[in] lt Turbulent length scale w/rt wall_dist 123831dbe9eSJames Wright * @param[in] h_node_sep Element lengths in coordinate directions 124ba6664aeSJames Wright * @param[in] nu Dynamic Viscosity; 125ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 126ba6664aeSJames Wright * @param[out] qn Spectrum coefficients, [nmodes] 127ba6664aeSJames Wright */ 128831dbe9eSJames Wright CEED_QFUNCTION_HELPER void CalcSpectrum(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h_node_sep[3], 129cbef7084SJames Wright const CeedScalar nu, CeedScalar qn[], const StgShur14Context stg_ctx) { 130ba6664aeSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 131ba6664aeSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 132e159aeacSJames Wright CeedScalar hmax, ke, keta, kcut, Ektot = 0.0; 1332b730f8bSJeremy L Thompson 134831dbe9eSJames Wright SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut); 135ba6664aeSJames Wright 136ba6664aeSJames Wright for (CeedInt n = 0; n < nmodes; n++) { 137e159aeacSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 138e159aeacSJames Wright qn[n] = Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 139ba6664aeSJames Wright Ektot += qn[n]; 140ba6664aeSJames Wright } 141ba6664aeSJames Wright 142961c9c98SJames Wright if (Ektot == 0) return; 143ba6664aeSJames Wright for (CeedInt n = 0; n < nmodes; n++) qn[n] /= Ektot; 144ba6664aeSJames Wright } 145ba6664aeSJames Wright 146ba6664aeSJames Wright /****************************************************** 147ba6664aeSJames Wright * @brief Calculate u(x,t) for STG inflow condition 148ba6664aeSJames Wright * 149ba6664aeSJames Wright * @param[in] X Location to evaluate u(X,t) 150ba6664aeSJames Wright * @param[in] t Time to evaluate u(X,t) 151ba6664aeSJames Wright * @param[in] ubar Mean velocity at X 152ba6664aeSJames Wright * @param[in] cij Cholesky decomposition at X 153ba6664aeSJames Wright * @param[in] qn Wavemode amplitudes at X, [nmodes] 154ba6664aeSJames Wright * @param[out] u Velocity at X and t 155ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 156ba6664aeSJames Wright */ 157cbef7084SJames Wright CEED_QFUNCTION_HELPER void StgShur14Calc(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 158cbef7084SJames Wright const CeedScalar qn[], CeedScalar u[3], const StgShur14Context stg_ctx) { 159ba6664aeSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 160ba6664aeSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 161ba6664aeSJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 162ba6664aeSJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 163ba6664aeSJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 164ba6664aeSJames Wright CeedScalar xdotd, vp[3] = {0.}; 165ba6664aeSJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 166ba6664aeSJames Wright 1672b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 168ba6664aeSJames Wright xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 169ba6664aeSJames Wright xdotd = 0.; 170ba6664aeSJames Wright for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 171ba6664aeSJames Wright const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 172961c9c98SJames Wright vp[0] += sqrt(qn[n]) * sigma[0 * nmodes + n] * cos_kxdp; 173961c9c98SJames Wright vp[1] += sqrt(qn[n]) * sigma[1 * nmodes + n] * cos_kxdp; 174961c9c98SJames Wright vp[2] += sqrt(qn[n]) * sigma[2 * nmodes + n] * cos_kxdp; 175ba6664aeSJames Wright } 176961c9c98SJames Wright for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 177ba6664aeSJames Wright 178ba6664aeSJames Wright u[0] = ubar[0] + cij[0] * vp[0]; 179ba6664aeSJames Wright u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 180ba6664aeSJames Wright u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 181ba6664aeSJames Wright } 182ba6664aeSJames Wright 183b277271eSJames Wright /****************************************************** 184b277271eSJames Wright * @brief Calculate u(x,t) for STG inflow condition 185b277271eSJames Wright * 186b277271eSJames Wright * @param[in] X Location to evaluate u(X,t) 187b277271eSJames Wright * @param[in] t Time to evaluate u(X,t) 188b277271eSJames Wright * @param[in] ubar Mean velocity at X 189b277271eSJames Wright * @param[in] cij Cholesky decomposition at X 190175f00a6SJames Wright * @param[in] Ektot Total spectrum energy at this location 191831dbe9eSJames Wright * @param[in] h_node_sep Element size in 3 directions 192175f00a6SJames Wright * @param[in] wall_dist Distance to closest wall 193175f00a6SJames Wright * @param[in] eps Turbulent dissipation 194175f00a6SJames Wright * @param[in] lt Turbulent length scale 195b277271eSJames Wright * @param[out] u Velocity at X and t 196b277271eSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 197b277271eSJames Wright */ 198cbef7084SJames Wright CEED_QFUNCTION_HELPER void StgShur14Calc_PrecompEktot(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 199831dbe9eSJames Wright const CeedScalar Ektot, const CeedScalar h_node_sep[3], const CeedScalar wall_dist, 200831dbe9eSJames Wright const CeedScalar eps, const CeedScalar lt, const CeedScalar nu, CeedScalar u[3], 201831dbe9eSJames Wright const StgShur14Context stg_ctx) { 202b277271eSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 203b277271eSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 204b277271eSJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 205b277271eSJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 206b277271eSJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 207b277271eSJames Wright CeedScalar hmax, ke, keta, kcut; 208831dbe9eSJames Wright SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut); 209b277271eSJames Wright CeedScalar xdotd, vp[3] = {0.}; 210b277271eSJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 211b277271eSJames Wright 2122b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 213b277271eSJames Wright xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 214b277271eSJames Wright xdotd = 0.; 215b277271eSJames Wright for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 216b277271eSJames Wright const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 217b277271eSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 218b277271eSJames Wright const CeedScalar qn = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot); 219b277271eSJames Wright vp[0] += sqrt(qn) * sigma[0 * nmodes + n] * cos_kxdp; 220b277271eSJames Wright vp[1] += sqrt(qn) * sigma[1 * nmodes + n] * cos_kxdp; 221b277271eSJames Wright vp[2] += sqrt(qn) * sigma[2 * nmodes + n] * cos_kxdp; 222b277271eSJames Wright } 223b277271eSJames Wright for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 224b277271eSJames Wright 225b277271eSJames Wright u[0] = ubar[0] + cij[0] * vp[0]; 226b277271eSJames Wright u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 227b277271eSJames Wright u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 228b277271eSJames Wright } 229b277271eSJames Wright 23062e628f8SJames Wright // Create preprocessed input for the stg calculation 23162e628f8SJames Wright // 23262e628f8SJames Wright // stg_data[0] = 1 / Ektot (inverse of total spectrum energy) 233cbef7084SJames Wright CEED_QFUNCTION(StgShur14Preprocess)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 23466170c20SJames Wright const CeedScalar(*dXdx_q)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0]; 23546603fc5SJames Wright const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 236b277271eSJames Wright 237b277271eSJames Wright CeedScalar(*stg_data) = (CeedScalar(*))out[0]; 238b277271eSJames Wright 239b277271eSJames Wright CeedScalar ubar[3], cij[6], eps, lt; 240cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 241b277271eSJames Wright const CeedScalar dx = stg_ctx->dx; 242b277271eSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 243b277271eSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 244b277271eSJames Wright const CeedScalar P0 = stg_ctx->P0; 24546603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 246b277271eSJames Wright const CeedScalar rho = P0 / (Rd * theta0); 247b277271eSJames Wright const CeedScalar nu = mu / rho; 248b277271eSJames Wright 249b277271eSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 250b277271eSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 2515dc40723SJames Wright CeedScalar hmax, ke, keta, kcut; 252b277271eSJames Wright 2532b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 254175f00a6SJames Wright const CeedScalar wall_dist = x[1][i]; 255b277271eSJames Wright const CeedScalar dXdx[2][3] = { 25666170c20SJames Wright {dXdx_q[0][0][i], dXdx_q[0][1][i], dXdx_q[0][2][i]}, 25766170c20SJames Wright {dXdx_q[1][0][i], dXdx_q[1][1][i], dXdx_q[1][2][i]}, 258b277271eSJames Wright }; 259b277271eSJames Wright 260831dbe9eSJames Wright CeedScalar h_node_sep[3]; 261831dbe9eSJames Wright h_node_sep[0] = dx; 262831dbe9eSJames 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]); 263831dbe9eSJames Wright ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3); 264b277271eSJames Wright 265175f00a6SJames Wright InterpolateProfile(wall_dist, ubar, cij, &eps, <, stg_ctx); 266831dbe9eSJames Wright SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut); 267b277271eSJames Wright 268b277271eSJames Wright // Calculate total TKE per spectrum 269d97dc904SJames Wright CeedScalar Ek_tot = 0; 2702b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 271b277271eSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 272d97dc904SJames Wright Ek_tot += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 273b277271eSJames Wright } 274d97dc904SJames Wright // avoid underflowed and poorly defined spectrum coefficients 275d97dc904SJames Wright stg_data[i] = Ek_tot != 0 ? 1 / Ek_tot : 0; 276b277271eSJames Wright } 277b277271eSJames Wright return 0; 278b277271eSJames Wright } 279b277271eSJames Wright 280b77c53c9SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition 281cbef7084SJames Wright CEED_QFUNCTION(ICsStg)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 28246603fc5SJames Wright const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 2831c299e57SJeremy L Thompson const CeedScalar(*J)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[1]; 284b77c53c9SJames Wright CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 285b77c53c9SJames Wright 286cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 287a2d72b6fSJames Wright const NewtonianIdealGasContext gas = &stg_ctx->newtonian_ctx; 28889060322SJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 28989060322SJames Wright const CeedScalar dx = stg_ctx->dx; 29089060322SJames Wright const CeedScalar time = stg_ctx->time; 291b77c53c9SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 292b77c53c9SJames Wright const CeedScalar P0 = stg_ctx->P0; 293a2d72b6fSJames Wright const CeedScalar rho = P0 / (GasConstant(gas) * theta0); 294a2d72b6fSJames Wright const CeedScalar nu = gas->mu / rho; 295b77c53c9SJames Wright 2962b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 29789060322SJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 2981c299e57SJeremy L Thompson CeedScalar dXdx[3][3]; 2998756a6ccSJames Wright InvertMappingJacobian_3D(Q, i, J, dXdx, NULL); 300831dbe9eSJames Wright CeedScalar h_node_sep[3]; 301831dbe9eSJames Wright h_node_sep[0] = dx; 302831dbe9eSJames 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])); 303831dbe9eSJames Wright ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3); 30489060322SJames Wright 30589060322SJames Wright InterpolateProfile(x_i[1], ubar, cij, &eps, <, stg_ctx); 30689060322SJames Wright if (stg_ctx->use_fluctuating_IC) { 307831dbe9eSJames Wright CalcSpectrum(x_i[1], eps, lt, h_node_sep, nu, qn, stg_ctx); 308cbef7084SJames Wright StgShur14Calc(x_i, time, ubar, cij, qn, u, stg_ctx); 30989060322SJames Wright } else { 31089060322SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 31189060322SJames Wright } 312b77c53c9SJames Wright 313a2d72b6fSJames Wright CeedScalar Y[5] = {P0, u[0], u[1], u[2], theta0}, q[5]; 314a2d72b6fSJames Wright State s = StateFromY(gas, Y); 315a2d72b6fSJames Wright StateToQ(gas, s, q, gas->state_var); 316a2d72b6fSJames Wright for (CeedInt j = 0; j < 5; j++) { 317a2d72b6fSJames Wright q0[j][i] = q[j]; 3187c4551aaSJames Wright } 319f0b01153SJames Wright } 320b77c53c9SJames Wright return 0; 321b77c53c9SJames Wright } 322b77c53c9SJames Wright 323ba6664aeSJames Wright /******************************************************************** 324ba6664aeSJames Wright * @brief QFunction to calculate the inflow boundary condition 325ba6664aeSJames Wright * 326ba6664aeSJames Wright * This will loop through quadrature points, calculate the wavemode amplitudes 327ba6664aeSJames Wright * at each location, then calculate the actual velocity. 328ba6664aeSJames Wright */ 329cbef7084SJames Wright CEED_QFUNCTION(StgShur14Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 33046603fc5SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 331f3e15844SJames Wright const CeedScalar(*q_data_sur) = in[2]; 33246603fc5SJames Wright const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 333ba6664aeSJames Wright 33446603fc5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 335f3e15844SJames Wright CeedScalar(*jac_data_sur) = out[1]; 336ba6664aeSJames Wright 337cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 338ba6664aeSJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 339ba6664aeSJames Wright const bool is_implicit = stg_ctx->is_implicit; 340ba6664aeSJames Wright const bool mean_only = stg_ctx->mean_only; 341ba6664aeSJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 342ba6664aeSJames Wright const CeedScalar dx = stg_ctx->dx; 343ba6664aeSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 344ba6664aeSJames Wright const CeedScalar time = stg_ctx->time; 345ba6664aeSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 346ba6664aeSJames Wright const CeedScalar P0 = stg_ctx->P0; 347ba6664aeSJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 34846603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 34946603fc5SJames Wright const CeedScalar gamma = HeatCapacityRatio(&stg_ctx->newtonian_ctx); 350ba6664aeSJames Wright 3512b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 352ba6664aeSJames Wright const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0); 353ba6664aeSJames Wright const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 354f3e15844SJames Wright CeedScalar wdetJb, dXdx[2][3], norm[3]; 355f3e15844SJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, dXdx, norm); 356f3e15844SJames Wright wdetJb *= is_implicit ? -1. : 1.; 357ba6664aeSJames Wright 358831dbe9eSJames Wright CeedScalar h_node_sep[3]; 359831dbe9eSJames Wright h_node_sep[0] = dx; 360831dbe9eSJames Wright for (CeedInt j = 1; j < 3; j++) h_node_sep[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 361831dbe9eSJames Wright ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3); 362ba6664aeSJames Wright 363ba6664aeSJames Wright InterpolateProfile(X[1][i], ubar, cij, &eps, <, stg_ctx); 364ba6664aeSJames Wright if (!mean_only) { 365831dbe9eSJames Wright CalcSpectrum(X[1][i], eps, lt, h_node_sep, mu / rho, qn, stg_ctx); 366cbef7084SJames Wright StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx); 367ba6664aeSJames Wright } else { 368ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 369ba6664aeSJames Wright } 370ba6664aeSJames Wright 3714dbab5e5SJames Wright const CeedScalar E_kinetic = .5 * rho * Dot3(u, u); 372ba6664aeSJames Wright CeedScalar E_internal, P; 373ba6664aeSJames Wright if (prescribe_T) { 374ba6664aeSJames Wright // Temperature is being set weakly (theta0) and for constant cv this sets E_internal 375ba6664aeSJames Wright E_internal = rho * cv * theta0; 376ba6664aeSJames Wright // Find pressure using 377ba6664aeSJames Wright P = rho * Rd * theta0; // interior rho with exterior T 378ba6664aeSJames Wright } else { 379ba6664aeSJames Wright E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution 380ba6664aeSJames Wright P = E_internal * (gamma - 1.); 381ba6664aeSJames Wright } 382ba6664aeSJames Wright 383ba6664aeSJames Wright const CeedScalar E = E_internal + E_kinetic; 384ba6664aeSJames Wright 385ba6664aeSJames Wright // Velocity normal to the boundary 3864dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, u); 3874dbab5e5SJames Wright 388ba6664aeSJames Wright // The Physics 389ba6664aeSJames Wright // Zero v so all future terms can safely sum into it 390ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 391ba6664aeSJames Wright 392ba6664aeSJames Wright // The Physics 393ba6664aeSJames Wright // -- Density 394ba6664aeSJames Wright v[0][i] -= wdetJb * rho * u_normal; 395ba6664aeSJames Wright 396ba6664aeSJames Wright // -- Momentum 3972b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 398ba6664aeSJames Wright 399ba6664aeSJames Wright // -- Total Energy Density 400ba6664aeSJames Wright v[4][i] -= wdetJb * u_normal * (E + P); 4014dbab5e5SJames Wright 402f3e15844SJames Wright const CeedScalar U[] = {rho, u[0], u[1], u[2], E}, kmstress[6] = {0.}; 403f3e15844SJames Wright StoredValuesPack(Q, i, 0, 5, U, jac_data_sur); 404f3e15844SJames Wright StoredValuesPack(Q, i, 5, 6, kmstress, jac_data_sur); 405ba6664aeSJames Wright } 406ba6664aeSJames Wright return 0; 407ba6664aeSJames Wright } 408ba6664aeSJames Wright 409cbef7084SJames Wright CEED_QFUNCTION(StgShur14Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 41046603fc5SJames Wright const CeedScalar(*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 41146603fc5SJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 41246603fc5SJames Wright const CeedScalar(*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 4134dbab5e5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 41446603fc5SJames Wright 415cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 4164dbab5e5SJames Wright const bool implicit = stg_ctx->is_implicit; 4174dbab5e5SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 41846603fc5SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 41946603fc5SJames Wright const CeedScalar gamma = HeatCapacityRatio(&stg_ctx->newtonian_ctx); 4204dbab5e5SJames Wright 4214dbab5e5SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 4224dbab5e5SJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 4234dbab5e5SJames Wright 424f0b01153SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 4254dbab5e5SJames Wright // Setup 4264dbab5e5SJames Wright // -- Interp-to-Interp q_data 4274dbab5e5SJames Wright // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 4284dbab5e5SJames Wright // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 4294dbab5e5SJames Wright // We can effect this by swapping the sign on this weight 4304dbab5e5SJames Wright const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 4314dbab5e5SJames Wright 4324dbab5e5SJames Wright // Calculate inflow values 4334dbab5e5SJames Wright CeedScalar velocity[3]; 4344dbab5e5SJames Wright for (CeedInt j = 0; j < 3; j++) velocity[j] = jac_data_sur[5 + j][i]; 435f3e15844SJames Wright // TODO This is almost certainly a bug. Velocity isn't stored here, only 0s. 4364dbab5e5SJames Wright 4374dbab5e5SJames Wright // enabling user to choose between weak T and weak rho inflow 4384dbab5e5SJames Wright CeedScalar drho, dE, dP; 4394dbab5e5SJames Wright if (prescribe_T) { 4404dbab5e5SJames Wright // rho should be from the current solution 4414dbab5e5SJames Wright drho = dq[0][i]; 4424dbab5e5SJames Wright CeedScalar dE_internal = drho * cv * theta0; 4434dbab5e5SJames Wright CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity); 4444dbab5e5SJames Wright dE = dE_internal + dE_kinetic; 4454dbab5e5SJames Wright dP = drho * Rd * theta0; // interior rho with exterior T 4464dbab5e5SJames Wright } else { // rho specified, E_internal from solution 4474dbab5e5SJames Wright drho = 0; 4484dbab5e5SJames Wright dE = dq[4][i]; 4494dbab5e5SJames Wright dP = dE * (gamma - 1.); 4504dbab5e5SJames Wright } 4512b730f8bSJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 4524dbab5e5SJames Wright 4534dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, velocity); 4544dbab5e5SJames Wright 4554dbab5e5SJames Wright v[0][i] = -wdetJb * drho * u_normal; 4562b730f8bSJeremy L Thompson for (int j = 0; j < 3; j++) v[j + 1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP); 4574dbab5e5SJames Wright v[4][i] = -wdetJb * u_normal * (dE + dP); 458f0b01153SJames Wright } 4594dbab5e5SJames Wright return 0; 4604dbab5e5SJames Wright } 4614dbab5e5SJames Wright 4620a6353c2SJames Wright /******************************************************************** 4630a6353c2SJames Wright * @brief QFunction to calculate the strongly enforce inflow BC 4640a6353c2SJames Wright * 4650a6353c2SJames Wright * This QF is for the strong application of STG via libCEED (rather than 4660a6353c2SJames Wright * through the native PETSc `DMAddBoundary` -> `bcFunc` method. 4670a6353c2SJames Wright */ 468cbef7084SJames Wright CEED_QFUNCTION(StgShur14InflowStrongQF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 46966170c20SJames Wright const CeedScalar(*dXdx_q)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0]; 47046603fc5SJames Wright const CeedScalar(*coords)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 47146603fc5SJames Wright const CeedScalar(*scale) = (const CeedScalar(*))in[2]; 472f8839eb4SJames Wright const CeedScalar(*inv_Ektotal) = (const CeedScalar(*))in[3]; 4730a6353c2SJames Wright CeedScalar(*bcval)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 4740a6353c2SJames Wright 475cbef7084SJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 476a2d72b6fSJames Wright const NewtonianIdealGasContext gas = &stg_ctx->newtonian_ctx; 47762e628f8SJames Wright CeedScalar u[3], ubar[3], cij[6], eps, lt; 4780a6353c2SJames Wright const bool mean_only = stg_ctx->mean_only; 4790a6353c2SJames Wright const CeedScalar dx = stg_ctx->dx; 4800a6353c2SJames Wright const CeedScalar time = stg_ctx->time; 4810a6353c2SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 4820a6353c2SJames Wright const CeedScalar P0 = stg_ctx->P0; 483a2d72b6fSJames Wright const CeedScalar rho = P0 / (GasConstant(gas) * theta0); 484a2d72b6fSJames Wright const CeedScalar nu = gas->mu / rho; 4850a6353c2SJames Wright 4862b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 4870a6353c2SJames Wright const CeedScalar x[] = {coords[0][i], coords[1][i], coords[2][i]}; 4880a6353c2SJames Wright const CeedScalar dXdx[2][3] = { 48966170c20SJames Wright {dXdx_q[0][0][i], dXdx_q[0][1][i], dXdx_q[0][2][i]}, 49066170c20SJames Wright {dXdx_q[1][0][i], dXdx_q[1][1][i], dXdx_q[1][2][i]}, 4910a6353c2SJames Wright }; 4920a6353c2SJames Wright 493831dbe9eSJames Wright CeedScalar h_node_sep[3]; 494831dbe9eSJames Wright h_node_sep[0] = dx; 495831dbe9eSJames Wright for (CeedInt j = 1; j < 3; j++) h_node_sep[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 496831dbe9eSJames Wright ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3); 4970a6353c2SJames Wright 4980a6353c2SJames Wright InterpolateProfile(coords[1][i], ubar, cij, &eps, <, stg_ctx); 4990a6353c2SJames Wright if (!mean_only) { 50062e628f8SJames Wright if (1) { 501831dbe9eSJames Wright StgShur14Calc_PrecompEktot(x, time, ubar, cij, inv_Ektotal[i], h_node_sep, x[1], eps, lt, nu, u, stg_ctx); 50262e628f8SJames Wright } else { // Original way 50362e628f8SJames Wright CeedScalar qn[STG_NMODES_MAX]; 504831dbe9eSJames Wright CalcSpectrum(coords[1][i], eps, lt, h_node_sep, nu, qn, stg_ctx); 505cbef7084SJames Wright StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx); 50662e628f8SJames Wright } 5070a6353c2SJames Wright } else { 5080a6353c2SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 5090a6353c2SJames Wright } 5100a6353c2SJames Wright 511a2d72b6fSJames Wright CeedScalar Y[5] = {P0, u[0], u[1], u[2], theta0}, q[5]; 512a2d72b6fSJames Wright State s = StateFromY(gas, Y); 513a2d72b6fSJames Wright StateToQ(gas, s, q, gas->state_var); 514a2d72b6fSJames Wright switch (gas->state_var) { 51597baf651SJames Wright case STATEVAR_CONSERVATIVE: 516a2d72b6fSJames Wright q[4] = 0.; // Don't set energy 51797baf651SJames Wright break; 51897baf651SJames Wright case STATEVAR_PRIMITIVE: 519a2d72b6fSJames Wright q[0] = 0; // Don't set pressure 52097baf651SJames Wright break; 521a2d72b6fSJames Wright case STATEVAR_ENTROPY: 522a2d72b6fSJames Wright q[0] = 0; // Don't set V_density 523a2d72b6fSJames Wright break; 524a2d72b6fSJames Wright } 525a2d72b6fSJames Wright for (CeedInt j = 0; j < 5; j++) { 526a2d72b6fSJames Wright bcval[j][i] = scale[i] * q[j]; 5270a6353c2SJames Wright } 5287c4551aaSJames Wright } 5290a6353c2SJames Wright return 0; 5300a6353c2SJames Wright } 531