xref: /libCEED/examples/fluids/qfunctions/stg_shur14.h (revision 831dbe9e1dadefb64122292c735c2cf2987406b2)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Implementation of the Synthetic Turbulence Generation (STG) algorithm
10 /// presented in Shur et al. 2014
11 //
12 /// SetupSTG_Rand reads in the input files and fills in STGShur14Context.
13 /// Then STGShur14_CalcQF is run over quadrature points.
14 /// Before the program exits, TearDownSTG is run to free the memory of the allocated arrays.
15 #include <ceed.h>
16 #include <math.h>
17 #include <stdlib.h>
18 
19 #include "newtonian_state.h"
20 #include "setupgeo_helpers.h"
21 #include "stg_shur14_type.h"
22 #include "utils.h"
23 
24 #define STG_NMODES_MAX 1024
25 
26 /*
27  * @brief Interpolate quantities from input profile to given location
28  *
29  * Assumed that prof_wd[i+1] > prof_wd[i] and prof_wd[0] = 0
30  * If wall_dist > prof_wd[-1], then the interpolation takes the values at prof_wd[-1]
31  *
32  * @param[in]  wall_dist Distance to the nearest wall
33  * @param[out] ubar      Mean velocity at wall_dist
34  * @param[out] cij       Cholesky decomposition at wall_dist
35  * @param[out] eps       Turbulent dissipation at wall_dist
36  * @param[out] lt        Turbulent length scale at wall_dist
37  * @param[in]  stg_ctx   STGShur14Context for the problem
38  */
39 CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar wall_dist, CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt,
40                                               const StgShur14Context stg_ctx) {
41   const CeedInt     nprofs    = stg_ctx->nprofs;
42   const CeedScalar *prof_wd   = &stg_ctx->data[stg_ctx->offsets.wall_dist];
43   const CeedScalar *prof_eps  = &stg_ctx->data[stg_ctx->offsets.eps];
44   const CeedScalar *prof_lt   = &stg_ctx->data[stg_ctx->offsets.lt];
45   const CeedScalar *prof_ubar = &stg_ctx->data[stg_ctx->offsets.ubar];
46   const CeedScalar *prof_cij  = &stg_ctx->data[stg_ctx->offsets.cij];
47   CeedInt           idx       = -1;
48 
49   for (CeedInt i = 0; i < nprofs; i++) {
50     if (wall_dist < prof_wd[i]) {
51       idx = i;
52       break;
53     }
54   }
55 
56   if (idx > 0) {  // y within the bounds of prof_wd
57     CeedScalar coeff = (wall_dist - prof_wd[idx - 1]) / (prof_wd[idx] - prof_wd[idx - 1]);
58 
59     ubar[0] = prof_ubar[0 * nprofs + idx - 1] + coeff * (prof_ubar[0 * nprofs + idx] - prof_ubar[0 * nprofs + idx - 1]);
60     ubar[1] = prof_ubar[1 * nprofs + idx - 1] + coeff * (prof_ubar[1 * nprofs + idx] - prof_ubar[1 * nprofs + idx - 1]);
61     ubar[2] = prof_ubar[2 * nprofs + idx - 1] + coeff * (prof_ubar[2 * nprofs + idx] - prof_ubar[2 * nprofs + idx - 1]);
62     cij[0]  = prof_cij[0 * nprofs + idx - 1] + coeff * (prof_cij[0 * nprofs + idx] - prof_cij[0 * nprofs + idx - 1]);
63     cij[1]  = prof_cij[1 * nprofs + idx - 1] + coeff * (prof_cij[1 * nprofs + idx] - prof_cij[1 * nprofs + idx - 1]);
64     cij[2]  = prof_cij[2 * nprofs + idx - 1] + coeff * (prof_cij[2 * nprofs + idx] - prof_cij[2 * nprofs + idx - 1]);
65     cij[3]  = prof_cij[3 * nprofs + idx - 1] + coeff * (prof_cij[3 * nprofs + idx] - prof_cij[3 * nprofs + idx - 1]);
66     cij[4]  = prof_cij[4 * nprofs + idx - 1] + coeff * (prof_cij[4 * nprofs + idx] - prof_cij[4 * nprofs + idx - 1]);
67     cij[5]  = prof_cij[5 * nprofs + idx - 1] + coeff * (prof_cij[5 * nprofs + idx] - prof_cij[5 * nprofs + idx - 1]);
68     *eps    = prof_eps[idx - 1] + coeff * (prof_eps[idx] - prof_eps[idx - 1]);
69     *lt     = prof_lt[idx - 1] + coeff * (prof_lt[idx] - prof_lt[idx - 1]);
70   } else {  // y outside bounds of prof_wd
71     ubar[0] = prof_ubar[1 * nprofs - 1];
72     ubar[1] = prof_ubar[2 * nprofs - 1];
73     ubar[2] = prof_ubar[3 * nprofs - 1];
74     cij[0]  = prof_cij[1 * nprofs - 1];
75     cij[1]  = prof_cij[2 * nprofs - 1];
76     cij[2]  = prof_cij[3 * nprofs - 1];
77     cij[3]  = prof_cij[4 * nprofs - 1];
78     cij[4]  = prof_cij[5 * nprofs - 1];
79     cij[5]  = prof_cij[6 * nprofs - 1];
80     *eps    = prof_eps[nprofs - 1];
81     *lt     = prof_lt[nprofs - 1];
82   }
83 }
84 
85 /*
86  * @brief Calculate spectrum coefficient, qn
87  *
88  * Calculates q_n at a given distance to the wall
89  *
90  * @param[in]  kappa     nth wavenumber
91  * @param[in]  dkappa    Difference between wavenumbers
92  * @param[in]  keta      Dissipation wavenumber
93  * @param[in]  kcut      Mesh-induced cutoff wavenumber
94  * @param[in]  ke        Energy-containing wavenumber
95  * @param[in]  Ektot_inv Inverse of total turbulent kinetic energy of spectrum
96  * @returns    qn        Spectrum coefficient
97  */
98 CEED_QFUNCTION_HELPER CeedScalar Calc_qn(const CeedScalar kappa, const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut,
99                                          const CeedScalar ke, const CeedScalar Ektot_inv) {
100   const CeedScalar feta_x_fcut = exp(-Square(12 * kappa / keta) - Cube(4 * Max(kappa - 0.9 * kcut, 0) / kcut));
101   return pow(kappa / ke, 4.) * pow(1 + 2.4 * Square(kappa / ke), -17. / 6) * feta_x_fcut * dkappa * Ektot_inv;
102 }
103 
104 // Calculate hmax, ke, keta, and kcut
105 CEED_QFUNCTION_HELPER void SpectrumConstants(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar hNodSep[3],
106                                              const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke, CeedScalar *keta, CeedScalar *kcut) {
107   *hmax = Max(Max(hNodSep[0], hNodSep[1]), hNodSep[2]);
108   *ke   = wall_dist == 0 ? 1e16 : 2 * M_PI / Min(2 * wall_dist, 3 * lt);
109   *keta = 2 * M_PI * pow(Cube(nu) / eps, -0.25);
110   *kcut = M_PI / Min(Max(Max(hNodSep[1], hNodSep[2]), 0.3 * (*hmax)) + 0.1 * wall_dist, *hmax);
111 }
112 
113 /*
114  * @brief Calculate spectrum coefficients for STG
115  *
116  * Calculates q_n at a given distance to the wall
117  *
118  * @param[in]  wall_dist  Distance to the nearest wall
119  * @param[in]  eps        Turbulent dissipation w/rt wall_dist
120  * @param[in]  lt         Turbulent length scale w/rt wall_dist
121  * @param[in]  h_node_sep Element lengths in coordinate directions
122  * @param[in]  nu         Dynamic Viscosity;
123  * @param[in]  stg_ctx    STGShur14Context for the problem
124  * @param[out] qn         Spectrum coefficients, [nmodes]
125  */
126 CEED_QFUNCTION_HELPER void CalcSpectrum(const CeedScalar wall_dist, const CeedScalar eps, const CeedScalar lt, const CeedScalar h_node_sep[3],
127                                         const CeedScalar nu, CeedScalar qn[], const StgShur14Context stg_ctx) {
128   const CeedInt     nmodes = stg_ctx->nmodes;
129   const CeedScalar *kappa  = &stg_ctx->data[stg_ctx->offsets.kappa];
130   CeedScalar        hmax, ke, keta, kcut, Ektot = 0.0;
131 
132   SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut);
133 
134   for (CeedInt n = 0; n < nmodes; n++) {
135     const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1];
136     qn[n]                   = Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0);
137     Ektot += qn[n];
138   }
139 
140   if (Ektot == 0) return;
141   for (CeedInt n = 0; n < nmodes; n++) qn[n] /= Ektot;
142 }
143 
144 /******************************************************
145  * @brief Calculate u(x,t) for STG inflow condition
146  *
147  * @param[in]  X       Location to evaluate u(X,t)
148  * @param[in]  t       Time to evaluate u(X,t)
149  * @param[in]  ubar    Mean velocity at X
150  * @param[in]  cij     Cholesky decomposition at X
151  * @param[in]  qn      Wavemode amplitudes at X, [nmodes]
152  * @param[out] u       Velocity at X and t
153  * @param[in]  stg_ctx STGShur14Context for the problem
154  */
155 CEED_QFUNCTION_HELPER void StgShur14Calc(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6],
156                                          const CeedScalar qn[], CeedScalar u[3], const StgShur14Context stg_ctx) {
157   const CeedInt     nmodes = stg_ctx->nmodes;
158   const CeedScalar *kappa  = &stg_ctx->data[stg_ctx->offsets.kappa];
159   const CeedScalar *phi    = &stg_ctx->data[stg_ctx->offsets.phi];
160   const CeedScalar *sigma  = &stg_ctx->data[stg_ctx->offsets.sigma];
161   const CeedScalar *d      = &stg_ctx->data[stg_ctx->offsets.d];
162   CeedScalar        xdotd, vp[3] = {0.};
163   CeedScalar        xhat[] = {0., X[1], X[2]};
164 
165   CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) {
166     xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1);
167     xdotd   = 0.;
168     for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i];
169     const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]);
170     vp[0] += sqrt(qn[n]) * sigma[0 * nmodes + n] * cos_kxdp;
171     vp[1] += sqrt(qn[n]) * sigma[1 * nmodes + n] * cos_kxdp;
172     vp[2] += sqrt(qn[n]) * sigma[2 * nmodes + n] * cos_kxdp;
173   }
174   for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5);
175 
176   u[0] = ubar[0] + cij[0] * vp[0];
177   u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1];
178   u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2];
179 }
180 
181 /******************************************************
182  * @brief Calculate u(x,t) for STG inflow condition
183  *
184  * @param[in]  X          Location to evaluate u(X,t)
185  * @param[in]  t          Time to evaluate u(X,t)
186  * @param[in]  ubar       Mean velocity at X
187  * @param[in]  cij        Cholesky decomposition at X
188  * @param[in]  Ektot      Total spectrum energy at this location
189  * @param[in]  h_node_sep Element size in 3 directions
190  * @param[in]  wall_dist  Distance to closest wall
191  * @param[in]  eps        Turbulent dissipation
192  * @param[in]  lt         Turbulent length scale
193  * @param[out] u          Velocity at X and t
194  * @param[in]  stg_ctx    STGShur14Context for the problem
195  */
196 CEED_QFUNCTION_HELPER void StgShur14Calc_PrecompEktot(const CeedScalar X[3], const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6],
197                                                       const CeedScalar Ektot, const CeedScalar h_node_sep[3], const CeedScalar wall_dist,
198                                                       const CeedScalar eps, const CeedScalar lt, const CeedScalar nu, CeedScalar u[3],
199                                                       const StgShur14Context stg_ctx) {
200   const CeedInt     nmodes = stg_ctx->nmodes;
201   const CeedScalar *kappa  = &stg_ctx->data[stg_ctx->offsets.kappa];
202   const CeedScalar *phi    = &stg_ctx->data[stg_ctx->offsets.phi];
203   const CeedScalar *sigma  = &stg_ctx->data[stg_ctx->offsets.sigma];
204   const CeedScalar *d      = &stg_ctx->data[stg_ctx->offsets.d];
205   CeedScalar        hmax, ke, keta, kcut;
206   SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut);
207   CeedScalar xdotd, vp[3] = {0.};
208   CeedScalar xhat[] = {0., X[1], X[2]};
209 
210   CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) {
211     xhat[0] = (X[0] - stg_ctx->u0 * t) * Max(2 * kappa[0] / kappa[n], 0.1);
212     xdotd   = 0.;
213     for (CeedInt i = 0; i < 3; i++) xdotd += d[i * nmodes + n] * xhat[i];
214     const CeedScalar cos_kxdp = cos(kappa[n] * xdotd + phi[n]);
215     const CeedScalar dkappa   = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1];
216     const CeedScalar qn       = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot);
217     vp[0] += sqrt(qn) * sigma[0 * nmodes + n] * cos_kxdp;
218     vp[1] += sqrt(qn) * sigma[1 * nmodes + n] * cos_kxdp;
219     vp[2] += sqrt(qn) * sigma[2 * nmodes + n] * cos_kxdp;
220   }
221   for (CeedInt i = 0; i < 3; i++) vp[i] *= 2 * sqrt(1.5);
222 
223   u[0] = ubar[0] + cij[0] * vp[0];
224   u[1] = ubar[1] + cij[3] * vp[0] + cij[1] * vp[1];
225   u[2] = ubar[2] + cij[4] * vp[0] + cij[5] * vp[1] + cij[2] * vp[2];
226 }
227 
228 // Create preprocessed input for the stg calculation
229 //
230 // stg_data[0] = 1 / Ektot (inverse of total spectrum energy)
231 CEED_QFUNCTION(StgShur14Preprocess)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
232   const CeedScalar(*dXdx_q)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0];
233   const CeedScalar(*x)[CEED_Q_VLA]         = (const CeedScalar(*)[CEED_Q_VLA])in[1];
234 
235   CeedScalar(*stg_data) = (CeedScalar(*))out[0];
236 
237   CeedScalar             ubar[3], cij[6], eps, lt;
238   const StgShur14Context stg_ctx = (StgShur14Context)ctx;
239   const CeedScalar       dx      = stg_ctx->dx;
240   const CeedScalar       mu      = stg_ctx->newtonian_ctx.mu;
241   const CeedScalar       theta0  = stg_ctx->theta0;
242   const CeedScalar       P0      = stg_ctx->P0;
243   const CeedScalar       Rd      = GasConstant(&stg_ctx->newtonian_ctx);
244   const CeedScalar       rho     = P0 / (Rd * theta0);
245   const CeedScalar       nu      = mu / rho;
246 
247   const CeedInt     nmodes = stg_ctx->nmodes;
248   const CeedScalar *kappa  = &stg_ctx->data[stg_ctx->offsets.kappa];
249   CeedScalar        hmax, ke, keta, kcut;
250 
251   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
252     const CeedScalar wall_dist  = x[1][i];
253     const CeedScalar dXdx[2][3] = {
254         {dXdx_q[0][0][i], dXdx_q[0][1][i], dXdx_q[0][2][i]},
255         {dXdx_q[1][0][i], dXdx_q[1][1][i], dXdx_q[1][2][i]},
256     };
257 
258     CeedScalar h_node_sep[3];
259     h_node_sep[0] = dx;
260     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]);
261     ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3);
262 
263     InterpolateProfile(wall_dist, ubar, cij, &eps, &lt, stg_ctx);
264     SpectrumConstants(wall_dist, eps, lt, h_node_sep, nu, &hmax, &ke, &keta, &kcut);
265 
266     // Calculate total TKE per spectrum
267     CeedScalar Ek_tot = 0;
268     CeedPragmaSIMD for (CeedInt n = 0; n < nmodes; n++) {
269       const CeedScalar dkappa = n == 0 ? kappa[0] : kappa[n] - kappa[n - 1];
270       Ek_tot += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0);
271     }
272     // avoid underflowed and poorly defined spectrum coefficients
273     stg_data[i] = Ek_tot != 0 ? 1 / Ek_tot : 0;
274   }
275   return 0;
276 }
277 
278 // Extrude the STGInflow profile through out the domain for an initial condition
279 CEED_QFUNCTION(ICsStg)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
280   const CeedScalar(*x)[CEED_Q_VLA]    = (const CeedScalar(*)[CEED_Q_VLA])in[0];
281   const CeedScalar(*J)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[1];
282   CeedScalar(*q0)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
283 
284   const StgShur14Context stg_ctx = (StgShur14Context)ctx;
285   CeedScalar             qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
286   const CeedScalar       dx     = stg_ctx->dx;
287   const CeedScalar       time   = stg_ctx->time;
288   const CeedScalar       theta0 = stg_ctx->theta0;
289   const CeedScalar       P0     = stg_ctx->P0;
290   const CeedScalar       cv     = stg_ctx->newtonian_ctx.cv;
291   const CeedScalar       rho    = P0 / (GasConstant(&stg_ctx->newtonian_ctx) * theta0);
292   const CeedScalar       nu     = stg_ctx->newtonian_ctx.mu / rho;
293 
294   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
295     const CeedScalar x_i[3] = {x[0][i], x[1][i], x[2][i]};
296     CeedScalar       dXdx[3][3];
297     InvertMappingJacobian_3D(Q, i, J, dXdx, NULL);
298     CeedScalar h_node_sep[3];
299     h_node_sep[0] = dx;
300     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]));
301     ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3);
302 
303     InterpolateProfile(x_i[1], ubar, cij, &eps, &lt, stg_ctx);
304     if (stg_ctx->use_fluctuating_IC) {
305       CalcSpectrum(x_i[1], eps, lt, h_node_sep, nu, qn, stg_ctx);
306       StgShur14Calc(x_i, time, ubar, cij, qn, u, stg_ctx);
307     } else {
308       for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
309     }
310 
311     switch (stg_ctx->newtonian_ctx.state_var) {
312       case STATEVAR_CONSERVATIVE:
313         q0[0][i] = rho;
314         q0[1][i] = u[0] * rho;
315         q0[2][i] = u[1] * rho;
316         q0[3][i] = u[2] * rho;
317         q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0);
318         break;
319 
320       case STATEVAR_PRIMITIVE:
321         q0[0][i] = P0;
322         q0[1][i] = u[0];
323         q0[2][i] = u[1];
324         q0[3][i] = u[2];
325         q0[4][i] = theta0;
326         break;
327     }
328   }
329   return 0;
330 }
331 
332 /********************************************************************
333  * @brief QFunction to calculate the inflow boundary condition
334  *
335  * This will loop through quadrature points, calculate the wavemode amplitudes
336  * at each location, then calculate the actual velocity.
337  */
338 CEED_QFUNCTION(StgShur14Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
339   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
340   const CeedScalar(*q_data_sur)    = in[2];
341   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
342 
343   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
344   CeedScalar(*jac_data_sur)  = out[1];
345 
346   const StgShur14Context stg_ctx = (StgShur14Context)ctx;
347   CeedScalar             qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
348   const bool             is_implicit = stg_ctx->is_implicit;
349   const bool             mean_only   = stg_ctx->mean_only;
350   const bool             prescribe_T = stg_ctx->prescribe_T;
351   const CeedScalar       dx          = stg_ctx->dx;
352   const CeedScalar       mu          = stg_ctx->newtonian_ctx.mu;
353   const CeedScalar       time        = stg_ctx->time;
354   const CeedScalar       theta0      = stg_ctx->theta0;
355   const CeedScalar       P0          = stg_ctx->P0;
356   const CeedScalar       cv          = stg_ctx->newtonian_ctx.cv;
357   const CeedScalar       Rd          = GasConstant(&stg_ctx->newtonian_ctx);
358   const CeedScalar       gamma       = HeatCapacityRatio(&stg_ctx->newtonian_ctx);
359 
360   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
361     const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0);
362     const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
363     CeedScalar       wdetJb, dXdx[2][3], norm[3];
364     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, dXdx, norm);
365     wdetJb *= is_implicit ? -1. : 1.;
366 
367     CeedScalar h_node_sep[3];
368     h_node_sep[0] = dx;
369     for (CeedInt j = 1; j < 3; j++) h_node_sep[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]));
370     ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3);
371 
372     InterpolateProfile(X[1][i], ubar, cij, &eps, &lt, stg_ctx);
373     if (!mean_only) {
374       CalcSpectrum(X[1][i], eps, lt, h_node_sep, mu / rho, qn, stg_ctx);
375       StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx);
376     } else {
377       for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
378     }
379 
380     const CeedScalar E_kinetic = .5 * rho * Dot3(u, u);
381     CeedScalar       E_internal, P;
382     if (prescribe_T) {
383       // Temperature is being set weakly (theta0) and for constant cv this sets E_internal
384       E_internal = rho * cv * theta0;
385       // Find pressure using
386       P = rho * Rd * theta0;  // interior rho with exterior T
387     } else {
388       E_internal = q[4][i] - E_kinetic;  // uses prescribed rho and u, E from solution
389       P          = E_internal * (gamma - 1.);
390     }
391 
392     const CeedScalar E = E_internal + E_kinetic;
393 
394     // Velocity normal to the boundary
395     const CeedScalar u_normal = Dot3(norm, u);
396 
397     // The Physics
398     // Zero v so all future terms can safely sum into it
399     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
400 
401     // The Physics
402     // -- Density
403     v[0][i] -= wdetJb * rho * u_normal;
404 
405     // -- Momentum
406     for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P);
407 
408     // -- Total Energy Density
409     v[4][i] -= wdetJb * u_normal * (E + P);
410 
411     const CeedScalar U[] = {rho, u[0], u[1], u[2], E}, kmstress[6] = {0.};
412     StoredValuesPack(Q, i, 0, 5, U, jac_data_sur);
413     StoredValuesPack(Q, i, 5, 6, kmstress, jac_data_sur);
414   }
415   return 0;
416 }
417 
418 CEED_QFUNCTION(StgShur14Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
419   const CeedScalar(*dq)[CEED_Q_VLA]           = (const CeedScalar(*)[CEED_Q_VLA])in[0];
420   const CeedScalar(*q_data_sur)[CEED_Q_VLA]   = (const CeedScalar(*)[CEED_Q_VLA])in[2];
421   const CeedScalar(*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4];
422   CeedScalar(*v)[CEED_Q_VLA]                  = (CeedScalar(*)[CEED_Q_VLA])out[0];
423 
424   const StgShur14Context stg_ctx  = (StgShur14Context)ctx;
425   const bool             implicit = stg_ctx->is_implicit;
426   const CeedScalar       cv       = stg_ctx->newtonian_ctx.cv;
427   const CeedScalar       Rd       = GasConstant(&stg_ctx->newtonian_ctx);
428   const CeedScalar       gamma    = HeatCapacityRatio(&stg_ctx->newtonian_ctx);
429 
430   const CeedScalar theta0      = stg_ctx->theta0;
431   const bool       prescribe_T = stg_ctx->prescribe_T;
432 
433   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
434     // Setup
435     // -- Interp-to-Interp q_data
436     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
437     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
438     // We can effect this by swapping the sign on this weight
439     const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i];
440 
441     // Calculate inflow values
442     CeedScalar velocity[3];
443     for (CeedInt j = 0; j < 3; j++) velocity[j] = jac_data_sur[5 + j][i];
444     // TODO This is almost certainly a bug. Velocity isn't stored here, only 0s.
445 
446     // enabling user to choose between weak T and weak rho inflow
447     CeedScalar drho, dE, dP;
448     if (prescribe_T) {
449       // rho should be from the current solution
450       drho                   = dq[0][i];
451       CeedScalar dE_internal = drho * cv * theta0;
452       CeedScalar dE_kinetic  = .5 * drho * Dot3(velocity, velocity);
453       dE                     = dE_internal + dE_kinetic;
454       dP                     = drho * Rd * theta0;  // interior rho with exterior T
455     } else {                                        // rho specified, E_internal from solution
456       drho = 0;
457       dE   = dq[4][i];
458       dP   = dE * (gamma - 1.);
459     }
460     const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
461 
462     const CeedScalar u_normal = Dot3(norm, velocity);
463 
464     v[0][i] = -wdetJb * drho * u_normal;
465     for (int j = 0; j < 3; j++) v[j + 1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP);
466     v[4][i] = -wdetJb * u_normal * (dE + dP);
467   }
468   return 0;
469 }
470 
471 /********************************************************************
472  * @brief QFunction to calculate the strongly enforce inflow BC
473  *
474  * This QF is for the strong application of STG via libCEED (rather than
475  * through the native PETSc `DMAddBoundary` -> `bcFunc` method.
476  */
477 CEED_QFUNCTION(StgShur14InflowStrongQF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
478   const CeedScalar(*dXdx_q)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0];
479   const CeedScalar(*coords)[CEED_Q_VLA]    = (const CeedScalar(*)[CEED_Q_VLA])in[1];
480   const CeedScalar(*scale)                 = (const CeedScalar(*))in[2];
481   const CeedScalar(*inv_Ektotal)           = (const CeedScalar(*))in[3];
482   CeedScalar(*bcval)[CEED_Q_VLA]           = (CeedScalar(*)[CEED_Q_VLA])out[0];
483 
484   const StgShur14Context stg_ctx = (StgShur14Context)ctx;
485   CeedScalar             u[3], ubar[3], cij[6], eps, lt;
486   const bool             mean_only = stg_ctx->mean_only;
487   const CeedScalar       dx        = stg_ctx->dx;
488   const CeedScalar       time      = stg_ctx->time;
489   const CeedScalar       theta0    = stg_ctx->theta0;
490   const CeedScalar       P0        = stg_ctx->P0;
491   const CeedScalar       rho       = P0 / (GasConstant(&stg_ctx->newtonian_ctx) * theta0);
492   const CeedScalar       nu        = stg_ctx->newtonian_ctx.mu / rho;
493 
494   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
495     const CeedScalar x[]        = {coords[0][i], coords[1][i], coords[2][i]};
496     const CeedScalar dXdx[2][3] = {
497         {dXdx_q[0][0][i], dXdx_q[0][1][i], dXdx_q[0][2][i]},
498         {dXdx_q[1][0][i], dXdx_q[1][1][i], dXdx_q[1][2][i]},
499     };
500 
501     CeedScalar h_node_sep[3];
502     h_node_sep[0] = dx;
503     for (CeedInt j = 1; j < 3; j++) h_node_sep[j] = 2 / sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]));
504     ScaleN(h_node_sep, stg_ctx->h_scale_factor, 3);
505 
506     InterpolateProfile(coords[1][i], ubar, cij, &eps, &lt, stg_ctx);
507     if (!mean_only) {
508       if (1) {
509         StgShur14Calc_PrecompEktot(x, time, ubar, cij, inv_Ektotal[i], h_node_sep, x[1], eps, lt, nu, u, stg_ctx);
510       } else {  // Original way
511         CeedScalar qn[STG_NMODES_MAX];
512         CalcSpectrum(coords[1][i], eps, lt, h_node_sep, nu, qn, stg_ctx);
513         StgShur14Calc(x, time, ubar, cij, qn, u, stg_ctx);
514       }
515     } else {
516       for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
517     }
518 
519     switch (stg_ctx->newtonian_ctx.state_var) {
520       case STATEVAR_CONSERVATIVE:
521         bcval[0][i] = scale[i] * rho;
522         bcval[1][i] = scale[i] * rho * u[0];
523         bcval[2][i] = scale[i] * rho * u[1];
524         bcval[3][i] = scale[i] * rho * u[2];
525         bcval[4][i] = 0.;
526         break;
527 
528       case STATEVAR_PRIMITIVE:
529         bcval[0][i] = 0;
530         bcval[1][i] = scale[i] * u[0];
531         bcval[2][i] = scale[i] * u[1];
532         bcval[3][i] = scale[i] * u[2];
533         bcval[4][i] = scale[i] * theta0;
534         break;
535     }
536   }
537   return 0;
538 }
539