xref: /honee/qfunctions/stg_shur14.h (revision a6e8f989aaa6b733abaa05a7879c1bbdb090d48c)
1493642f1SJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2493642f1SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3493642f1SJames Wright //
4493642f1SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5493642f1SJames Wright //
6493642f1SJames Wright // This file is part of CEED:  http://github.com/ceed
7493642f1SJames Wright 
8493642f1SJames Wright /// @file
9493642f1SJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm
10493642f1SJames Wright /// presented in Shur et al. 2014
11493642f1SJames Wright //
12493642f1SJames Wright /// SetupSTG_Rand reads in the input files and fills in STGShur14Context. Then
13493642f1SJames Wright /// STGShur14_CalcQF is run over quadrature points. Before the program exits,
14493642f1SJames Wright /// TearDownSTG is run to free the memory of the allocated arrays.
15493642f1SJames Wright 
16493642f1SJames Wright #ifndef stg_shur14_h
17493642f1SJames Wright #define stg_shur14_h
18493642f1SJames Wright 
19493642f1SJames Wright #include <math.h>
20493642f1SJames Wright #include <ceed.h>
21493642f1SJames Wright #include <stdlib.h>
22493642f1SJames Wright #include "stg_shur14_type.h"
23493642f1SJames Wright 
24493642f1SJames Wright #ifndef M_PI
25493642f1SJames Wright #define M_PI    3.14159265358979323846
26493642f1SJames Wright #endif
27493642f1SJames Wright 
28493642f1SJames Wright #define STG_NMODES_MAX 1024
29493642f1SJames Wright 
30493642f1SJames Wright CEED_QFUNCTION_HELPER CeedScalar Max(CeedScalar a, CeedScalar b) { return a < b ? b : a; }
31493642f1SJames Wright CEED_QFUNCTION_HELPER CeedScalar Min(CeedScalar a, CeedScalar b) { return a < b ? a : b; }
32493642f1SJames Wright 
33493642f1SJames Wright /*
34493642f1SJames Wright  * @brief Interpolate quantities from input profile to given location
35493642f1SJames Wright  *
36493642f1SJames Wright  * Assumed that prof_dw[i+1] > prof_dw[i] and prof_dw[0] = 0
37493642f1SJames Wright  * If dw > prof_dw[-1], then the interpolation takes the values at prof_dw[-1]
38493642f1SJames Wright  *
39493642f1SJames Wright  * @param[in]  dw      Distance to the nearest wall
40493642f1SJames Wright  * @param[out] ubar    Mean velocity at dw
41493642f1SJames Wright  * @param[out] cij     Cholesky decomposition at dw
42493642f1SJames Wright  * @param[out] eps     Turbulent dissipation at dw
43493642f1SJames Wright  * @param[out] lt      Turbulent length scale at dw
44493642f1SJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
45493642f1SJames Wright  */
46493642f1SJames Wright CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar dw,
47493642f1SJames Wright     CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt,
48493642f1SJames Wright     const STGShur14Context stg_ctx) {
49493642f1SJames Wright 
50493642f1SJames Wright   const CeedInt    nprofs    = stg_ctx->nprofs;
51493642f1SJames Wright   const CeedScalar *prof_dw  = &stg_ctx->data[stg_ctx->offsets.prof_dw];
52493642f1SJames Wright   const CeedScalar *prof_eps = &stg_ctx->data[stg_ctx->offsets.eps];
53493642f1SJames Wright   const CeedScalar *prof_lt  = &stg_ctx->data[stg_ctx->offsets.lt];
54493642f1SJames Wright   const CeedScalar *prof_ubar = &stg_ctx->data[stg_ctx->offsets.ubar];
55493642f1SJames Wright   const CeedScalar *prof_cij  = &stg_ctx->data[stg_ctx->offsets.cij];
56493642f1SJames Wright   CeedInt idx=-1;
57493642f1SJames Wright 
58493642f1SJames Wright   for(CeedInt i=0; i<nprofs; i++) {
59493642f1SJames Wright     if (dw < prof_dw[i]) {
60493642f1SJames Wright       idx = i;
61493642f1SJames Wright       break;
62493642f1SJames Wright     }
63493642f1SJames Wright   }
64493642f1SJames Wright 
65493642f1SJames Wright   if (idx > 0) { // y within the bounds of prof_dw
66493642f1SJames Wright     CeedScalar coeff = (dw - prof_dw[idx-1]) / (prof_dw[idx] - prof_dw[idx-1]);
67493642f1SJames Wright 
68493642f1SJames Wright     //*INDENT-OFF*
69493642f1SJames Wright     ubar[0] = prof_ubar[0*nprofs+idx-1] + coeff*( prof_ubar[0*nprofs+idx] - prof_ubar[0*nprofs+idx-1] );
70493642f1SJames Wright     ubar[1] = prof_ubar[1*nprofs+idx-1] + coeff*( prof_ubar[1*nprofs+idx] - prof_ubar[1*nprofs+idx-1] );
71493642f1SJames Wright     ubar[2] = prof_ubar[2*nprofs+idx-1] + coeff*( prof_ubar[2*nprofs+idx] - prof_ubar[2*nprofs+idx-1] );
72493642f1SJames Wright     cij[0]  = prof_cij[0*nprofs+idx-1]  + coeff*( prof_cij[0*nprofs+idx]  - prof_cij[0*nprofs+idx-1] );
73493642f1SJames Wright     cij[1]  = prof_cij[1*nprofs+idx-1]  + coeff*( prof_cij[1*nprofs+idx]  - prof_cij[1*nprofs+idx-1] );
74493642f1SJames Wright     cij[2]  = prof_cij[2*nprofs+idx-1]  + coeff*( prof_cij[2*nprofs+idx]  - prof_cij[2*nprofs+idx-1] );
75493642f1SJames Wright     cij[3]  = prof_cij[3*nprofs+idx-1]  + coeff*( prof_cij[3*nprofs+idx]  - prof_cij[3*nprofs+idx-1] );
76493642f1SJames Wright     cij[4]  = prof_cij[4*nprofs+idx-1]  + coeff*( prof_cij[4*nprofs+idx]  - prof_cij[4*nprofs+idx-1] );
77493642f1SJames Wright     cij[5]  = prof_cij[5*nprofs+idx-1]  + coeff*( prof_cij[5*nprofs+idx]  - prof_cij[5*nprofs+idx-1] );
78493642f1SJames Wright     *eps    = prof_eps[idx-1]           + coeff*( prof_eps[idx]           - prof_eps[idx-1] );
79493642f1SJames Wright     *lt     = prof_lt[idx-1]            + coeff*( prof_lt[idx]            - prof_lt[idx-1] );
80493642f1SJames Wright     //*INDENT-ON*
81493642f1SJames Wright   } else { // y outside bounds of prof_dw
82493642f1SJames Wright     ubar[0] = prof_ubar[1*nprofs-1];
83493642f1SJames Wright     ubar[1] = prof_ubar[2*nprofs-1];
84493642f1SJames Wright     ubar[2] = prof_ubar[3*nprofs-1];
85493642f1SJames Wright     cij[0]  = prof_cij[1*nprofs-1];
86493642f1SJames Wright     cij[1]  = prof_cij[2*nprofs-1];
87493642f1SJames Wright     cij[2]  = prof_cij[3*nprofs-1];
88493642f1SJames Wright     cij[3]  = prof_cij[4*nprofs-1];
89493642f1SJames Wright     cij[4]  = prof_cij[5*nprofs-1];
90493642f1SJames Wright     cij[5]  = prof_cij[6*nprofs-1];
91493642f1SJames Wright     *eps    = prof_eps[nprofs-1];
92493642f1SJames Wright     *lt     = prof_lt[nprofs-1];
93493642f1SJames Wright   }
94493642f1SJames Wright }
95493642f1SJames Wright 
96493642f1SJames Wright /*
97493642f1SJames Wright  * @brief Calculate spectrum coefficients for STG
98493642f1SJames Wright  *
99493642f1SJames Wright  * Calculates q_n at a given distance to the wall
100493642f1SJames Wright  *
101493642f1SJames Wright  * @param[in]  dw      Distance to the nearest wall
102493642f1SJames Wright  * @param[in]  eps     Turbulent dissipation w/rt dw
103493642f1SJames Wright  * @param[in]  lt      Turbulent length scale w/rt dw
104493642f1SJames Wright  * @param[in]  h       Element lengths in coordinate directions
105493642f1SJames Wright  * @param[in]  nu      Dynamic Viscosity;
106493642f1SJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
107493642f1SJames Wright  * @param[out] qn      Spectrum coefficients, [nmodes]
108493642f1SJames Wright  */
109493642f1SJames Wright void CEED_QFUNCTION_HELPER(CalcSpectrum)(const CeedScalar dw,
110493642f1SJames Wright     const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3],
111493642f1SJames Wright     const CeedScalar nu, CeedScalar qn[], const STGShur14Context stg_ctx) {
112493642f1SJames Wright 
113493642f1SJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
114493642f1SJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
115493642f1SJames Wright 
116493642f1SJames Wright   const CeedScalar hmax = Max( Max(h[0], h[1]), h[2]);
1172be05fbdSJames Wright   const CeedScalar ke   = dw==0 ? 1e16 : 2*M_PI/Min(2*dw, 3*lt);
118493642f1SJames Wright   const CeedScalar keta = 2*M_PI*pow(pow(nu,3.0)/eps, -0.25);
119493642f1SJames Wright   const CeedScalar kcut =
120493642f1SJames Wright     M_PI/ Min( Max(Max(h[1], h[2]), 0.3*hmax) + 0.1*dw, hmax );
121493642f1SJames Wright   CeedScalar fcut, feta, Ektot=0.0;
122493642f1SJames Wright 
123493642f1SJames Wright   for(CeedInt n=0; n<nmodes; n++) {
124493642f1SJames Wright     feta   = exp(-Square(12*kappa[n]/keta));
125493642f1SJames Wright     fcut   = exp( -pow(4*Max(kappa[n] - 0.9*kcut, 0)/kcut, 3.) );
126493642f1SJames Wright     qn[n]  = pow(kappa[n]/ke, 4.)
127493642f1SJames Wright              * pow(1 + 2.4*Square(kappa[n]/ke),-17./6)*feta*fcut;
128493642f1SJames Wright     qn[n] *= n==0 ? kappa[0] : kappa[n] - kappa[n-1];
129493642f1SJames Wright     Ektot += qn[n];
130493642f1SJames Wright   }
131493642f1SJames Wright 
1320a8dc919SJames Wright   if (Ektot == 0) return;
133493642f1SJames Wright   for(CeedInt n=0; n<nmodes; n++) qn[n] /= Ektot;
134493642f1SJames Wright }
135493642f1SJames Wright 
136493642f1SJames Wright /******************************************************
137493642f1SJames Wright  * @brief Calculate u(x,t) for STG inflow condition
138493642f1SJames Wright  *
139493642f1SJames Wright  * @param[in]  X       Location to evaluate u(X,t)
140493642f1SJames Wright  * @param[in]  t       Time to evaluate u(X,t)
141493642f1SJames Wright  * @param[in]  ubar    Mean velocity at X
142493642f1SJames Wright  * @param[in]  cij     Cholesky decomposition at X
143493642f1SJames Wright  * @param[in]  qn      Wavemode amplitudes at X, [nmodes]
144493642f1SJames Wright  * @param[out] u       Velocity at X and t
145493642f1SJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
146493642f1SJames Wright  */
147493642f1SJames Wright void CEED_QFUNCTION_HELPER(STGShur14_Calc)(const CeedScalar X[3],
148493642f1SJames Wright     const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6],
149493642f1SJames Wright     const CeedScalar qn[], CeedScalar u[3],
150493642f1SJames Wright     const STGShur14Context stg_ctx) {
151493642f1SJames Wright 
152493642f1SJames Wright   //*INDENT-OFF*
153493642f1SJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
154493642f1SJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
155493642f1SJames Wright   const CeedScalar *phi   = &stg_ctx->data[stg_ctx->offsets.phi];
156493642f1SJames Wright   const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma];
157493642f1SJames Wright   const CeedScalar *d     = &stg_ctx->data[stg_ctx->offsets.d];
158493642f1SJames Wright   //*INDENT-ON*
159493642f1SJames Wright   CeedScalar xdotd, vp[3] = {0.};
160493642f1SJames Wright   CeedScalar xhat[] = {0., X[1], X[2]};
161493642f1SJames Wright 
162493642f1SJames Wright   CeedPragmaSIMD
163493642f1SJames Wright   for(CeedInt n=0; n<nmodes; n++) {
164493642f1SJames Wright     xhat[0] = (X[0] - stg_ctx->u0*t)*Max(2*kappa[0]/kappa[n], 0.1);
165493642f1SJames Wright     xdotd = 0.;
166493642f1SJames Wright     for(CeedInt i=0; i<3; i++) xdotd += d[i*nmodes+n]*xhat[i];
167493642f1SJames Wright     const CeedScalar cos_kxdp = cos(kappa[n]*xdotd + phi[n]);
1680a8dc919SJames Wright     vp[0] += sqrt(qn[n])*sigma[0*nmodes+n] * cos_kxdp;
1690a8dc919SJames Wright     vp[1] += sqrt(qn[n])*sigma[1*nmodes+n] * cos_kxdp;
1700a8dc919SJames Wright     vp[2] += sqrt(qn[n])*sigma[2*nmodes+n] * cos_kxdp;
171493642f1SJames Wright   }
1720a8dc919SJames Wright   for(CeedInt i=0; i<3; i++) vp[i] *= 2*sqrt(1.5);
173493642f1SJames Wright 
174493642f1SJames Wright   u[0] = ubar[0] + cij[0]*vp[0];
175493642f1SJames Wright   u[1] = ubar[1] + cij[3]*vp[0] + cij[1]*vp[1];
176493642f1SJames Wright   u[2] = ubar[2] + cij[4]*vp[0] + cij[5]*vp[1] + cij[2]*vp[2];
177493642f1SJames Wright }
178493642f1SJames Wright 
17943bff553SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition
18043bff553SJames Wright CEED_QFUNCTION(ICsSTG)(void *ctx, CeedInt Q,
18143bff553SJames Wright                        const CeedScalar *const *in, CeedScalar *const *out) {
18243bff553SJames Wright   // Inputs
18343bff553SJames Wright   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
18443bff553SJames Wright 
18543bff553SJames Wright   // Outputs
18643bff553SJames Wright   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
18743bff553SJames Wright 
18843bff553SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
18943bff553SJames Wright   CeedScalar u[3], cij[6], eps, lt;
19043bff553SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
19143bff553SJames Wright   const CeedScalar P0     = stg_ctx->P0;
19243bff553SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
19343bff553SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
19443bff553SJames Wright   const CeedScalar Rd     = cp - cv;
19543bff553SJames Wright   const CeedScalar rho = P0 / (Rd * theta0);
19643bff553SJames Wright 
19743bff553SJames Wright   CeedPragmaSIMD
19843bff553SJames Wright   for(CeedInt i=0; i<Q; i++) {
19943bff553SJames Wright     InterpolateProfile(X[1][i], u, cij, &eps, &lt, stg_ctx);
20043bff553SJames Wright 
20143bff553SJames Wright     q0[0][i] = rho;
20243bff553SJames Wright     q0[1][i] = u[0] * rho;
20343bff553SJames Wright     q0[2][i] = u[1] * rho;
20443bff553SJames Wright     q0[3][i] = u[2] * rho;
20543bff553SJames Wright     q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0);
20643bff553SJames Wright   } // End of Quadrature Point Loop
20743bff553SJames Wright   return 0;
20843bff553SJames Wright }
20943bff553SJames Wright 
210493642f1SJames Wright /********************************************************************
211493642f1SJames Wright  * @brief QFunction to calculate the inflow boundary condition
212493642f1SJames Wright  *
213493642f1SJames Wright  * This will loop through quadrature points, calculate the wavemode amplitudes
214493642f1SJames Wright  * at each location, then calculate the actual velocity.
215493642f1SJames Wright  */
216493642f1SJames Wright CEED_QFUNCTION(STGShur14_Inflow)(void *ctx, CeedInt Q,
217493642f1SJames Wright                                  const CeedScalar *const *in,
218493642f1SJames Wright                                  CeedScalar *const *out) {
219493642f1SJames Wright 
220493642f1SJames Wright   //*INDENT-OFF*
221493642f1SJames Wright   const CeedScalar (*q)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
222dd64951cSJames Wright                    (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[2],
223dd64951cSJames Wright                    (*X)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[3];
224493642f1SJames Wright 
225*a6e8f989SJames Wright   CeedScalar(*v)[CEED_Q_VLA]            = (CeedScalar(*)[CEED_Q_VLA]) out[0],
226*a6e8f989SJames Wright             (*jac_data_sur)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[1];
227493642f1SJames Wright 
228493642f1SJames Wright   //*INDENT-ON*
229493642f1SJames Wright 
230493642f1SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
231493642f1SJames Wright   CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
232493642f1SJames Wright   const bool is_implicit  = stg_ctx->is_implicit;
233493642f1SJames Wright   const bool mean_only    = stg_ctx->mean_only;
234493642f1SJames Wright   const bool prescribe_T  = stg_ctx->prescribe_T;
235493642f1SJames Wright   const CeedScalar dx     = stg_ctx->dx;
236493642f1SJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
237493642f1SJames Wright   const CeedScalar time   = stg_ctx->time;
238493642f1SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
239493642f1SJames Wright   const CeedScalar P0     = stg_ctx->P0;
240493642f1SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
241493642f1SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
242493642f1SJames Wright   const CeedScalar Rd     = cp - cv;
243493642f1SJames Wright   const CeedScalar gamma  = cp/cv;
244493642f1SJames Wright 
245493642f1SJames Wright   CeedPragmaSIMD
246493642f1SJames Wright   for(CeedInt i=0; i<Q; i++) {
247493642f1SJames Wright     const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0);
248493642f1SJames Wright     const CeedScalar x[] = { X[0][i], X[1][i], X[2][i] };
249493642f1SJames Wright     const CeedScalar dXdx[2][3] = {
250493642f1SJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
251493642f1SJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
252493642f1SJames Wright     };
253493642f1SJames Wright 
254493642f1SJames Wright     CeedScalar h[3];
255493642f1SJames Wright     for (CeedInt j=0; j<3; j++)
256493642f1SJames Wright       h[j] = 2/sqrt(dXdx[0][j]*dXdx[0][j] + dXdx[1][j]*dXdx[1][j]);
257493642f1SJames Wright     h[0] = dx;
258493642f1SJames Wright 
259493642f1SJames Wright     InterpolateProfile(X[1][i], ubar, cij, &eps, &lt, stg_ctx);
260493642f1SJames Wright     if (!mean_only) {
261493642f1SJames Wright       CalcSpectrum(X[1][i], eps, lt, h, mu/rho, qn, stg_ctx);
262493642f1SJames Wright       STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
263493642f1SJames Wright     } else {
264493642f1SJames Wright       for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
265493642f1SJames Wright     }
266493642f1SJames Wright 
267*a6e8f989SJames Wright     const CeedScalar E_kinetic = .5 * rho * Dot3(u, u);
268493642f1SJames Wright     CeedScalar E_internal, P;
269493642f1SJames Wright     if (prescribe_T) {
270493642f1SJames Wright       // Temperature is being set weakly (theta0) and for constant cv this sets E_internal
271493642f1SJames Wright       E_internal = rho * cv * theta0;
272493642f1SJames Wright       // Find pressure using
273493642f1SJames Wright       P = rho * Rd * theta0; // interior rho with exterior T
274493642f1SJames Wright     } else {
275493642f1SJames Wright       E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution
276493642f1SJames Wright       P = E_internal * (gamma - 1.);
277493642f1SJames Wright     }
278493642f1SJames Wright 
279493642f1SJames Wright     const CeedScalar wdetJb  = (is_implicit ? -1. : 1.) * q_data_sur[0][i];
280493642f1SJames Wright     // ---- Normal vect
281493642f1SJames Wright     const CeedScalar norm[3] = {q_data_sur[1][i],
282493642f1SJames Wright                                 q_data_sur[2][i],
283493642f1SJames Wright                                 q_data_sur[3][i]
284493642f1SJames Wright                                };
285493642f1SJames Wright 
286493642f1SJames Wright     const CeedScalar E = E_internal + E_kinetic;
287493642f1SJames Wright 
288493642f1SJames Wright     // Velocity normal to the boundary
289*a6e8f989SJames Wright     const CeedScalar u_normal = Dot3(norm, u);
290*a6e8f989SJames Wright 
291493642f1SJames Wright     // The Physics
292493642f1SJames Wright     // Zero v so all future terms can safely sum into it
293493642f1SJames Wright     for (CeedInt j=0; j<5; j++) v[j][i] = 0.;
294493642f1SJames Wright 
295493642f1SJames Wright     // The Physics
296493642f1SJames Wright     // -- Density
297493642f1SJames Wright     v[0][i] -= wdetJb * rho * u_normal;
298493642f1SJames Wright 
299493642f1SJames Wright     // -- Momentum
300493642f1SJames Wright     for (CeedInt j=0; j<3; j++)
301493642f1SJames Wright       v[j+1][i] -= wdetJb *(rho * u_normal * u[j] +
302493642f1SJames Wright                             norm[j] * P);
303493642f1SJames Wright 
304493642f1SJames Wright     // -- Total Energy Density
305493642f1SJames Wright     v[4][i] -= wdetJb * u_normal * (E + P);
306*a6e8f989SJames Wright 
307*a6e8f989SJames Wright     jac_data_sur[0][i] = rho;
308*a6e8f989SJames Wright     jac_data_sur[1][i] = u[0];
309*a6e8f989SJames Wright     jac_data_sur[2][i] = u[1];
310*a6e8f989SJames Wright     jac_data_sur[3][i] = u[2];
311*a6e8f989SJames Wright     jac_data_sur[4][i] = E;
312*a6e8f989SJames Wright     for (int j=0; j<6; j++) jac_data_sur[5+j][i] = 0.;
313493642f1SJames Wright   }
314493642f1SJames Wright   return 0;
315493642f1SJames Wright }
316493642f1SJames Wright 
317*a6e8f989SJames Wright CEED_QFUNCTION(STGShur14_Inflow_Jacobian)(void *ctx, CeedInt Q,
318*a6e8f989SJames Wright     const CeedScalar *const *in,
319*a6e8f989SJames Wright     CeedScalar *const *out) {
320*a6e8f989SJames Wright   // *INDENT-OFF*
321*a6e8f989SJames Wright   // Inputs
322*a6e8f989SJames Wright   const CeedScalar (*dq)[CEED_Q_VLA]           = (const CeedScalar(*)[CEED_Q_VLA])in[0],
323*a6e8f989SJames Wright                    (*q_data_sur)[CEED_Q_VLA]   = (const CeedScalar(*)[CEED_Q_VLA])in[2],
324*a6e8f989SJames Wright                    (*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4];
325*a6e8f989SJames Wright   // Outputs
326*a6e8f989SJames Wright   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
327*a6e8f989SJames Wright   // *INDENT-ON*
328*a6e8f989SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
329*a6e8f989SJames Wright   const bool implicit     = stg_ctx->is_implicit;
330*a6e8f989SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
331*a6e8f989SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
332*a6e8f989SJames Wright   const CeedScalar Rd     = cp - cv;
333*a6e8f989SJames Wright   const CeedScalar gamma  = cp/cv;
334*a6e8f989SJames Wright 
335*a6e8f989SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
336*a6e8f989SJames Wright   const bool prescribe_T  = stg_ctx->prescribe_T;
337*a6e8f989SJames Wright 
338*a6e8f989SJames Wright   CeedPragmaSIMD
339*a6e8f989SJames Wright   // Quadrature Point Loop
340*a6e8f989SJames Wright   for (CeedInt i=0; i<Q; i++) {
341*a6e8f989SJames Wright     // Setup
342*a6e8f989SJames Wright     // Setup
343*a6e8f989SJames Wright     // -- Interp-to-Interp q_data
344*a6e8f989SJames Wright     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
345*a6e8f989SJames Wright     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
346*a6e8f989SJames Wright     // We can effect this by swapping the sign on this weight
347*a6e8f989SJames Wright     const CeedScalar wdetJb  = (implicit ? -1. : 1.) * q_data_sur[0][i];
348*a6e8f989SJames Wright 
349*a6e8f989SJames Wright     // Calculate inflow values
350*a6e8f989SJames Wright     CeedScalar velocity[3];
351*a6e8f989SJames Wright     for (CeedInt j=0; j<3; j++) velocity[j] = jac_data_sur[5+j][i];
352*a6e8f989SJames Wright 
353*a6e8f989SJames Wright     // enabling user to choose between weak T and weak rho inflow
354*a6e8f989SJames Wright     CeedScalar drho, dE, dP;
355*a6e8f989SJames Wright     if (prescribe_T) {
356*a6e8f989SJames Wright       // rho should be from the current solution
357*a6e8f989SJames Wright       drho = dq[0][i];
358*a6e8f989SJames Wright       CeedScalar dE_internal = drho * cv * theta0;
359*a6e8f989SJames Wright       CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity);
360*a6e8f989SJames Wright       dE = dE_internal + dE_kinetic;
361*a6e8f989SJames Wright       dP = drho * Rd * theta0; // interior rho with exterior T
362*a6e8f989SJames Wright     } else { // rho specified, E_internal from solution
363*a6e8f989SJames Wright       drho = 0;
364*a6e8f989SJames Wright       dE = dq[4][i];
365*a6e8f989SJames Wright       dP = dE * (gamma - 1.);
366*a6e8f989SJames Wright     }
367*a6e8f989SJames Wright     const CeedScalar norm[3] = {q_data_sur[1][i],
368*a6e8f989SJames Wright                                 q_data_sur[2][i],
369*a6e8f989SJames Wright                                 q_data_sur[3][i]
370*a6e8f989SJames Wright                                };
371*a6e8f989SJames Wright 
372*a6e8f989SJames Wright     const CeedScalar u_normal = Dot3(norm, velocity);
373*a6e8f989SJames Wright 
374*a6e8f989SJames Wright     v[0][i] = - wdetJb * drho * u_normal;
375*a6e8f989SJames Wright     for (int j=0; j<3; j++)
376*a6e8f989SJames Wright       v[j+1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP);
377*a6e8f989SJames Wright     v[4][i] = - wdetJb * u_normal * (dE + dP);
378*a6e8f989SJames Wright   } // End Quadrature Point Loop
379*a6e8f989SJames Wright   return 0;
380*a6e8f989SJames Wright }
381*a6e8f989SJames Wright 
382493642f1SJames Wright #endif // stg_shur14_h
383