xref: /libCEED/examples/fluids/qfunctions/stg_shur14.h (revision e159aeac15e8e644cec5ea034206454e5074ed03)
1ba6664aeSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2ba6664aeSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3ba6664aeSJames Wright //
4ba6664aeSJames Wright // SPDX-License-Identifier: BSD-2-Clause
5ba6664aeSJames Wright //
6ba6664aeSJames Wright // This file is part of CEED:  http://github.com/ceed
7ba6664aeSJames Wright 
8ba6664aeSJames Wright /// @file
9ba6664aeSJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm
10ba6664aeSJames Wright /// presented in Shur et al. 2014
11ba6664aeSJames Wright //
12ba6664aeSJames Wright /// SetupSTG_Rand reads in the input files and fills in STGShur14Context. Then
13ba6664aeSJames Wright /// STGShur14_CalcQF is run over quadrature points. Before the program exits,
14ba6664aeSJames Wright /// TearDownSTG is run to free the memory of the allocated arrays.
15ba6664aeSJames Wright 
16ba6664aeSJames Wright #ifndef stg_shur14_h
17ba6664aeSJames Wright #define stg_shur14_h
18ba6664aeSJames Wright 
19ba6664aeSJames Wright #include <math.h>
20ba6664aeSJames Wright #include <ceed.h>
21ba6664aeSJames Wright #include <stdlib.h>
22ba6664aeSJames Wright #include "stg_shur14_type.h"
2313fa47b2SJames Wright #include "utils.h"
24ba6664aeSJames Wright 
25ba6664aeSJames Wright #define STG_NMODES_MAX 1024
26ba6664aeSJames Wright 
27ba6664aeSJames Wright /*
28ba6664aeSJames Wright  * @brief Interpolate quantities from input profile to given location
29ba6664aeSJames Wright  *
30ba6664aeSJames Wright  * Assumed that prof_dw[i+1] > prof_dw[i] and prof_dw[0] = 0
31ba6664aeSJames Wright  * If dw > prof_dw[-1], then the interpolation takes the values at prof_dw[-1]
32ba6664aeSJames Wright  *
33ba6664aeSJames Wright  * @param[in]  dw      Distance to the nearest wall
34ba6664aeSJames Wright  * @param[out] ubar    Mean velocity at dw
35ba6664aeSJames Wright  * @param[out] cij     Cholesky decomposition at dw
36ba6664aeSJames Wright  * @param[out] eps     Turbulent dissipation at dw
37ba6664aeSJames Wright  * @param[out] lt      Turbulent length scale at dw
38ba6664aeSJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
39ba6664aeSJames Wright  */
40ba6664aeSJames Wright CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar dw,
41ba6664aeSJames Wright     CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt,
42ba6664aeSJames Wright     const STGShur14Context stg_ctx) {
43ba6664aeSJames Wright 
44ba6664aeSJames Wright   const CeedInt    nprofs    = stg_ctx->nprofs;
45ba6664aeSJames Wright   const CeedScalar *prof_dw  = &stg_ctx->data[stg_ctx->offsets.prof_dw];
46ba6664aeSJames Wright   const CeedScalar *prof_eps = &stg_ctx->data[stg_ctx->offsets.eps];
47ba6664aeSJames Wright   const CeedScalar *prof_lt  = &stg_ctx->data[stg_ctx->offsets.lt];
48ba6664aeSJames Wright   const CeedScalar *prof_ubar = &stg_ctx->data[stg_ctx->offsets.ubar];
49ba6664aeSJames Wright   const CeedScalar *prof_cij  = &stg_ctx->data[stg_ctx->offsets.cij];
50ba6664aeSJames Wright   CeedInt idx=-1;
51ba6664aeSJames Wright 
52ba6664aeSJames Wright   for(CeedInt i=0; i<nprofs; i++) {
53ba6664aeSJames Wright     if (dw < prof_dw[i]) {
54ba6664aeSJames Wright       idx = i;
55ba6664aeSJames Wright       break;
56ba6664aeSJames Wright     }
57ba6664aeSJames Wright   }
58ba6664aeSJames Wright 
59ba6664aeSJames Wright   if (idx > 0) { // y within the bounds of prof_dw
60ba6664aeSJames Wright     CeedScalar coeff = (dw - prof_dw[idx-1]) / (prof_dw[idx] - prof_dw[idx-1]);
61ba6664aeSJames Wright 
62ba6664aeSJames Wright     //*INDENT-OFF*
63ba6664aeSJames Wright     ubar[0] = prof_ubar[0*nprofs+idx-1] + coeff*( prof_ubar[0*nprofs+idx] - prof_ubar[0*nprofs+idx-1] );
64ba6664aeSJames Wright     ubar[1] = prof_ubar[1*nprofs+idx-1] + coeff*( prof_ubar[1*nprofs+idx] - prof_ubar[1*nprofs+idx-1] );
65ba6664aeSJames Wright     ubar[2] = prof_ubar[2*nprofs+idx-1] + coeff*( prof_ubar[2*nprofs+idx] - prof_ubar[2*nprofs+idx-1] );
66ba6664aeSJames Wright     cij[0]  = prof_cij[0*nprofs+idx-1]  + coeff*( prof_cij[0*nprofs+idx]  - prof_cij[0*nprofs+idx-1] );
67ba6664aeSJames Wright     cij[1]  = prof_cij[1*nprofs+idx-1]  + coeff*( prof_cij[1*nprofs+idx]  - prof_cij[1*nprofs+idx-1] );
68ba6664aeSJames Wright     cij[2]  = prof_cij[2*nprofs+idx-1]  + coeff*( prof_cij[2*nprofs+idx]  - prof_cij[2*nprofs+idx-1] );
69ba6664aeSJames Wright     cij[3]  = prof_cij[3*nprofs+idx-1]  + coeff*( prof_cij[3*nprofs+idx]  - prof_cij[3*nprofs+idx-1] );
70ba6664aeSJames Wright     cij[4]  = prof_cij[4*nprofs+idx-1]  + coeff*( prof_cij[4*nprofs+idx]  - prof_cij[4*nprofs+idx-1] );
71ba6664aeSJames Wright     cij[5]  = prof_cij[5*nprofs+idx-1]  + coeff*( prof_cij[5*nprofs+idx]  - prof_cij[5*nprofs+idx-1] );
72ba6664aeSJames Wright     *eps    = prof_eps[idx-1]           + coeff*( prof_eps[idx]           - prof_eps[idx-1] );
73ba6664aeSJames Wright     *lt     = prof_lt[idx-1]            + coeff*( prof_lt[idx]            - prof_lt[idx-1] );
74ba6664aeSJames Wright     //*INDENT-ON*
75ba6664aeSJames Wright   } else { // y outside bounds of prof_dw
76ba6664aeSJames Wright     ubar[0] = prof_ubar[1*nprofs-1];
77ba6664aeSJames Wright     ubar[1] = prof_ubar[2*nprofs-1];
78ba6664aeSJames Wright     ubar[2] = prof_ubar[3*nprofs-1];
79ba6664aeSJames Wright     cij[0]  = prof_cij[1*nprofs-1];
80ba6664aeSJames Wright     cij[1]  = prof_cij[2*nprofs-1];
81ba6664aeSJames Wright     cij[2]  = prof_cij[3*nprofs-1];
82ba6664aeSJames Wright     cij[3]  = prof_cij[4*nprofs-1];
83ba6664aeSJames Wright     cij[4]  = prof_cij[5*nprofs-1];
84ba6664aeSJames Wright     cij[5]  = prof_cij[6*nprofs-1];
85ba6664aeSJames Wright     *eps    = prof_eps[nprofs-1];
86ba6664aeSJames Wright     *lt     = prof_lt[nprofs-1];
87ba6664aeSJames Wright   }
88ba6664aeSJames Wright }
89ba6664aeSJames Wright 
90ba6664aeSJames Wright /*
91*e159aeacSJames Wright  * @brief Calculate spectrum coefficient, qn
92*e159aeacSJames Wright  *
93*e159aeacSJames Wright  * Calculates q_n at a given distance to the wall
94*e159aeacSJames Wright  *
95*e159aeacSJames Wright  * @param[in]  kappa  nth wavenumber
96*e159aeacSJames Wright  * @param[in]  dkappa Difference between wavenumbers
97*e159aeacSJames Wright  * @param[in]  keta   Dissipation wavenumber
98*e159aeacSJames Wright  * @param[in]  kcut   Mesh-induced cutoff wavenumber
99*e159aeacSJames Wright  * @param[in]  ke     Energy-containing wavenumber
100*e159aeacSJames Wright  * @param[in]  Ektot  Total turbulent kinetic energy of spectrum
101*e159aeacSJames Wright  * @returns    qn     Spectrum coefficient
102*e159aeacSJames Wright  */
103*e159aeacSJames Wright CeedScalar CEED_QFUNCTION_HELPER(Calc_qn)(const CeedScalar kappa,
104*e159aeacSJames Wright     const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut,
105*e159aeacSJames Wright     const CeedScalar ke, const CeedScalar Ektot) {
106*e159aeacSJames Wright   const CeedScalar feta_x_fcut   = exp(-Square(12*kappa/keta)
107*e159aeacSJames Wright                                        -Cube(4*Max(kappa - 0.9*kcut, 0)/kcut) );
108*e159aeacSJames Wright   return pow(kappa/ke, 4.) * pow(1 + 2.4*Square(kappa/ke),-17./6)
109*e159aeacSJames Wright          *feta_x_fcut*dkappa/Ektot;
110*e159aeacSJames Wright }
111*e159aeacSJames Wright 
112*e159aeacSJames Wright // Calculate hmax, ke, keta, and kcut
113*e159aeacSJames Wright void CEED_QFUNCTION_HELPER(SpectrumConstants)(const CeedScalar dw,
114*e159aeacSJames Wright     const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3],
115*e159aeacSJames Wright     const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke,
116*e159aeacSJames Wright     CeedScalar *keta, CeedScalar *kcut) {
117*e159aeacSJames Wright   *hmax = Max( Max(h[0], h[1]), h[2]);
118*e159aeacSJames Wright   *ke   = dw==0 ? 1e16 : 2*M_PI/Min(2*dw, 3*lt);
119*e159aeacSJames Wright   *keta = 2*M_PI*pow(Cube(nu)/eps, -0.25);
120*e159aeacSJames Wright   *kcut = M_PI/ Min( Max(Max(h[1], h[2]), 0.3*(*hmax)) + 0.1*dw, *hmax );
121*e159aeacSJames Wright }
122*e159aeacSJames Wright 
123*e159aeacSJames Wright /*
124ba6664aeSJames Wright  * @brief Calculate spectrum coefficients for STG
125ba6664aeSJames Wright  *
126ba6664aeSJames Wright  * Calculates q_n at a given distance to the wall
127ba6664aeSJames Wright  *
128ba6664aeSJames Wright  * @param[in]  dw      Distance to the nearest wall
129ba6664aeSJames Wright  * @param[in]  eps     Turbulent dissipation w/rt dw
130ba6664aeSJames Wright  * @param[in]  lt      Turbulent length scale w/rt dw
131ba6664aeSJames Wright  * @param[in]  h       Element lengths in coordinate directions
132ba6664aeSJames Wright  * @param[in]  nu      Dynamic Viscosity;
133ba6664aeSJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
134ba6664aeSJames Wright  * @param[out] qn      Spectrum coefficients, [nmodes]
135ba6664aeSJames Wright  */
136ba6664aeSJames Wright void CEED_QFUNCTION_HELPER(CalcSpectrum)(const CeedScalar dw,
137ba6664aeSJames Wright     const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3],
138ba6664aeSJames Wright     const CeedScalar nu, CeedScalar qn[], const STGShur14Context stg_ctx) {
139ba6664aeSJames Wright 
140ba6664aeSJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
141ba6664aeSJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
142*e159aeacSJames Wright   CeedScalar hmax, ke, keta, kcut, Ektot=0.0;
143*e159aeacSJames Wright   SpectrumConstants(dw, eps, lt, h, nu, &hmax, &ke, &keta, &kcut);
144ba6664aeSJames Wright 
145ba6664aeSJames Wright   for(CeedInt n=0; n<nmodes; n++) {
146*e159aeacSJames Wright     const CeedScalar dkappa = n==0 ? kappa[0] : kappa[n] - kappa[n-1];
147*e159aeacSJames Wright     qn[n] = Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0);
148ba6664aeSJames Wright     Ektot += qn[n];
149ba6664aeSJames Wright   }
150ba6664aeSJames Wright 
151961c9c98SJames Wright   if (Ektot == 0) return;
152ba6664aeSJames Wright   for(CeedInt n=0; n<nmodes; n++) qn[n] /= Ektot;
153ba6664aeSJames Wright }
154ba6664aeSJames Wright 
155ba6664aeSJames Wright /******************************************************
156ba6664aeSJames Wright  * @brief Calculate u(x,t) for STG inflow condition
157ba6664aeSJames Wright  *
158ba6664aeSJames Wright  * @param[in]  X       Location to evaluate u(X,t)
159ba6664aeSJames Wright  * @param[in]  t       Time to evaluate u(X,t)
160ba6664aeSJames Wright  * @param[in]  ubar    Mean velocity at X
161ba6664aeSJames Wright  * @param[in]  cij     Cholesky decomposition at X
162ba6664aeSJames Wright  * @param[in]  qn      Wavemode amplitudes at X, [nmodes]
163ba6664aeSJames Wright  * @param[out] u       Velocity at X and t
164ba6664aeSJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
165ba6664aeSJames Wright  */
166ba6664aeSJames Wright void CEED_QFUNCTION_HELPER(STGShur14_Calc)(const CeedScalar X[3],
167ba6664aeSJames Wright     const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6],
168ba6664aeSJames Wright     const CeedScalar qn[], CeedScalar u[3],
169ba6664aeSJames Wright     const STGShur14Context stg_ctx) {
170ba6664aeSJames Wright 
171ba6664aeSJames Wright   //*INDENT-OFF*
172ba6664aeSJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
173ba6664aeSJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
174ba6664aeSJames Wright   const CeedScalar *phi   = &stg_ctx->data[stg_ctx->offsets.phi];
175ba6664aeSJames Wright   const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma];
176ba6664aeSJames Wright   const CeedScalar *d     = &stg_ctx->data[stg_ctx->offsets.d];
177ba6664aeSJames Wright   //*INDENT-ON*
178ba6664aeSJames Wright   CeedScalar xdotd, vp[3] = {0.};
179ba6664aeSJames Wright   CeedScalar xhat[] = {0., X[1], X[2]};
180ba6664aeSJames Wright 
181ba6664aeSJames Wright   CeedPragmaSIMD
182ba6664aeSJames Wright   for(CeedInt n=0; n<nmodes; n++) {
183ba6664aeSJames Wright     xhat[0] = (X[0] - stg_ctx->u0*t)*Max(2*kappa[0]/kappa[n], 0.1);
184ba6664aeSJames Wright     xdotd = 0.;
185ba6664aeSJames Wright     for(CeedInt i=0; i<3; i++) xdotd += d[i*nmodes+n]*xhat[i];
186ba6664aeSJames Wright     const CeedScalar cos_kxdp = cos(kappa[n]*xdotd + phi[n]);
187961c9c98SJames Wright     vp[0] += sqrt(qn[n])*sigma[0*nmodes+n] * cos_kxdp;
188961c9c98SJames Wright     vp[1] += sqrt(qn[n])*sigma[1*nmodes+n] * cos_kxdp;
189961c9c98SJames Wright     vp[2] += sqrt(qn[n])*sigma[2*nmodes+n] * cos_kxdp;
190ba6664aeSJames Wright   }
191961c9c98SJames Wright   for(CeedInt i=0; i<3; i++) vp[i] *= 2*sqrt(1.5);
192ba6664aeSJames Wright 
193ba6664aeSJames Wright   u[0] = ubar[0] + cij[0]*vp[0];
194ba6664aeSJames Wright   u[1] = ubar[1] + cij[3]*vp[0] + cij[1]*vp[1];
195ba6664aeSJames Wright   u[2] = ubar[2] + cij[4]*vp[0] + cij[5]*vp[1] + cij[2]*vp[2];
196ba6664aeSJames Wright }
197ba6664aeSJames Wright 
198b77c53c9SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition
199b77c53c9SJames Wright CEED_QFUNCTION(ICsSTG)(void *ctx, CeedInt Q,
200b77c53c9SJames Wright                        const CeedScalar *const *in, CeedScalar *const *out) {
201b77c53c9SJames Wright   // Inputs
202b77c53c9SJames Wright   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
203b77c53c9SJames Wright 
204b77c53c9SJames Wright   // Outputs
205b77c53c9SJames Wright   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
206b77c53c9SJames Wright 
207b77c53c9SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
208b77c53c9SJames Wright   CeedScalar u[3], cij[6], eps, lt;
209b77c53c9SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
210b77c53c9SJames Wright   const CeedScalar P0     = stg_ctx->P0;
211b77c53c9SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
212b77c53c9SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
213b77c53c9SJames Wright   const CeedScalar Rd     = cp - cv;
214b77c53c9SJames Wright   const CeedScalar rho = P0 / (Rd * theta0);
215b77c53c9SJames Wright 
216b77c53c9SJames Wright   CeedPragmaSIMD
217b77c53c9SJames Wright   for(CeedInt i=0; i<Q; i++) {
218b77c53c9SJames Wright     InterpolateProfile(X[1][i], u, cij, &eps, &lt, stg_ctx);
219b77c53c9SJames Wright 
220b77c53c9SJames Wright     q0[0][i] = rho;
221b77c53c9SJames Wright     q0[1][i] = u[0] * rho;
222b77c53c9SJames Wright     q0[2][i] = u[1] * rho;
223b77c53c9SJames Wright     q0[3][i] = u[2] * rho;
224b77c53c9SJames Wright     q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0);
225b77c53c9SJames Wright   } // End of Quadrature Point Loop
226b77c53c9SJames Wright   return 0;
227b77c53c9SJames Wright }
228b77c53c9SJames Wright 
229ba6664aeSJames Wright /********************************************************************
230ba6664aeSJames Wright  * @brief QFunction to calculate the inflow boundary condition
231ba6664aeSJames Wright  *
232ba6664aeSJames Wright  * This will loop through quadrature points, calculate the wavemode amplitudes
233ba6664aeSJames Wright  * at each location, then calculate the actual velocity.
234ba6664aeSJames Wright  */
235ba6664aeSJames Wright CEED_QFUNCTION(STGShur14_Inflow)(void *ctx, CeedInt Q,
236ba6664aeSJames Wright                                  const CeedScalar *const *in,
237ba6664aeSJames Wright                                  CeedScalar *const *out) {
238ba6664aeSJames Wright 
239ba6664aeSJames Wright   //*INDENT-OFF*
240ba6664aeSJames Wright   const CeedScalar (*q)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
241e8b03feeSJames Wright                    (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[2],
242e8b03feeSJames Wright                    (*X)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[3];
243ba6664aeSJames Wright 
2444dbab5e5SJames Wright   CeedScalar(*v)[CEED_Q_VLA]            = (CeedScalar(*)[CEED_Q_VLA]) out[0],
2454dbab5e5SJames Wright             (*jac_data_sur)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[1];
246ba6664aeSJames Wright 
247ba6664aeSJames Wright   //*INDENT-ON*
248ba6664aeSJames Wright 
249ba6664aeSJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
250ba6664aeSJames Wright   CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
251ba6664aeSJames Wright   const bool is_implicit  = stg_ctx->is_implicit;
252ba6664aeSJames Wright   const bool mean_only    = stg_ctx->mean_only;
253ba6664aeSJames Wright   const bool prescribe_T  = stg_ctx->prescribe_T;
254ba6664aeSJames Wright   const CeedScalar dx     = stg_ctx->dx;
255ba6664aeSJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
256ba6664aeSJames Wright   const CeedScalar time   = stg_ctx->time;
257ba6664aeSJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
258ba6664aeSJames Wright   const CeedScalar P0     = stg_ctx->P0;
259ba6664aeSJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
260ba6664aeSJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
261ba6664aeSJames Wright   const CeedScalar Rd     = cp - cv;
262ba6664aeSJames Wright   const CeedScalar gamma  = cp/cv;
263ba6664aeSJames Wright 
264ba6664aeSJames Wright   CeedPragmaSIMD
265ba6664aeSJames Wright   for(CeedInt i=0; i<Q; i++) {
266ba6664aeSJames Wright     const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0);
267ba6664aeSJames Wright     const CeedScalar x[] = { X[0][i], X[1][i], X[2][i] };
268ba6664aeSJames Wright     const CeedScalar dXdx[2][3] = {
269ba6664aeSJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
270ba6664aeSJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
271ba6664aeSJames Wright     };
272ba6664aeSJames Wright 
273ba6664aeSJames Wright     CeedScalar h[3];
274ba6664aeSJames Wright     h[0] = dx;
275a939fbabSJames Wright     for (CeedInt j=1; j<3; j++)
276a939fbabSJames Wright       h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]));
277ba6664aeSJames Wright 
278ba6664aeSJames Wright     InterpolateProfile(X[1][i], ubar, cij, &eps, &lt, stg_ctx);
279ba6664aeSJames Wright     if (!mean_only) {
280ba6664aeSJames Wright       CalcSpectrum(X[1][i], eps, lt, h, mu/rho, qn, stg_ctx);
281ba6664aeSJames Wright       STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
282ba6664aeSJames Wright     } else {
283ba6664aeSJames Wright       for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
284ba6664aeSJames Wright     }
285ba6664aeSJames Wright 
2864dbab5e5SJames Wright     const CeedScalar E_kinetic = .5 * rho * Dot3(u, u);
287ba6664aeSJames Wright     CeedScalar E_internal, P;
288ba6664aeSJames Wright     if (prescribe_T) {
289ba6664aeSJames Wright       // Temperature is being set weakly (theta0) and for constant cv this sets E_internal
290ba6664aeSJames Wright       E_internal = rho * cv * theta0;
291ba6664aeSJames Wright       // Find pressure using
292ba6664aeSJames Wright       P = rho * Rd * theta0; // interior rho with exterior T
293ba6664aeSJames Wright     } else {
294ba6664aeSJames Wright       E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution
295ba6664aeSJames Wright       P = E_internal * (gamma - 1.);
296ba6664aeSJames Wright     }
297ba6664aeSJames Wright 
298ba6664aeSJames Wright     const CeedScalar wdetJb  = (is_implicit ? -1. : 1.) * q_data_sur[0][i];
299ba6664aeSJames Wright     // ---- Normal vect
300ba6664aeSJames Wright     const CeedScalar norm[3] = {q_data_sur[1][i],
301ba6664aeSJames Wright                                 q_data_sur[2][i],
302ba6664aeSJames Wright                                 q_data_sur[3][i]
303ba6664aeSJames Wright                                };
304ba6664aeSJames Wright 
305ba6664aeSJames Wright     const CeedScalar E = E_internal + E_kinetic;
306ba6664aeSJames Wright 
307ba6664aeSJames Wright     // Velocity normal to the boundary
3084dbab5e5SJames Wright     const CeedScalar u_normal = Dot3(norm, u);
3094dbab5e5SJames Wright 
310ba6664aeSJames Wright     // The Physics
311ba6664aeSJames Wright     // Zero v so all future terms can safely sum into it
312ba6664aeSJames Wright     for (CeedInt j=0; j<5; j++) v[j][i] = 0.;
313ba6664aeSJames Wright 
314ba6664aeSJames Wright     // The Physics
315ba6664aeSJames Wright     // -- Density
316ba6664aeSJames Wright     v[0][i] -= wdetJb * rho * u_normal;
317ba6664aeSJames Wright 
318ba6664aeSJames Wright     // -- Momentum
319ba6664aeSJames Wright     for (CeedInt j=0; j<3; j++)
320ba6664aeSJames Wright       v[j+1][i] -= wdetJb *(rho * u_normal * u[j] +
321ba6664aeSJames Wright                             norm[j] * P);
322ba6664aeSJames Wright 
323ba6664aeSJames Wright     // -- Total Energy Density
324ba6664aeSJames Wright     v[4][i] -= wdetJb * u_normal * (E + P);
3254dbab5e5SJames Wright 
3264dbab5e5SJames Wright     jac_data_sur[0][i] = rho;
3274dbab5e5SJames Wright     jac_data_sur[1][i] = u[0];
3284dbab5e5SJames Wright     jac_data_sur[2][i] = u[1];
3294dbab5e5SJames Wright     jac_data_sur[3][i] = u[2];
3304dbab5e5SJames Wright     jac_data_sur[4][i] = E;
3314dbab5e5SJames Wright     for (int j=0; j<6; j++) jac_data_sur[5+j][i] = 0.;
332ba6664aeSJames Wright   }
333ba6664aeSJames Wright   return 0;
334ba6664aeSJames Wright }
335ba6664aeSJames Wright 
3364dbab5e5SJames Wright CEED_QFUNCTION(STGShur14_Inflow_Jacobian)(void *ctx, CeedInt Q,
3374dbab5e5SJames Wright     const CeedScalar *const *in,
3384dbab5e5SJames Wright     CeedScalar *const *out) {
3394dbab5e5SJames Wright   // *INDENT-OFF*
3404dbab5e5SJames Wright   // Inputs
3414dbab5e5SJames Wright   const CeedScalar (*dq)[CEED_Q_VLA]           = (const CeedScalar(*)[CEED_Q_VLA])in[0],
3424dbab5e5SJames Wright                    (*q_data_sur)[CEED_Q_VLA]   = (const CeedScalar(*)[CEED_Q_VLA])in[2],
3434dbab5e5SJames Wright                    (*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4];
3444dbab5e5SJames Wright   // Outputs
3454dbab5e5SJames Wright   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
3464dbab5e5SJames Wright   // *INDENT-ON*
3474dbab5e5SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
3484dbab5e5SJames Wright   const bool implicit     = stg_ctx->is_implicit;
3494dbab5e5SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
3504dbab5e5SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
3514dbab5e5SJames Wright   const CeedScalar Rd     = cp - cv;
3524dbab5e5SJames Wright   const CeedScalar gamma  = cp/cv;
3534dbab5e5SJames Wright 
3544dbab5e5SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
3554dbab5e5SJames Wright   const bool prescribe_T  = stg_ctx->prescribe_T;
3564dbab5e5SJames Wright 
3574dbab5e5SJames Wright   CeedPragmaSIMD
3584dbab5e5SJames Wright   // Quadrature Point Loop
3594dbab5e5SJames Wright   for (CeedInt i=0; i<Q; i++) {
3604dbab5e5SJames Wright     // Setup
3614dbab5e5SJames Wright     // Setup
3624dbab5e5SJames Wright     // -- Interp-to-Interp q_data
3634dbab5e5SJames Wright     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
3644dbab5e5SJames Wright     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
3654dbab5e5SJames Wright     // We can effect this by swapping the sign on this weight
3664dbab5e5SJames Wright     const CeedScalar wdetJb  = (implicit ? -1. : 1.) * q_data_sur[0][i];
3674dbab5e5SJames Wright 
3684dbab5e5SJames Wright     // Calculate inflow values
3694dbab5e5SJames Wright     CeedScalar velocity[3];
3704dbab5e5SJames Wright     for (CeedInt j=0; j<3; j++) velocity[j] = jac_data_sur[5+j][i];
3714dbab5e5SJames Wright 
3724dbab5e5SJames Wright     // enabling user to choose between weak T and weak rho inflow
3734dbab5e5SJames Wright     CeedScalar drho, dE, dP;
3744dbab5e5SJames Wright     if (prescribe_T) {
3754dbab5e5SJames Wright       // rho should be from the current solution
3764dbab5e5SJames Wright       drho = dq[0][i];
3774dbab5e5SJames Wright       CeedScalar dE_internal = drho * cv * theta0;
3784dbab5e5SJames Wright       CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity);
3794dbab5e5SJames Wright       dE = dE_internal + dE_kinetic;
3804dbab5e5SJames Wright       dP = drho * Rd * theta0; // interior rho with exterior T
3814dbab5e5SJames Wright     } else { // rho specified, E_internal from solution
3824dbab5e5SJames Wright       drho = 0;
3834dbab5e5SJames Wright       dE = dq[4][i];
3844dbab5e5SJames Wright       dP = dE * (gamma - 1.);
3854dbab5e5SJames Wright     }
3864dbab5e5SJames Wright     const CeedScalar norm[3] = {q_data_sur[1][i],
3874dbab5e5SJames Wright                                 q_data_sur[2][i],
3884dbab5e5SJames Wright                                 q_data_sur[3][i]
3894dbab5e5SJames Wright                                };
3904dbab5e5SJames Wright 
3914dbab5e5SJames Wright     const CeedScalar u_normal = Dot3(norm, velocity);
3924dbab5e5SJames Wright 
3934dbab5e5SJames Wright     v[0][i] = - wdetJb * drho * u_normal;
3944dbab5e5SJames Wright     for (int j=0; j<3; j++)
3954dbab5e5SJames Wright       v[j+1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP);
3964dbab5e5SJames Wright     v[4][i] = - wdetJb * u_normal * (dE + dP);
3974dbab5e5SJames Wright   } // End Quadrature Point Loop
3984dbab5e5SJames Wright   return 0;
3994dbab5e5SJames Wright }
4004dbab5e5SJames Wright 
4010a6353c2SJames Wright /********************************************************************
4020a6353c2SJames Wright  * @brief QFunction to calculate the strongly enforce inflow BC
4030a6353c2SJames Wright  *
4040a6353c2SJames Wright  * This QF is for the strong application of STG via libCEED (rather than
4050a6353c2SJames Wright  * through the native PETSc `DMAddBoundary` -> `bcFunc` method.
4060a6353c2SJames Wright  */
4070a6353c2SJames Wright CEED_QFUNCTION(STGShur14_Inflow_StrongQF)(void *ctx, CeedInt Q,
4080a6353c2SJames Wright     const CeedScalar *const *in, CeedScalar *const *out) {
4090a6353c2SJames Wright 
4100a6353c2SJames Wright   //*INDENT-OFF*
4110a6353c2SJames Wright   const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
4120a6353c2SJames Wright                    (*coords)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA]) in[1],
4130a6353c2SJames Wright                    (*scale)                  = (const CeedScalar(*)) in[2];
4140a6353c2SJames Wright 
4150a6353c2SJames Wright   CeedScalar(*bcval)[CEED_Q_VLA]            = (CeedScalar(*)[CEED_Q_VLA]) out[0];
4160a6353c2SJames Wright   //*INDENT-ON*
4170a6353c2SJames Wright 
4180a6353c2SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
4190a6353c2SJames Wright   CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
4200a6353c2SJames Wright   const bool mean_only    = stg_ctx->mean_only;
4210a6353c2SJames Wright   const CeedScalar dx     = stg_ctx->dx;
4220a6353c2SJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
4230a6353c2SJames Wright   const CeedScalar time   = stg_ctx->time;
4240a6353c2SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
4250a6353c2SJames Wright   const CeedScalar P0     = stg_ctx->P0;
4260a6353c2SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
4270a6353c2SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
4280a6353c2SJames Wright   const CeedScalar Rd     = cp - cv;
4290a6353c2SJames Wright   const CeedScalar rho    = P0 / (Rd * theta0);
4300a6353c2SJames Wright 
4310a6353c2SJames Wright   CeedPragmaSIMD
4320a6353c2SJames Wright   for(CeedInt i=0; i<Q; i++) {
4330a6353c2SJames Wright     const CeedScalar x[] = { coords[0][i], coords[1][i], coords[2][i] };
4340a6353c2SJames Wright     const CeedScalar dXdx[2][3] = {
4350a6353c2SJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
4360a6353c2SJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
4370a6353c2SJames Wright     };
4380a6353c2SJames Wright 
4390a6353c2SJames Wright     CeedScalar h[3];
4400a6353c2SJames Wright     h[0] = dx;
441a939fbabSJames Wright     for (CeedInt j=1; j<3; j++)
442a939fbabSJames Wright       h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]));
4430a6353c2SJames Wright 
4440a6353c2SJames Wright     InterpolateProfile(coords[1][i], ubar, cij, &eps, &lt, stg_ctx);
4450a6353c2SJames Wright     if (!mean_only) {
4460a6353c2SJames Wright       CalcSpectrum(coords[1][i], eps, lt, h, mu/rho, qn, stg_ctx);
4470a6353c2SJames Wright       STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
4480a6353c2SJames Wright     } else {
4490a6353c2SJames Wright       for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
4500a6353c2SJames Wright     }
4510a6353c2SJames Wright 
4520a6353c2SJames Wright     bcval[0][i] = scale[i] * rho;
4530a6353c2SJames Wright     bcval[1][i] = scale[i] * rho * u[0];
4540a6353c2SJames Wright     bcval[2][i] = scale[i] * rho * u[1];
4550a6353c2SJames Wright     bcval[3][i] = scale[i] * rho * u[2];
456cf3d54ffSJames Wright     bcval[4][i] = 0.;
4570a6353c2SJames Wright   }
4580a6353c2SJames Wright   return 0;
4590a6353c2SJames Wright }
4600a6353c2SJames Wright 
461ba6664aeSJames Wright #endif // stg_shur14_h
462