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