1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2493642f1SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3493642f1SJames Wright // 4493642f1SJames Wright // SPDX-License-Identifier: BSD-2-Clause 5493642f1SJames Wright // 6493642f1SJames Wright // This file is part of CEED: http://github.com/ceed 7493642f1SJames Wright 8493642f1SJames Wright /// @file 9493642f1SJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm 10493642f1SJames Wright /// presented in Shur et al. 2014 11493642f1SJames Wright // 1204e40bb6SJeremy L Thompson /// SetupSTG_Rand reads in the input files and fills in STGShur14Context. 1304e40bb6SJeremy L Thompson /// Then STGShur14_CalcQF is run over quadrature points. 1404e40bb6SJeremy L Thompson /// Before the program exits, TearDownSTG is run to free the memory of the allocated arrays. 15493642f1SJames Wright #include <ceed.h> 16d0cce58aSJeremy L Thompson #include <math.h> 17493642f1SJames Wright #include <stdlib.h> 182b916ea7SJeremy L Thompson 193d65b166SJames Wright #include "newtonian_state.h" 201a74fa30SJames Wright #include "setupgeo_helpers.h" 21493642f1SJames Wright #include "stg_shur14_type.h" 22704b8bbeSJames Wright #include "utils.h" 23493642f1SJames Wright 24493642f1SJames Wright #define STG_NMODES_MAX 1024 25493642f1SJames Wright 26493642f1SJames Wright /* 27493642f1SJames Wright * @brief Interpolate quantities from input profile to given location 28493642f1SJames Wright * 29c77f3192SJames Wright * Assumed that prof_wd[i+1] > prof_wd[i] and prof_wd[0] = 0 30c77f3192SJames Wright * If wall_dist > prof_wd[-1], then the interpolation takes the values at prof_wd[-1] 31493642f1SJames Wright * 32c77f3192SJames Wright * @param[in] wall_dist Distance to the nearest wall 33c77f3192SJames Wright * @param[out] ubar Mean velocity at wall_dist 34c77f3192SJames Wright * @param[out] cij Cholesky decomposition at wall_dist 35c77f3192SJames Wright * @param[out] eps Turbulent dissipation at wall_dist 36c77f3192SJames Wright * @param[out] lt Turbulent length scale at wall_dist 37493642f1SJames Wright * @param[in] stg_ctx STGShur14Context for the problem 38493642f1SJames Wright */ 392b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar wall_dist, CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt, 4042454adaSJames Wright const StgShur14Context stg_ctx) { 41493642f1SJames Wright const CeedInt nprofs = stg_ctx->nprofs; 42c77f3192SJames Wright const CeedScalar *prof_wd = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 43493642f1SJames Wright const CeedScalar *prof_eps = &stg_ctx->data[stg_ctx->offsets.eps]; 44493642f1SJames Wright const CeedScalar *prof_lt = &stg_ctx->data[stg_ctx->offsets.lt]; 45493642f1SJames Wright const CeedScalar *prof_ubar = &stg_ctx->data[stg_ctx->offsets.ubar]; 46493642f1SJames Wright const CeedScalar *prof_cij = &stg_ctx->data[stg_ctx->offsets.cij]; 47493642f1SJames Wright CeedInt idx = -1; 48493642f1SJames Wright 49493642f1SJames Wright for (CeedInt i = 0; i < nprofs; i++) { 50c77f3192SJames Wright if (wall_dist < prof_wd[i]) { 51493642f1SJames Wright idx = i; 52493642f1SJames Wright break; 53493642f1SJames Wright } 54493642f1SJames Wright } 55493642f1SJames Wright 56c77f3192SJames Wright if (idx > 0) { // y within the bounds of prof_wd 57c77f3192SJames Wright CeedScalar coeff = (wall_dist - prof_wd[idx - 1]) / (prof_wd[idx] - prof_wd[idx - 1]); 58c77f3192SJames Wright 59493642f1SJames Wright ubar[0] = prof_ubar[0 * nprofs + idx - 1] + coeff * (prof_ubar[0 * nprofs + idx] - prof_ubar[0 * nprofs + idx - 1]); 60493642f1SJames Wright ubar[1] = prof_ubar[1 * nprofs + idx - 1] + coeff * (prof_ubar[1 * nprofs + idx] - prof_ubar[1 * nprofs + idx - 1]); 61493642f1SJames Wright ubar[2] = prof_ubar[2 * nprofs + idx - 1] + coeff * (prof_ubar[2 * nprofs + idx] - prof_ubar[2 * nprofs + idx - 1]); 62493642f1SJames Wright cij[0] = prof_cij[0 * nprofs + idx - 1] + coeff * (prof_cij[0 * nprofs + idx] - prof_cij[0 * nprofs + idx - 1]); 63493642f1SJames Wright cij[1] = prof_cij[1 * nprofs + idx - 1] + coeff * (prof_cij[1 * nprofs + idx] - prof_cij[1 * nprofs + idx - 1]); 64493642f1SJames Wright cij[2] = prof_cij[2 * nprofs + idx - 1] + coeff * (prof_cij[2 * nprofs + idx] - prof_cij[2 * nprofs + idx - 1]); 65493642f1SJames Wright cij[3] = prof_cij[3 * nprofs + idx - 1] + coeff * (prof_cij[3 * nprofs + idx] - prof_cij[3 * nprofs + idx - 1]); 66493642f1SJames Wright cij[4] = prof_cij[4 * nprofs + idx - 1] + coeff * (prof_cij[4 * nprofs + idx] - prof_cij[4 * nprofs + idx - 1]); 67493642f1SJames Wright cij[5] = prof_cij[5 * nprofs + idx - 1] + coeff * (prof_cij[5 * nprofs + idx] - prof_cij[5 * nprofs + idx - 1]); 68493642f1SJames Wright *eps = prof_eps[idx - 1] + coeff * (prof_eps[idx] - prof_eps[idx - 1]); 69493642f1SJames Wright *lt = prof_lt[idx - 1] + coeff * (prof_lt[idx] - prof_lt[idx - 1]); 70c77f3192SJames Wright } else { // y outside bounds of prof_wd 71493642f1SJames Wright ubar[0] = prof_ubar[1 * nprofs - 1]; 72493642f1SJames Wright ubar[1] = prof_ubar[2 * nprofs - 1]; 73493642f1SJames Wright ubar[2] = prof_ubar[3 * nprofs - 1]; 74493642f1SJames Wright cij[0] = prof_cij[1 * nprofs - 1]; 75493642f1SJames Wright cij[1] = prof_cij[2 * nprofs - 1]; 76493642f1SJames Wright cij[2] = prof_cij[3 * nprofs - 1]; 77493642f1SJames Wright cij[3] = prof_cij[4 * nprofs - 1]; 78493642f1SJames Wright cij[4] = prof_cij[5 * nprofs - 1]; 79493642f1SJames Wright cij[5] = prof_cij[6 * nprofs - 1]; 80493642f1SJames Wright *eps = prof_eps[nprofs - 1]; 81493642f1SJames Wright *lt = prof_lt[nprofs - 1]; 82493642f1SJames Wright } 83493642f1SJames Wright } 84493642f1SJames Wright 85493642f1SJames Wright /* 8671cd6200SJames Wright * @brief Calculate spectrum coefficient, qn 8771cd6200SJames Wright * 8871cd6200SJames Wright * Calculates q_n at a given distance to the wall 8971cd6200SJames Wright * 9071cd6200SJames Wright * @param[in] kappa nth wavenumber 9171cd6200SJames Wright * @param[in] dkappa Difference between wavenumbers 9271cd6200SJames Wright * @param[in] keta Dissipation wavenumber 9371cd6200SJames Wright * @param[in] kcut Mesh-induced cutoff wavenumber 9471cd6200SJames Wright * @param[in] ke Energy-containing wavenumber 959ef62cddSJames Wright * @param[in] Ektot_inv Inverse of total turbulent kinetic energy of spectrum 9671cd6200SJames Wright * @returns qn Spectrum coefficient 9771cd6200SJames Wright */ 982b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER CeedScalar Calc_qn(const CeedScalar kappa, const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut, 9970b0cb14SJames Wright const CeedScalar ke, const CeedScalar Ektot_inv) { 1002b916ea7SJeremy L Thompson const CeedScalar feta_x_fcut = exp(-Square(12 * kappa / keta) - Cube(4 * Max(kappa - 0.9 * kcut, 0) / kcut)); 1012b916ea7SJeremy L Thompson return pow(kappa / ke, 4.) * pow(1 + 2.4 * Square(kappa / ke), -17. / 6) * feta_x_fcut * dkappa * Ektot_inv; 10271cd6200SJames Wright } 10371cd6200SJames Wright 10471cd6200SJames Wright // Calculate hmax, ke, keta, and kcut 1052b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void SpectrumConstants(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 1062b916ea7SJeremy L Thompson const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke, CeedScalar *keta, CeedScalar *kcut) { 10771cd6200SJames Wright *hmax = Max(Max(h[0], h[1]), h[2]); 108c77f3192SJames Wright *ke = wall_dist == 0 ? 1e16 : 2 * M_PI / Min(2 * wall_dist, 3 * lt); 10971cd6200SJames Wright *keta = 2 * M_PI * pow(Cube(nu) / eps, -0.25); 110c77f3192SJames Wright *kcut = M_PI / Min(Max(Max(h[1], h[2]), 0.3 * (*hmax)) + 0.1 * wall_dist, *hmax); 11171cd6200SJames Wright } 11271cd6200SJames Wright 11371cd6200SJames Wright /* 114493642f1SJames Wright * @brief Calculate spectrum coefficients for STG 115493642f1SJames Wright * 116493642f1SJames Wright * Calculates q_n at a given distance to the wall 117493642f1SJames Wright * 118c77f3192SJames Wright * @param[in] wall_dist Distance to the nearest wall 119c77f3192SJames Wright * @param[in] eps Turbulent dissipation w/rt wall_dist 120c77f3192SJames Wright * @param[in] lt Turbulent length scale w/rt wall_dist 121493642f1SJames Wright * @param[in] h Element lengths in coordinate directions 122493642f1SJames Wright * @param[in] nu Dynamic Viscosity; 123493642f1SJames Wright * @param[in] stg_ctx STGShur14Context for the problem 124493642f1SJames Wright * @param[out] qn Spectrum coefficients, [nmodes] 125493642f1SJames Wright */ 1262b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void CalcSpectrum(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 12742454adaSJames Wright const CeedScalar nu, CeedScalar qn[], const StgShur14Context stg_ctx) { 128493642f1SJames Wright const CeedInt nmodes = stg_ctx->nmodes; 129493642f1SJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 13071cd6200SJames Wright CeedScalar hmax, ke, keta, kcut, Ektot = 0.0; 1312b916ea7SJeremy L Thompson 132c77f3192SJames Wright SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 133493642f1SJames Wright 134493642f1SJames Wright for (CeedInt n = 0; n < nmodes; n++) { 13571cd6200SJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 13671cd6200SJames Wright qn[n] = Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 137493642f1SJames Wright Ektot += qn[n]; 138493642f1SJames Wright } 139493642f1SJames Wright 1400a8dc919SJames Wright if (Ektot == 0) return; 141493642f1SJames Wright for (CeedInt n = 0; n < nmodes; n++) qn[n] /= Ektot; 142493642f1SJames Wright } 143493642f1SJames Wright 144493642f1SJames Wright /****************************************************** 145493642f1SJames Wright * @brief Calculate u(x,t) for STG inflow condition 146493642f1SJames Wright * 147493642f1SJames Wright * @param[in] X Location to evaluate u(X,t) 148493642f1SJames Wright * @param[in] t Time to evaluate u(X,t) 149493642f1SJames Wright * @param[in] ubar Mean velocity at X 150493642f1SJames Wright * @param[in] cij Cholesky decomposition at X 151493642f1SJames Wright * @param[in] qn Wavemode amplitudes at X, [nmodes] 152493642f1SJames Wright * @param[out] u Velocity at X and t 153493642f1SJames Wright * @param[in] stg_ctx STGShur14Context for the problem 154493642f1SJames Wright */ 15542454adaSJames Wright CEED_QFUNCTION_HELPER void StgShur14Calc(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 15642454adaSJames Wright const CeedScalar qn[], CeedScalar u[3], const StgShur14Context stg_ctx) { 157493642f1SJames Wright const CeedInt nmodes = stg_ctx->nmodes; 158493642f1SJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 159493642f1SJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 160493642f1SJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 161493642f1SJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 162493642f1SJames Wright CeedScalar xdotd, vp[3] = {0.}; 163493642f1SJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 164493642f1SJames Wright 1652b916ea7SJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 166493642f1SJames Wright xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 167493642f1SJames Wright xdotd = 0.; 168493642f1SJames Wright for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 169493642f1SJames Wright const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 1700a8dc919SJames Wright vp[0] += sqrt(qn[n]) * sigma[0 * nmodes + n] * cos_kxdp; 1710a8dc919SJames Wright vp[1] += sqrt(qn[n]) * sigma[1 * nmodes + n] * cos_kxdp; 1720a8dc919SJames Wright vp[2] += sqrt(qn[n]) * sigma[2 * nmodes + n] * cos_kxdp; 173493642f1SJames Wright } 1740a8dc919SJames Wright for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 175493642f1SJames Wright 176493642f1SJames Wright u[0] = ubar[0] + cij[0] * vp[0]; 177493642f1SJames Wright u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 178493642f1SJames Wright u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 179493642f1SJames Wright } 180493642f1SJames Wright 1818eea80fcSJames Wright /****************************************************** 1828eea80fcSJames Wright * @brief Calculate u(x,t) for STG inflow condition 1838eea80fcSJames Wright * 1848eea80fcSJames Wright * @param[in] X Location to evaluate u(X,t) 1858eea80fcSJames Wright * @param[in] t Time to evaluate u(X,t) 1868eea80fcSJames Wright * @param[in] ubar Mean velocity at X 1878eea80fcSJames Wright * @param[in] cij Cholesky decomposition at X 188c77f3192SJames Wright * @param[in] Ektot Total spectrum energy at this location 189c77f3192SJames Wright * @param[in] h Element size in 3 directions 190c77f3192SJames Wright * @param[in] wall_dist Distance to closest wall 191c77f3192SJames Wright * @param[in] eps Turbulent dissipation 192c77f3192SJames Wright * @param[in] lt Turbulent length scale 1938eea80fcSJames Wright * @param[out] u Velocity at X and t 1948eea80fcSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 1958eea80fcSJames Wright */ 19642454adaSJames Wright CEED_QFUNCTION_HELPER void StgShur14Calc_PrecompEktot(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 19742454adaSJames Wright const CeedScalar Ektot, const CeedScalar h[3], const CeedScalar wall_dist, const CeedScalar eps, 19842454adaSJames Wright const CeedScalar lt, const CeedScalar nu, CeedScalar u[3], const StgShur14Context stg_ctx) { 1998eea80fcSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 2008eea80fcSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 2018eea80fcSJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 2028eea80fcSJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 2038eea80fcSJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 2048eea80fcSJames Wright CeedScalar hmax, ke, keta, kcut; 205c77f3192SJames Wright SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 2068eea80fcSJames Wright CeedScalar xdotd, vp[3] = {0.}; 2078eea80fcSJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 2088eea80fcSJames Wright 2092b916ea7SJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 2108eea80fcSJames Wright xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1); 2118eea80fcSJames Wright xdotd = 0.; 2128eea80fcSJames Wright for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i]; 2138eea80fcSJames Wright const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]); 2148eea80fcSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 2158eea80fcSJames Wright const CeedScalar qn = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot); 2168eea80fcSJames Wright vp[0] += sqrt(qn) * sigma[0 * nmodes + n] * cos_kxdp; 2178eea80fcSJames Wright vp[1] += sqrt(qn) * sigma[1 * nmodes + n] * cos_kxdp; 2188eea80fcSJames Wright vp[2] += sqrt(qn) * sigma[2 * nmodes + n] * cos_kxdp; 2198eea80fcSJames Wright } 2208eea80fcSJames Wright for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5); 2218eea80fcSJames Wright 2228eea80fcSJames Wright u[0] = ubar[0] + cij[0] * vp[0]; 2238eea80fcSJames Wright u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1]; 2248eea80fcSJames Wright u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2]; 2258eea80fcSJames Wright } 2268eea80fcSJames Wright 22770b0cb14SJames Wright // Create preprocessed input for the stg calculation 22870b0cb14SJames Wright // 22970b0cb14SJames Wright // stg_data[0] = 1 / Ektot (inverse of total spectrum energy) 23042454adaSJames Wright CEED_QFUNCTION(StgShur14Preprocess)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 2316f188493SJames Wright const CeedScalar(*dXdx_q)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0]; 2323d65b166SJames Wright const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 2338eea80fcSJames Wright 2348eea80fcSJames Wright CeedScalar(*stg_data) = (CeedScalar(*))out[0]; 2358eea80fcSJames Wright 2368eea80fcSJames Wright CeedScalar ubar[3], cij[6], eps, lt; 23742454adaSJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 2388eea80fcSJames Wright const CeedScalar dx = stg_ctx->dx; 2398eea80fcSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 2408eea80fcSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 2418eea80fcSJames Wright const CeedScalar P0 = stg_ctx->P0; 2423d65b166SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 2438eea80fcSJames Wright const CeedScalar rho = P0 / (Rd * theta0); 2448eea80fcSJames Wright const CeedScalar nu = mu / rho; 2458eea80fcSJames Wright 2468eea80fcSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 2478eea80fcSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 2489eeef72bSJames Wright CeedScalar hmax, ke, keta, kcut; 2498eea80fcSJames Wright 2502b916ea7SJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 251c77f3192SJames Wright const CeedScalar wall_dist = x[1][i]; 2528eea80fcSJames Wright const CeedScalar dXdx[2][3] = { 2536f188493SJames Wright {dXdx_q[0][0][i], dXdx_q[0][1][i], dXdx_q[0][2][i]}, 2546f188493SJames Wright {dXdx_q[1][0][i], dXdx_q[1][1][i], dXdx_q[1][2][i]}, 2558eea80fcSJames Wright }; 2568eea80fcSJames Wright 2578eea80fcSJames Wright CeedScalar h[3]; 2588eea80fcSJames Wright h[0] = dx; 2592b916ea7SJeremy 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]); 2608eea80fcSJames Wright 261c77f3192SJames Wright InterpolateProfile(wall_dist, ubar, cij, &eps, <, stg_ctx); 262c77f3192SJames Wright SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 2638eea80fcSJames Wright 2648eea80fcSJames Wright // Calculate total TKE per spectrum 2652f638ed2SJames Wright CeedScalar Ek_tot = 0; 2662b916ea7SJeremy L Thompson CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) { 2678eea80fcSJames Wright const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1]; 2682f638ed2SJames Wright Ek_tot += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 2698eea80fcSJames Wright } 2702f638ed2SJames Wright // avoid underflowed and poorly defined spectrum coefficients 2712f638ed2SJames Wright stg_data[i] = Ek_tot != 0 ? 1 / Ek_tot : 0; 2728eea80fcSJames Wright } 2738eea80fcSJames Wright return 0; 2748eea80fcSJames Wright } 2758eea80fcSJames Wright 27643bff553SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition 27742454adaSJames Wright CEED_QFUNCTION(ICsStg)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 2783d65b166SJames Wright const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 2794f0244d1SJeremy L Thompson const CeedScalar(*J)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[1]; 28043bff553SJames Wright CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 28143bff553SJames Wright 28242454adaSJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 283d4e0f297SJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 284d4e0f297SJames Wright const CeedScalar dx = stg_ctx->dx; 285d4e0f297SJames Wright const CeedScalar time = stg_ctx->time; 28643bff553SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 28743bff553SJames Wright const CeedScalar P0 = stg_ctx->P0; 28843bff553SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 2893d65b166SJames Wright const CeedScalar rho = P0 / (GasConstant(&stg_ctx->newtonian_ctx) * theta0); 2903d65b166SJames Wright const CeedScalar nu = stg_ctx->newtonian_ctx.mu / rho; 29143bff553SJames Wright 2922b916ea7SJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 293d4e0f297SJames Wright const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]}; 2944f0244d1SJeremy L Thompson CeedScalar dXdx[3][3]; 2951a74fa30SJames Wright InvertMappingJacobian_3D(Q, i, J, dXdx, NULL); 296d4e0f297SJames Wright CeedScalar h[3]; 297d4e0f297SJames Wright h[0] = dx; 2982b916ea7SJeremy 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])); 299d4e0f297SJames Wright 300d4e0f297SJames Wright InterpolateProfile(x_i[1], ubar, cij, &eps, <, stg_ctx); 301d4e0f297SJames Wright if (stg_ctx->use_fluctuating_IC) { 3023d65b166SJames Wright CalcSpectrum(x_i[1], eps, lt, h, nu, qn, stg_ctx); 30342454adaSJames Wright StgShur14Calc(x_i, time, ubar, cij, qn, u, stg_ctx); 304d4e0f297SJames Wright } else { 305d4e0f297SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 306d4e0f297SJames Wright } 30743bff553SJames Wright 3083636f6a4SJames Wright switch (stg_ctx->newtonian_ctx.state_var) { 3093636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 31043bff553SJames Wright q0[0][i] = rho; 31143bff553SJames Wright q0[1][i] = u[0] * rho; 31243bff553SJames Wright q0[2][i] = u[1] * rho; 31343bff553SJames Wright q0[3][i] = u[2] * rho; 31443bff553SJames Wright q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0); 3153636f6a4SJames Wright break; 3163636f6a4SJames Wright 3173636f6a4SJames Wright case STATEVAR_PRIMITIVE: 3183636f6a4SJames Wright q0[0][i] = P0; 3193636f6a4SJames Wright q0[1][i] = u[0]; 3203636f6a4SJames Wright q0[2][i] = u[1]; 3213636f6a4SJames Wright q0[3][i] = u[2]; 3223636f6a4SJames Wright q0[4][i] = theta0; 3233636f6a4SJames Wright break; 32488243482SJames Wright } 325*b193fadcSJames Wright } 32643bff553SJames Wright return 0; 32743bff553SJames Wright } 32843bff553SJames Wright 329493642f1SJames Wright /******************************************************************** 330493642f1SJames Wright * @brief QFunction to calculate the inflow boundary condition 331493642f1SJames Wright * 332493642f1SJames Wright * This will loop through quadrature points, calculate the wavemode amplitudes 333493642f1SJames Wright * at each location, then calculate the actual velocity. 334493642f1SJames Wright */ 33542454adaSJames Wright CEED_QFUNCTION(StgShur14Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 3363d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 337ade49511SJames Wright const CeedScalar(*q_data_sur) = in[2]; 3383d65b166SJames Wright const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 339493642f1SJames Wright 3403d65b166SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 341ade49511SJames Wright CeedScalar(*jac_data_sur) = out[1]; 342493642f1SJames Wright 34342454adaSJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 344493642f1SJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 345493642f1SJames Wright const bool is_implicit = stg_ctx->is_implicit; 346493642f1SJames Wright const bool mean_only = stg_ctx->mean_only; 347493642f1SJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 348493642f1SJames Wright const CeedScalar dx = stg_ctx->dx; 349493642f1SJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 350493642f1SJames Wright const CeedScalar time = stg_ctx->time; 351493642f1SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 352493642f1SJames Wright const CeedScalar P0 = stg_ctx->P0; 353493642f1SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 3543d65b166SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 3553d65b166SJames Wright const CeedScalar gamma = HeatCapacityRatio(&stg_ctx->newtonian_ctx); 356493642f1SJames Wright 3572b916ea7SJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 358493642f1SJames Wright const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0); 359493642f1SJames Wright const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 360ade49511SJames Wright CeedScalar wdetJb, dXdx[2][3], norm[3]; 361ade49511SJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, dXdx, norm); 362ade49511SJames Wright wdetJb *= is_implicit ? -1. : 1.; 363493642f1SJames Wright 364493642f1SJames Wright CeedScalar h[3]; 365493642f1SJames Wright h[0] = dx; 3662b916ea7SJeremy L Thompson for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 367493642f1SJames Wright 368493642f1SJames Wright InterpolateProfile(X[1][i], ubar, cij, &eps, <, stg_ctx); 369493642f1SJames Wright if (!mean_only) { 370493642f1SJames Wright CalcSpectrum(X[1][i], eps, lt, h, mu / rho, qn, stg_ctx); 37142454adaSJames Wright StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx); 372493642f1SJames Wright } else { 373493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 374493642f1SJames Wright } 375493642f1SJames Wright 376a6e8f989SJames Wright const CeedScalar E_kinetic = .5 * rho * Dot3(u, u); 377493642f1SJames Wright CeedScalar E_internal, P; 378493642f1SJames Wright if (prescribe_T) { 379493642f1SJames Wright // Temperature is being set weakly (theta0) and for constant cv this sets E_internal 380493642f1SJames Wright E_internal = rho * cv * theta0; 381493642f1SJames Wright // Find pressure using 382493642f1SJames Wright P = rho * Rd * theta0; // interior rho with exterior T 383493642f1SJames Wright } else { 384493642f1SJames Wright E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution 385493642f1SJames Wright P = E_internal * (gamma - 1.); 386493642f1SJames Wright } 387493642f1SJames Wright 388493642f1SJames Wright const CeedScalar E = E_internal + E_kinetic; 389493642f1SJames Wright 390493642f1SJames Wright // Velocity normal to the boundary 391a6e8f989SJames Wright const CeedScalar u_normal = Dot3(norm, u); 392a6e8f989SJames Wright 393493642f1SJames Wright // The Physics 394493642f1SJames Wright // Zero v so all future terms can safely sum into it 395493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 396493642f1SJames Wright 397493642f1SJames Wright // The Physics 398493642f1SJames Wright // -- Density 399493642f1SJames Wright v[0][i] -= wdetJb * rho * u_normal; 400493642f1SJames Wright 401493642f1SJames Wright // -- Momentum 4022b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 403493642f1SJames Wright 404493642f1SJames Wright // -- Total Energy Density 405493642f1SJames Wright v[4][i] -= wdetJb * u_normal * (E + P); 406a6e8f989SJames Wright 407ade49511SJames Wright const CeedScalar U[] = {rho, u[0], u[1], u[2], E}, kmstress[6] = {0.}; 408ade49511SJames Wright StoredValuesPack(Q, i, 0, 5, U, jac_data_sur); 409ade49511SJames Wright StoredValuesPack(Q, i, 5, 6, kmstress, jac_data_sur); 410493642f1SJames Wright } 411493642f1SJames Wright return 0; 412493642f1SJames Wright } 413493642f1SJames Wright 41442454adaSJames Wright CEED_QFUNCTION(StgShur14Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 4153d65b166SJames Wright const CeedScalar(*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 4163d65b166SJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 4173d65b166SJames Wright const CeedScalar(*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 418a6e8f989SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 4193d65b166SJames Wright 42042454adaSJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 421a6e8f989SJames Wright const bool implicit = stg_ctx->is_implicit; 422a6e8f989SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 4233d65b166SJames Wright const CeedScalar Rd = GasConstant(&stg_ctx->newtonian_ctx); 4243d65b166SJames Wright const CeedScalar gamma = HeatCapacityRatio(&stg_ctx->newtonian_ctx); 425a6e8f989SJames Wright 426a6e8f989SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 427a6e8f989SJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 428a6e8f989SJames Wright 429*b193fadcSJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 430a6e8f989SJames Wright // Setup 431a6e8f989SJames Wright // -- Interp-to-Interp q_data 432a6e8f989SJames Wright // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 433a6e8f989SJames Wright // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 434a6e8f989SJames Wright // We can effect this by swapping the sign on this weight 435a6e8f989SJames Wright const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 436a6e8f989SJames Wright 437a6e8f989SJames Wright // Calculate inflow values 438a6e8f989SJames Wright CeedScalar velocity[3]; 439a6e8f989SJames Wright for (CeedInt j = 0; j < 3; j++) velocity[j] = jac_data_sur[5 + j][i]; 440ade49511SJames Wright // TODO This is almost certainly a bug. Velocity isn't stored here, only 0s. 441a6e8f989SJames Wright 442a6e8f989SJames Wright // enabling user to choose between weak T and weak rho inflow 443a6e8f989SJames Wright CeedScalar drho, dE, dP; 444a6e8f989SJames Wright if (prescribe_T) { 445a6e8f989SJames Wright // rho should be from the current solution 446a6e8f989SJames Wright drho = dq[0][i]; 447a6e8f989SJames Wright CeedScalar dE_internal = drho * cv * theta0; 448a6e8f989SJames Wright CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity); 449a6e8f989SJames Wright dE = dE_internal + dE_kinetic; 450a6e8f989SJames Wright dP = drho * Rd * theta0; // interior rho with exterior T 451a6e8f989SJames Wright } else { // rho specified, E_internal from solution 452a6e8f989SJames Wright drho = 0; 453a6e8f989SJames Wright dE = dq[4][i]; 454a6e8f989SJames Wright dP = dE * (gamma - 1.); 455a6e8f989SJames Wright } 4562b916ea7SJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 457a6e8f989SJames Wright 458a6e8f989SJames Wright const CeedScalar u_normal = Dot3(norm, velocity); 459a6e8f989SJames Wright 460a6e8f989SJames Wright v[0][i] = -wdetJb * drho * u_normal; 4612b916ea7SJeremy L Thompson for (int j = 0; j < 3; j++) v[j + 1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP); 462a6e8f989SJames Wright v[4][i] = -wdetJb * u_normal * (dE + dP); 463*b193fadcSJames Wright } 464a6e8f989SJames Wright return 0; 465a6e8f989SJames Wright } 466a6e8f989SJames Wright 467b7190ff7SJames Wright /******************************************************************** 468b7190ff7SJames Wright * @brief QFunction to calculate the strongly enforce inflow BC 469b7190ff7SJames Wright * 470b7190ff7SJames Wright * This QF is for the strong application of STG via libCEED (rather than 471b7190ff7SJames Wright * through the native PETSc `DMAddBoundary` -> `bcFunc` method. 472b7190ff7SJames Wright */ 47342454adaSJames Wright CEED_QFUNCTION(StgShur14InflowStrongQF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 4746f188493SJames Wright const CeedScalar(*dXdx_q)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0]; 4753d65b166SJames Wright const CeedScalar(*coords)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 4763d65b166SJames Wright const CeedScalar(*scale) = (const CeedScalar(*))in[2]; 4779ef62cddSJames Wright const CeedScalar(*inv_Ektotal) = (const CeedScalar(*))in[3]; 478b7190ff7SJames Wright CeedScalar(*bcval)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 479b7190ff7SJames Wright 48042454adaSJames Wright const StgShur14Context stg_ctx = (StgShur14Context)ctx; 48170b0cb14SJames Wright CeedScalar u[3], ubar[3], cij[6], eps, lt; 482b7190ff7SJames Wright const bool mean_only = stg_ctx->mean_only; 483b7190ff7SJames Wright const CeedScalar dx = stg_ctx->dx; 484b7190ff7SJames Wright const CeedScalar time = stg_ctx->time; 485b7190ff7SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 486b7190ff7SJames Wright const CeedScalar P0 = stg_ctx->P0; 4873d65b166SJames Wright const CeedScalar rho = P0 / (GasConstant(&stg_ctx->newtonian_ctx) * theta0); 4883d65b166SJames Wright const CeedScalar nu = stg_ctx->newtonian_ctx.mu / rho; 489b7190ff7SJames Wright 4902b916ea7SJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 491b7190ff7SJames Wright const CeedScalar x[] = {coords[0][i], coords[1][i], coords[2][i]}; 492b7190ff7SJames Wright const CeedScalar dXdx[2][3] = { 4936f188493SJames Wright {dXdx_q[0][0][i], dXdx_q[0][1][i], dXdx_q[0][2][i]}, 4946f188493SJames Wright {dXdx_q[1][0][i], dXdx_q[1][1][i], dXdx_q[1][2][i]}, 495b7190ff7SJames Wright }; 496b7190ff7SJames Wright 497b7190ff7SJames Wright CeedScalar h[3]; 498b7190ff7SJames Wright h[0] = dx; 4992b916ea7SJeremy L Thompson for (CeedInt j = 1; j < 3; j++) h[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 500b7190ff7SJames Wright 501b7190ff7SJames Wright InterpolateProfile(coords[1][i], ubar, cij, &eps, <, stg_ctx); 502b7190ff7SJames Wright if (!mean_only) { 50370b0cb14SJames Wright if (1) { 50442454adaSJames Wright StgShur14Calc_PrecompEktot(x, time, ubar, cij, inv_Ektotal[i], h, x[1], eps, lt, nu, u, stg_ctx); 50570b0cb14SJames Wright } else { // Original way 50670b0cb14SJames Wright CeedScalar qn[STG_NMODES_MAX]; 5073d65b166SJames Wright CalcSpectrum(coords[1][i], eps, lt, h, nu, qn, stg_ctx); 50842454adaSJames Wright StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx); 50970b0cb14SJames Wright } 510b7190ff7SJames Wright } else { 511b7190ff7SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 512b7190ff7SJames Wright } 513b7190ff7SJames Wright 5143636f6a4SJames Wright switch (stg_ctx->newtonian_ctx.state_var) { 5153636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 516b7190ff7SJames Wright bcval[0][i] = scale[i] * rho; 517b7190ff7SJames Wright bcval[1][i] = scale[i] * rho * u[0]; 518b7190ff7SJames Wright bcval[2][i] = scale[i] * rho * u[1]; 519b7190ff7SJames Wright bcval[3][i] = scale[i] * rho * u[2]; 52066531c8bSJames Wright bcval[4][i] = 0.; 5213636f6a4SJames Wright break; 5223636f6a4SJames Wright 5233636f6a4SJames Wright case STATEVAR_PRIMITIVE: 5243636f6a4SJames Wright bcval[0][i] = 0; 5253636f6a4SJames Wright bcval[1][i] = scale[i] * u[0]; 5263636f6a4SJames Wright bcval[2][i] = scale[i] * u[1]; 5273636f6a4SJames Wright bcval[3][i] = scale[i] * u[2]; 5283636f6a4SJames Wright bcval[4][i] = scale[i] * theta0; 5293636f6a4SJames Wright break; 530b7190ff7SJames Wright } 53188243482SJames Wright } 532b7190ff7SJames Wright return 0; 533b7190ff7SJames Wright } 534