xref: /honee/qfunctions/stg_shur14.h (revision 475f0cac5d40259768f4556cf888e8f2448554cb)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3493642f1SJames Wright 
4493642f1SJames Wright /// @file
5493642f1SJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm
6493642f1SJames Wright /// presented in Shur et al. 2014
7493642f1SJames Wright //
804e40bb6SJeremy L Thompson /// SetupSTG_Rand reads in the input files and fills in STGShur14Context.
904e40bb6SJeremy L Thompson /// Then STGShur14_CalcQF is run over quadrature points.
1004e40bb6SJeremy L Thompson /// Before the program exits, TearDownSTG is run to free the memory of the allocated arrays.
113e17a7a1SJames Wright #include <ceed/types.h>
123e17a7a1SJames Wright #ifndef CEED_RUNNING_JIT_PASS
13d0cce58aSJeremy L Thompson #include <math.h>
143e17a7a1SJames Wright #endif
152b916ea7SJeremy L Thompson 
163d65b166SJames Wright #include "newtonian_state.h"
171a74fa30SJames Wright #include "setupgeo_helpers.h"
18493642f1SJames Wright #include "stg_shur14_type.h"
19704b8bbeSJames Wright #include "utils.h"
20493642f1SJames Wright 
21493642f1SJames Wright #define STG_NMODES_MAX 1024
22493642f1SJames Wright 
23493642f1SJames Wright /*
24493642f1SJames Wright  * @brief Interpolate quantities from input profile to given location
25493642f1SJames Wright  *
26c77f3192SJames Wright  * Assumed that prof_wd[i+1] > prof_wd[i] and prof_wd[0] = 0
27c77f3192SJames Wright  * If wall_dist > prof_wd[-1], then the interpolation takes the values at prof_wd[-1]
28493642f1SJames Wright  *
29c77f3192SJames Wright  * @param[in]  wall_dist Distance to the nearest wall
30c77f3192SJames Wright  * @param[out] ubar      Mean velocity at wall_dist
31c77f3192SJames Wright  * @param[out] cij       Cholesky decomposition at wall_dist
32c77f3192SJames Wright  * @param[out] eps       Turbulent dissipation at wall_dist
33c77f3192SJames Wright  * @param[out] lt        Turbulent length scale at wall_dist
34493642f1SJames Wright  * @param[in]  stg_ctx   STGShur14Context for the problem
35493642f1SJames Wright  */
InterpolateProfile(const CeedScalar wall_dist,CeedScalar ubar[3],CeedScalar cij[6],CeedScalar * eps,CeedScalar * lt,const StgShur14Context stg_ctx)362b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar wall_dist, CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt,
3742454adaSJames Wright                                               const StgShur14Context stg_ctx) {
38493642f1SJames Wright   const CeedInt     nprofs    = stg_ctx->nprofs;
39c77f3192SJames Wright   const CeedScalar *prof_wd   = &stg_ctx->data[stg_ctx->offsets.wall_dist];
40493642f1SJames Wright   const CeedScalar *prof_eps  = &stg_ctx->data[stg_ctx->offsets.eps];
41493642f1SJames Wright   const CeedScalar *prof_lt   = &stg_ctx->data[stg_ctx->offsets.lt];
42493642f1SJames Wright   const CeedScalar *prof_ubar = &stg_ctx->data[stg_ctx->offsets.ubar];
43493642f1SJames Wright   const CeedScalar *prof_cij  = &stg_ctx->data[stg_ctx->offsets.cij];
44493642f1SJames Wright   CeedInt           idx       = -1;
45493642f1SJames Wright 
46493642f1SJames Wright   for (CeedInt i = 0; i < nprofs; i++) {
47c77f3192SJames Wright     if (wall_dist < prof_wd[i]) {
48493642f1SJames Wright       idx = i;
49493642f1SJames Wright       break;
50493642f1SJames Wright     }
51493642f1SJames Wright   }
52493642f1SJames Wright 
53c77f3192SJames Wright   if (idx > 0) {  // y within the bounds of prof_wd
54c77f3192SJames Wright     CeedScalar coeff = (wall_dist - prof_wd[idx - 1]) / (prof_wd[idx] - prof_wd[idx - 1]);
55c77f3192SJames Wright 
56493642f1SJames Wright     ubar[0] = prof_ubar[0 * nprofs + idx - 1] + coeff * (prof_ubar[0 * nprofs + idx] - prof_ubar[0 * nprofs + idx - 1]);
57493642f1SJames Wright     ubar[1] = prof_ubar[1 * nprofs + idx - 1] + coeff * (prof_ubar[1 * nprofs + idx] - prof_ubar[1 * nprofs + idx - 1]);
58493642f1SJames Wright     ubar[2] = prof_ubar[2 * nprofs + idx - 1] + coeff * (prof_ubar[2 * nprofs + idx] - prof_ubar[2 * nprofs + idx - 1]);
59493642f1SJames Wright     cij[0]  = prof_cij[0 * nprofs + idx - 1] + coeff * (prof_cij[0 * nprofs + idx] - prof_cij[0 * nprofs + idx - 1]);
60493642f1SJames Wright     cij[1]  = prof_cij[1 * nprofs + idx - 1] + coeff * (prof_cij[1 * nprofs + idx] - prof_cij[1 * nprofs + idx - 1]);
61493642f1SJames Wright     cij[2]  = prof_cij[2 * nprofs + idx - 1] + coeff * (prof_cij[2 * nprofs + idx] - prof_cij[2 * nprofs + idx - 1]);
62493642f1SJames Wright     cij[3]  = prof_cij[3 * nprofs + idx - 1] + coeff * (prof_cij[3 * nprofs + idx] - prof_cij[3 * nprofs + idx - 1]);
63493642f1SJames Wright     cij[4]  = prof_cij[4 * nprofs + idx - 1] + coeff * (prof_cij[4 * nprofs + idx] - prof_cij[4 * nprofs + idx - 1]);
64493642f1SJames Wright     cij[5]  = prof_cij[5 * nprofs + idx - 1] + coeff * (prof_cij[5 * nprofs + idx] - prof_cij[5 * nprofs + idx - 1]);
65493642f1SJames Wright     *eps    = prof_eps[idx - 1] + coeff * (prof_eps[idx] - prof_eps[idx - 1]);
66493642f1SJames Wright     *lt     = prof_lt[idx - 1] + coeff * (prof_lt[idx] - prof_lt[idx - 1]);
67c77f3192SJames Wright   } else {  // y outside bounds of prof_wd
68493642f1SJames Wright     ubar[0] = prof_ubar[1 * nprofs - 1];
69493642f1SJames Wright     ubar[1] = prof_ubar[2 * nprofs - 1];
70493642f1SJames Wright     ubar[2] = prof_ubar[3 * nprofs - 1];
71493642f1SJames Wright     cij[0]  = prof_cij[1 * nprofs - 1];
72493642f1SJames Wright     cij[1]  = prof_cij[2 * nprofs - 1];
73493642f1SJames Wright     cij[2]  = prof_cij[3 * nprofs - 1];
74493642f1SJames Wright     cij[3]  = prof_cij[4 * nprofs - 1];
75493642f1SJames Wright     cij[4]  = prof_cij[5 * nprofs - 1];
76493642f1SJames Wright     cij[5]  = prof_cij[6 * nprofs - 1];
77493642f1SJames Wright     *eps    = prof_eps[nprofs - 1];
78493642f1SJames Wright     *lt     = prof_lt[nprofs - 1];
79493642f1SJames Wright   }
80493642f1SJames Wright }
81493642f1SJames Wright 
82493642f1SJames Wright /*
8371cd6200SJames Wright  * @brief Calculate spectrum coefficient, qn
8471cd6200SJames Wright  *
8571cd6200SJames Wright  * Calculates q_n at a given distance to the wall
8671cd6200SJames Wright  *
8771cd6200SJames Wright  * @param[in]  kappa     nth wavenumber
8871cd6200SJames Wright  * @param[in]  dkappa    Difference between wavenumbers
8971cd6200SJames Wright  * @param[in]  keta      Dissipation wavenumber
9071cd6200SJames Wright  * @param[in]  kcut      Mesh-induced cutoff wavenumber
9171cd6200SJames Wright  * @param[in]  ke        Energy-containing wavenumber
929ef62cddSJames Wright  * @param[in]  Ektot_inv Inverse of total turbulent kinetic energy of spectrum
9371cd6200SJames Wright  * @returns    qn        Spectrum coefficient
9471cd6200SJames Wright  */
Calc_qn(const CeedScalar kappa,const CeedScalar dkappa,const CeedScalar keta,const CeedScalar kcut,const CeedScalar ke,const CeedScalar Ektot_inv)952b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER CeedScalar Calc_qn(const CeedScalar kappa, const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut,
9670b0cb14SJames Wright                                          const CeedScalar ke, const CeedScalar Ektot_inv) {
972b916ea7SJeremy L Thompson   const CeedScalar feta_x_fcut = exp(-Square(12 * kappa / keta) - Cube(4 * Max(kappa - 0.9 * kcut, 0) / kcut));
982b916ea7SJeremy L Thompson   return pow(kappa / ke, 4.) * pow(1 + 2.4 * Square(kappa / ke), -17. / 6) * feta_x_fcut * dkappa * Ektot_inv;
9971cd6200SJames Wright }
10071cd6200SJames Wright 
10171cd6200SJames Wright // Calculate hmax, ke, keta, and kcut
SpectrumConstants(const CeedScalar wall_dist,const CeedScalar eps,const CeedScalar lt,const CeedScalar hNodSep[3],const CeedScalar nu,CeedScalar * hmax,CeedScalar * ke,CeedScalar * keta,CeedScalar * kcut)10284b557acSJames Wright CEED_QFUNCTION_HELPER void SpectrumConstants(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar hNodSep[3],
1032b916ea7SJeremy L Thompson                                              const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke, CeedScalar *keta, CeedScalar *kcut) {
10484b557acSJames Wright   *hmax = Max(Max(hNodSep[0], hNodSep[1]), hNodSep[2]);
105c77f3192SJames Wright   *ke   = wall_dist == 0 ? 1e16 : 2 * M_PI / Min(2 * wall_dist, 3 * lt);
10671cd6200SJames Wright   *keta = 2 * M_PI * pow(Cube(nu) / eps, -0.25);
10784b557acSJames Wright   *kcut = M_PI / Min(Max(Max(hNodSep[1], hNodSep[2]), 0.3 * (*hmax)) + 0.1 * wall_dist, *hmax);
10871cd6200SJames Wright }
10971cd6200SJames Wright 
11071cd6200SJames Wright /*
111493642f1SJames Wright  * @brief Calculate spectrum coefficients for STG
112493642f1SJames Wright  *
113493642f1SJames Wright  * Calculates q_n at a given distance to the wall
114493642f1SJames Wright  *
115c77f3192SJames Wright  * @param[in]  wall_dist  Distance to the nearest wall
116c77f3192SJames Wright  * @param[in]  eps        Turbulent dissipation w/rt wall_dist
117c77f3192SJames Wright  * @param[in]  lt         Turbulent length scale w/rt wall_dist
11884b557acSJames Wright  * @param[in]  h_node_sep Element lengths in coordinate directions
119493642f1SJames Wright  * @param[in]  nu         Dynamic Viscosity;
120493642f1SJames Wright  * @param[in]  stg_ctx    STGShur14Context for the problem
121493642f1SJames Wright  * @param[out] qn         Spectrum coefficients, [nmodes]
122493642f1SJames Wright  */
CalcSpectrum(const CeedScalar wall_dist,const CeedScalar eps,const CeedScalar lt,const CeedScalar h_node_sep[3],const CeedScalar nu,CeedScalar qn[],const StgShur14Context stg_ctx)12384b557acSJames Wright CEED_QFUNCTION_HELPER void CalcSpectrum(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h_node_sep[3],
12442454adaSJames Wright                                         const CeedScalar nu, CeedScalar qn[], const StgShur14Context stg_ctx) {
125493642f1SJames Wright   const CeedInt     nmodes = stg_ctx->nmodes;
126493642f1SJames Wright   const CeedScalar *kappa  = &stg_ctx->data[stg_ctx->offsets.kappa];
12771cd6200SJames Wright   CeedScalar        hmax, ke, keta, kcut, Ektot = 0.0;
1282b916ea7SJeremy L Thompson 
12984b557acSJames Wright   SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut);
130493642f1SJames Wright 
131493642f1SJames Wright   for (CeedInt n = 0; n < nmodes; n++) {
13271cd6200SJames Wright     const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1];
13371cd6200SJames Wright     qn[n]                   = Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0);
134493642f1SJames Wright     Ektot += qn[n];
135493642f1SJames Wright   }
136493642f1SJames Wright 
1370a8dc919SJames Wright   if (Ektot == 0) return;
138493642f1SJames Wright   for (CeedInt n = 0; n < nmodes; n++) qn[n] /= Ektot;
139493642f1SJames Wright }
140493642f1SJames Wright 
141493642f1SJames Wright /******************************************************
142493642f1SJames Wright  * @brief Calculate u(x,t) for STG inflow condition
143493642f1SJames Wright  *
144493642f1SJames Wright  * @param[in]  X       Location to evaluate u(X,t)
145493642f1SJames Wright  * @param[in]  t       Time to evaluate u(X,t)
146493642f1SJames Wright  * @param[in]  ubar    Mean velocity at X
147493642f1SJames Wright  * @param[in]  cij     Cholesky decomposition at X
148493642f1SJames Wright  * @param[in]  qn      Wavemode amplitudes at X, [nmodes]
149493642f1SJames Wright  * @param[out] u       Velocity at X and t
150493642f1SJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
151493642f1SJames Wright  */
StgShur14Calc(const CeedScalar X[3],const CeedScalar t,const CeedScalar ubar[3],const CeedScalar cij[6],const CeedScalar qn[],CeedScalar u[3],const StgShur14Context stg_ctx)15242454adaSJames Wright CEED_QFUNCTION_HELPER void StgShur14Calc(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6],
15342454adaSJames Wright                                          const CeedScalar qn[], CeedScalar u[3], const StgShur14Context stg_ctx) {
154493642f1SJames Wright   const CeedInt     nmodes = stg_ctx->nmodes;
155493642f1SJames Wright   const CeedScalar *kappa  = &stg_ctx->data[stg_ctx->offsets.kappa];
156493642f1SJames Wright   const CeedScalar *phi    = &stg_ctx->data[stg_ctx->offsets.phi];
157493642f1SJames Wright   const CeedScalar *sigma  = &stg_ctx->data[stg_ctx->offsets.sigma];
158493642f1SJames Wright   const CeedScalar *d      = &stg_ctx->data[stg_ctx->offsets.d];
159493642f1SJames Wright   CeedScalar        xdotd, vp[3] = {0.};
160493642f1SJames Wright   CeedScalar        xhat[] = {0., X[1], X[2]};
161493642f1SJames Wright 
1622b916ea7SJeremy L Thompson   CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) {
163493642f1SJames Wright     xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1);
164493642f1SJames Wright     xdotd   = 0.;
165493642f1SJames Wright     for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i];
166493642f1SJames Wright     const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]);
1670a8dc919SJames Wright     vp[0] += sqrt(qn[n]) * sigma[0 * nmodes + n] * cos_kxdp;
1680a8dc919SJames Wright     vp[1] += sqrt(qn[n]) * sigma[1 * nmodes + n] * cos_kxdp;
1690a8dc919SJames Wright     vp[2] += sqrt(qn[n]) * sigma[2 * nmodes + n] * cos_kxdp;
170493642f1SJames Wright   }
1710a8dc919SJames Wright   for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5);
172493642f1SJames Wright 
173493642f1SJames Wright   u[0] = ubar[0] + cij[0] * vp[0];
174493642f1SJames Wright   u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1];
175493642f1SJames Wright   u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2];
176493642f1SJames Wright }
177493642f1SJames Wright 
1788eea80fcSJames Wright /******************************************************
1798eea80fcSJames Wright  * @brief Calculate u(x,t) for STG inflow condition
1808eea80fcSJames Wright  *
1818eea80fcSJames Wright  * @param[in]  X          Location to evaluate u(X,t)
1828eea80fcSJames Wright  * @param[in]  t          Time to evaluate u(X,t)
1838eea80fcSJames Wright  * @param[in]  ubar       Mean velocity at X
1848eea80fcSJames Wright  * @param[in]  cij        Cholesky decomposition at X
185c77f3192SJames Wright  * @param[in]  Ektot      Total spectrum energy at this location
18684b557acSJames Wright  * @param[in]  h_node_sep Element size in 3 directions
187c77f3192SJames Wright  * @param[in]  wall_dist  Distance to closest wall
188c77f3192SJames Wright  * @param[in]  eps        Turbulent dissipation
189c77f3192SJames Wright  * @param[in]  lt         Turbulent length scale
1908eea80fcSJames Wright  * @param[out] u          Velocity at X and t
1918eea80fcSJames Wright  * @param[in]  stg_ctx    STGShur14Context for the problem
1928eea80fcSJames Wright  */
StgShur14Calc_PrecompEktot(const CeedScalar X[3],const CeedScalar t,const CeedScalar ubar[3],const CeedScalar cij[6],const CeedScalar Ektot,const CeedScalar h_node_sep[3],const CeedScalar wall_dist,const CeedScalar eps,const CeedScalar lt,const CeedScalar nu,CeedScalar u[3],const StgShur14Context stg_ctx)19342454adaSJames Wright CEED_QFUNCTION_HELPER void StgShur14Calc_PrecompEktot(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6],
19484b557acSJames Wright                                                       const CeedScalar Ektot, const CeedScalar h_node_sep[3], const CeedScalar wall_dist,
19584b557acSJames Wright                                                       const CeedScalar eps, const CeedScalar lt, const CeedScalar nu, CeedScalar u[3],
19684b557acSJames Wright                                                       const StgShur14Context stg_ctx) {
1978eea80fcSJames Wright   const CeedInt     nmodes = stg_ctx->nmodes;
1988eea80fcSJames Wright   const CeedScalar *kappa  = &stg_ctx->data[stg_ctx->offsets.kappa];
1998eea80fcSJames Wright   const CeedScalar *phi    = &stg_ctx->data[stg_ctx->offsets.phi];
2008eea80fcSJames Wright   const CeedScalar *sigma  = &stg_ctx->data[stg_ctx->offsets.sigma];
2018eea80fcSJames Wright   const CeedScalar *d      = &stg_ctx->data[stg_ctx->offsets.d];
2028eea80fcSJames Wright   CeedScalar        hmax, ke, keta, kcut;
20384b557acSJames Wright   SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut);
2048eea80fcSJames Wright   CeedScalar xdotd, vp[3] = {0.};
2058eea80fcSJames Wright   CeedScalar xhat[] = {0., X[1], X[2]};
2068eea80fcSJames Wright 
2072b916ea7SJeremy L Thompson   CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) {
2088eea80fcSJames Wright     xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1);
2098eea80fcSJames Wright     xdotd   = 0.;
2108eea80fcSJames Wright     for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i];
2118eea80fcSJames Wright     const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]);
2128eea80fcSJames Wright     const CeedScalar dkappa   = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1];
2138eea80fcSJames Wright     const CeedScalar qn       = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot);
2148eea80fcSJames Wright     vp[0] += sqrt(qn) * sigma[0 * nmodes + n] * cos_kxdp;
2158eea80fcSJames Wright     vp[1] += sqrt(qn) * sigma[1 * nmodes + n] * cos_kxdp;
2168eea80fcSJames Wright     vp[2] += sqrt(qn) * sigma[2 * nmodes + n] * cos_kxdp;
2178eea80fcSJames Wright   }
2188eea80fcSJames Wright   for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5);
2198eea80fcSJames Wright 
2208eea80fcSJames Wright   u[0] = ubar[0] + cij[0] * vp[0];
2218eea80fcSJames Wright   u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1];
2228eea80fcSJames Wright   u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2];
2238eea80fcSJames Wright }
2248eea80fcSJames Wright 
2258424c23aSJames Wright /**
2268424c23aSJames Wright    @brief Calculate the element length scales based on dXdx
2278424c23aSJames Wright 
2288424c23aSJames Wright    WARNING: This assumes the reference domain is [-1,1], which is not true for tetrahedral elements
2298424c23aSJames Wright 
2308424c23aSJames Wright    @param[in]  dXdx    Inverse mapping Jacobian, d\xi/dx
2318424c23aSJames Wright    @param[in]  scale   Scale factor for the element lengths
2328424c23aSJames Wright    @param[out] lengths The element lengths in each cartesian direction
2338424c23aSJames Wright **/
CalculateElementLengths(CeedScalar dXdx[3][3],CeedScalar scale,CeedScalar lengths[3])2348424c23aSJames Wright CEED_QFUNCTION_HELPER void CalculateElementLengths(CeedScalar dXdx[3][3], CeedScalar scale, CeedScalar lengths[3]) {
2358424c23aSJames Wright   for (CeedInt j = 0; j < 3; j++) lengths[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]) + Square(dXdx[2][j]));
2368424c23aSJames Wright   ScaleN(lengths, scale, 3);
2378424c23aSJames Wright }
2388424c23aSJames Wright 
23970b0cb14SJames Wright // Create preprocessed input for the stg calculation
24070b0cb14SJames Wright //
24170b0cb14SJames Wright // stg_data[0] = 1 / Ektot (inverse of total spectrum energy)
StgShur14Preprocess(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)24242454adaSJames Wright CEED_QFUNCTION(StgShur14Preprocess)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
24321ba7ba4SJames Wright   const CeedScalar *dXdx_q         = in[0];
2443d65b166SJames Wright   const CeedScalar(*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
2458eea80fcSJames Wright 
2468eea80fcSJames Wright   CeedScalar(*stg_data) = (CeedScalar(*))out[0];
2478eea80fcSJames Wright 
2488eea80fcSJames Wright   CeedScalar             ubar[3], cij[6], eps, lt;
24942454adaSJames Wright   const StgShur14Context stg_ctx = (StgShur14Context)ctx;
250*cde3d787SJames Wright   const CeedScalar       mu      = stg_ctx->newt_ctx.gas.mu;
2518eea80fcSJames Wright   const CeedScalar       theta0  = stg_ctx->theta0;
2528eea80fcSJames Wright   const CeedScalar       P0      = stg_ctx->P0;
253*cde3d787SJames Wright   const CeedScalar       Rd      = GasConstant(stg_ctx->newt_ctx.gas);
2548eea80fcSJames Wright   const CeedScalar       rho     = P0 / (Rd * theta0);
2558eea80fcSJames Wright   const CeedScalar       nu      = mu / rho;
2568eea80fcSJames Wright 
2578eea80fcSJames Wright   const CeedInt     nmodes = stg_ctx->nmodes;
2588eea80fcSJames Wright   const CeedScalar *kappa  = &stg_ctx->data[stg_ctx->offsets.kappa];
2599eeef72bSJames Wright   CeedScalar        hmax, ke, keta, kcut;
2608eea80fcSJames Wright 
2612b916ea7SJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
262c77f3192SJames Wright     const CeedScalar wall_dist = x[1][i];
2638424c23aSJames Wright     CeedScalar       dXdx[3][3], h_node_sep[3];
26421ba7ba4SJames Wright     StoredValuesUnpack(Q, i, 0, 9, dXdx_q, (CeedScalar *)dXdx);
2658eea80fcSJames Wright 
2668424c23aSJames Wright     CalculateElementLengths(dXdx, stg_ctx->h_scale_factor, h_node_sep);
267c77f3192SJames Wright     InterpolateProfile(wall_dist, ubar, cij, &eps, &lt, stg_ctx);
26884b557acSJames Wright     SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut);
2698eea80fcSJames Wright 
2708eea80fcSJames Wright     // Calculate total TKE per spectrum
2712f638ed2SJames Wright     CeedScalar Ek_tot = 0;
2722b916ea7SJeremy L Thompson     CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) {
2738eea80fcSJames Wright       const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1];
2742f638ed2SJames Wright       Ek_tot += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0);
2758eea80fcSJames Wright     }
2762f638ed2SJames Wright     // avoid underflowed and poorly defined spectrum coefficients
2772f638ed2SJames Wright     stg_data[i] = Ek_tot != 0 ? 1 / Ek_tot : 0;
2788eea80fcSJames Wright   }
2798eea80fcSJames Wright   return 0;
2808eea80fcSJames Wright }
2818eea80fcSJames Wright 
28243bff553SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition
ICsStg(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)28342454adaSJames Wright CEED_QFUNCTION(ICsStg)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
2843d65b166SJames Wright   const CeedScalar(*x)[CEED_Q_VLA]    = (const CeedScalar(*)[CEED_Q_VLA])in[0];
2854f0244d1SJeremy L Thompson   const CeedScalar(*J)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[1];
28643bff553SJames Wright   CeedScalar(*q0)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
28743bff553SJames Wright 
28842454adaSJames Wright   const StgShur14Context         stg_ctx = (StgShur14Context)ctx;
289*cde3d787SJames Wright   const NewtonianIdealGasContext context = &stg_ctx->newt_ctx;
290*cde3d787SJames Wright   const NewtonianIGProperties    gas     = context->gas;
291d4e0f297SJames Wright   CeedScalar                     qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
292d4e0f297SJames Wright   const CeedScalar               dx     = stg_ctx->dx;
293d4e0f297SJames Wright   const CeedScalar               time   = stg_ctx->time;
29443bff553SJames Wright   const CeedScalar               theta0 = stg_ctx->theta0;
29543bff553SJames Wright   const CeedScalar               P0     = stg_ctx->P0;
2969b103f75SJames Wright   const CeedScalar               rho    = P0 / (GasConstant(gas) * theta0);
297*cde3d787SJames Wright   const CeedScalar               nu     = gas.mu / rho;
29843bff553SJames Wright 
2992b916ea7SJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
300d4e0f297SJames Wright     const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]};
3014f0244d1SJeremy L Thompson     CeedScalar       dXdx[3][3];
30284b557acSJames Wright     CeedScalar       h_node_sep[3];
30328a2a646SJames Wright 
30428a2a646SJames Wright     InvertMappingJacobian_3D(Q, i, J, dXdx, NULL);
30528a2a646SJames Wright     CalculateElementLengths(dXdx, stg_ctx->h_scale_factor, h_node_sep);
30628a2a646SJames Wright     h_node_sep[0] = dx * stg_ctx->h_scale_factor;
307d4e0f297SJames Wright 
308d4e0f297SJames Wright     InterpolateProfile(x_i[1], ubar, cij, &eps, &lt, stg_ctx);
309d4e0f297SJames Wright     if (stg_ctx->use_fluctuating_IC) {
31084b557acSJames Wright       CalcSpectrum(x_i[1], eps, lt, h_node_sep, nu, qn, stg_ctx);
31142454adaSJames Wright       StgShur14Calc(x_i, time, ubar, cij, qn, u, stg_ctx);
312d4e0f297SJames Wright     } else {
313d4e0f297SJames Wright       for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
314d4e0f297SJames Wright     }
31543bff553SJames Wright 
3169b103f75SJames Wright     CeedScalar Y[5] = {P0, u[0], u[1], u[2], theta0}, q[5];
3179b103f75SJames Wright     State      s    = StateFromY(gas, Y);
318*cde3d787SJames Wright     StateToQ(gas, s, q, context->state_var);
3199b103f75SJames Wright     for (CeedInt j = 0; j < 5; j++) {
3209b103f75SJames Wright       q0[j][i] = q[j];
32188243482SJames Wright     }
322b193fadcSJames Wright   }
32343bff553SJames Wright   return 0;
32443bff553SJames Wright }
32543bff553SJames Wright 
326493642f1SJames Wright /********************************************************************
327493642f1SJames Wright  * @brief QFunction to calculate the inflow boundary condition
328493642f1SJames Wright  *
329493642f1SJames Wright  * This will loop through quadrature points, calculate the wavemode amplitudes
330493642f1SJames Wright  * at each location, then calculate the actual velocity.
331493642f1SJames Wright  */
StgShur14Inflow(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)33242454adaSJames Wright CEED_QFUNCTION(StgShur14Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
3333d65b166SJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
334ade49511SJames Wright   const CeedScalar(*q_data_sur)    = in[2];
3353d65b166SJames Wright   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
336493642f1SJames Wright 
3373d65b166SJames Wright   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
338ade49511SJames Wright   CeedScalar(*jac_data_sur)  = out[1];
339493642f1SJames Wright 
34042454adaSJames Wright   const StgShur14Context      stg_ctx = (StgShur14Context)ctx;
341*cde3d787SJames Wright   const NewtonianIGProperties gas     = stg_ctx->newt_ctx.gas;
342493642f1SJames Wright   CeedScalar                  qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
343493642f1SJames Wright   const bool                  is_implicit = stg_ctx->is_implicit;
344493642f1SJames Wright   const bool                  mean_only   = stg_ctx->mean_only;
345493642f1SJames Wright   const bool                  prescribe_T = stg_ctx->prescribe_T;
346493642f1SJames Wright   const CeedScalar            dx          = stg_ctx->dx;
347*cde3d787SJames Wright   const CeedScalar            mu          = gas.mu;
348493642f1SJames Wright   const CeedScalar            time        = stg_ctx->time;
349493642f1SJames Wright   const CeedScalar            theta0      = stg_ctx->theta0;
350493642f1SJames Wright   const CeedScalar            P0          = stg_ctx->P0;
351*cde3d787SJames Wright   const CeedScalar            cv          = gas.cv;
352*cde3d787SJames Wright   const CeedScalar            Rd          = GasConstant(gas);
353*cde3d787SJames Wright   const CeedScalar            gamma       = HeatCapacityRatio(gas);
354493642f1SJames Wright 
3552b916ea7SJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
356493642f1SJames Wright     const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0);
357493642f1SJames Wright     const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
35878e8b7daSJames Wright     CeedScalar       wdetJb, dXdx[2][3], normal[3];
35978e8b7daSJames Wright     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, dXdx, normal);
360ade49511SJames Wright     wdetJb *= is_implicit ? -1. : 1.;
361493642f1SJames Wright 
36284b557acSJames Wright     CeedScalar h_node_sep[3];
36384b557acSJames Wright     h_node_sep[0] = dx;
36484b557acSJames Wright     for (CeedInt j = 1; j < 3; j++) h_node_sep[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]));
36584b557acSJames Wright     ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3);
366493642f1SJames Wright 
367493642f1SJames Wright     InterpolateProfile(X[1][i], ubar, cij, &eps, &lt, stg_ctx);
368493642f1SJames Wright     if (!mean_only) {
36984b557acSJames Wright       CalcSpectrum(X[1][i], eps, lt, h_node_sep, mu / rho, qn, stg_ctx);
37042454adaSJames Wright       StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx);
371493642f1SJames Wright     } else {
372493642f1SJames Wright       for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
373493642f1SJames Wright     }
374493642f1SJames Wright 
375a6e8f989SJames Wright     const CeedScalar E_kinetic = .5 * rho * Dot3(u, u);
376493642f1SJames Wright     CeedScalar       E_internal, P;
377493642f1SJames Wright     if (prescribe_T) {
378493642f1SJames Wright       // Temperature is being set weakly (theta0) and for constant cv this sets E_internal
379493642f1SJames Wright       E_internal = rho * cv * theta0;
380493642f1SJames Wright       // Find pressure using
381493642f1SJames Wright       P = rho * Rd * theta0;  // interior rho with exterior T
382493642f1SJames Wright     } else {
383493642f1SJames Wright       E_internal = q[4][i] - E_kinetic;  // uses prescribed rho and u, E from solution
384493642f1SJames Wright       P          = E_internal * (gamma - 1.);
385493642f1SJames Wright     }
386493642f1SJames Wright 
387493642f1SJames Wright     const CeedScalar E = E_internal + E_kinetic;
388493642f1SJames Wright 
389493642f1SJames Wright     // Velocity normal to the boundary
39078e8b7daSJames Wright     const CeedScalar u_normal = Dot3(normal, u);
391a6e8f989SJames Wright 
392493642f1SJames Wright     // The Physics
393493642f1SJames Wright     // Zero v so all future terms can safely sum into it
394493642f1SJames Wright     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
395493642f1SJames Wright 
396493642f1SJames Wright     // The Physics
397493642f1SJames Wright     // -- Density
398493642f1SJames Wright     v[0][i] -= wdetJb * rho * u_normal;
399493642f1SJames Wright 
400493642f1SJames Wright     // -- Momentum
40178e8b7daSJames Wright     for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + normal[j] * P);
402493642f1SJames Wright 
403493642f1SJames Wright     // -- Total Energy Density
404493642f1SJames Wright     v[4][i] -= wdetJb * u_normal * (E + P);
405a6e8f989SJames Wright 
406ade49511SJames Wright     const CeedScalar U[] = {rho, u[0], u[1], u[2], E}, kmstress[6] = {0.};
407ade49511SJames Wright     StoredValuesPack(Q, i, 0, 5, U, jac_data_sur);
408ade49511SJames Wright     StoredValuesPack(Q, i, 5, 6, kmstress, jac_data_sur);
409493642f1SJames Wright   }
410493642f1SJames Wright   return 0;
411493642f1SJames Wright }
412493642f1SJames Wright 
StgShur14Inflow_Jacobian(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)41342454adaSJames Wright CEED_QFUNCTION(StgShur14Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
4143d65b166SJames Wright   const CeedScalar(*dq)[CEED_Q_VLA]           = (const CeedScalar(*)[CEED_Q_VLA])in[0];
4153d65b166SJames Wright   const CeedScalar(*q_data_sur)[CEED_Q_VLA]   = (const CeedScalar(*)[CEED_Q_VLA])in[2];
4163d65b166SJames Wright   const CeedScalar(*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4];
417a6e8f989SJames Wright   CeedScalar(*v)[CEED_Q_VLA]                  = (CeedScalar(*)[CEED_Q_VLA])out[0];
4183d65b166SJames Wright 
41942454adaSJames Wright   const StgShur14Context      stg_ctx  = (StgShur14Context)ctx;
420*cde3d787SJames Wright   const NewtonianIGProperties gas      = stg_ctx->newt_ctx.gas;
421a6e8f989SJames Wright   const bool                  implicit = stg_ctx->is_implicit;
422*cde3d787SJames Wright   const CeedScalar            cv       = gas.cv;
423*cde3d787SJames Wright   const CeedScalar            Rd       = GasConstant(gas);
424*cde3d787SJames Wright   const CeedScalar            gamma    = HeatCapacityRatio(gas);
425a6e8f989SJames Wright 
426a6e8f989SJames Wright   const CeedScalar theta0      = stg_ctx->theta0;
427a6e8f989SJames Wright   const bool       prescribe_T = stg_ctx->prescribe_T;
428a6e8f989SJames Wright 
429b193fadcSJames 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     }
45678e8b7daSJames Wright     const CeedScalar normal[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
457a6e8f989SJames Wright 
45878e8b7daSJames Wright     const CeedScalar u_normal = Dot3(normal, velocity);
459a6e8f989SJames Wright 
460a6e8f989SJames Wright     v[0][i] = -wdetJb * drho * u_normal;
46178e8b7daSJames Wright     for (int j = 0; j < 3; j++) v[j + 1][i] = -wdetJb * (drho * u_normal * velocity[j] + normal[j] * dP);
462a6e8f989SJames Wright     v[4][i] = -wdetJb * u_normal * (dE + dP);
463b193fadcSJames 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  */
StgShur14InflowStrongQF(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)47342454adaSJames Wright CEED_QFUNCTION(StgShur14InflowStrongQF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
47421ba7ba4SJames Wright   const CeedScalar *dXdx_q              = 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;
481*cde3d787SJames Wright   const NewtonianIGProperties gas     = stg_ctx->newt_ctx.gas;
48270b0cb14SJames Wright   CeedScalar                  u[3], ubar[3], cij[6], eps, lt;
483b7190ff7SJames Wright   const bool                  mean_only = stg_ctx->mean_only;
484b7190ff7SJames Wright   const CeedScalar            time      = stg_ctx->time;
485b7190ff7SJames Wright   const CeedScalar            theta0    = stg_ctx->theta0;
486b7190ff7SJames Wright   const CeedScalar            P0        = stg_ctx->P0;
4879b103f75SJames Wright   const CeedScalar            rho       = P0 / (GasConstant(gas) * theta0);
488*cde3d787SJames Wright   const CeedScalar            nu        = gas.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]};
4928424c23aSJames Wright     CeedScalar       dXdx[3][3], h_node_sep[3];
49321ba7ba4SJames Wright     StoredValuesUnpack(Q, i, 0, 9, dXdx_q, (CeedScalar *)dXdx);
494b7190ff7SJames Wright 
4958424c23aSJames Wright     CalculateElementLengths(dXdx, stg_ctx->h_scale_factor, h_node_sep);
496b7190ff7SJames Wright     InterpolateProfile(coords[1][i], ubar, cij, &eps, &lt, stg_ctx);
497b7190ff7SJames Wright     if (!mean_only) {
49870b0cb14SJames Wright       if (1) {
49984b557acSJames Wright         StgShur14Calc_PrecompEktot(x, time, ubar, cij, inv_Ektotal[i], h_node_sep, x[1], eps, lt, nu, u, stg_ctx);
50070b0cb14SJames Wright       } else {  // Original way
50170b0cb14SJames Wright         CeedScalar qn[STG_NMODES_MAX];
50284b557acSJames Wright         CalcSpectrum(coords[1][i], eps, lt, h_node_sep, nu, qn, stg_ctx);
50342454adaSJames Wright         StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx);
50470b0cb14SJames Wright       }
505b7190ff7SJames Wright     } else {
506b7190ff7SJames Wright       for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
507b7190ff7SJames Wright     }
508b7190ff7SJames Wright 
5099b103f75SJames Wright     CeedScalar Y[5] = {P0, u[0], u[1], u[2], theta0}, q[5];
5109b103f75SJames Wright     State      s    = StateFromY(gas, Y);
511*cde3d787SJames Wright     StateToQ(gas, s, q, stg_ctx->newt_ctx.state_var);
512*cde3d787SJames Wright     switch (stg_ctx->newt_ctx.state_var) {
5133636f6a4SJames Wright       case STATEVAR_CONSERVATIVE:
5149b103f75SJames Wright         q[4] = 0.;  // Don't set energy
5153636f6a4SJames Wright         break;
5163636f6a4SJames Wright       case STATEVAR_PRIMITIVE:
5179b103f75SJames Wright         q[0] = 0;  // Don't set pressure
5183636f6a4SJames Wright         break;
5199b103f75SJames Wright       case STATEVAR_ENTROPY:
5209b103f75SJames Wright         q[0] = 0;  // Don't set V_density
5219b103f75SJames Wright         break;
5229b103f75SJames Wright     }
5239b103f75SJames Wright     for (CeedInt j = 0; j < 5; j++) {
5249b103f75SJames Wright       bcval[j][i] = scale[i] * q[j];
525b7190ff7SJames Wright     }
52688243482SJames Wright   }
527b7190ff7SJames Wright   return 0;
528b7190ff7SJames Wright }
529