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