xref: /libCEED/examples/fluids/qfunctions/stg_shur14.h (revision 175f00a62957e4f3f8aaa73c4289f1b68f0117f5)
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 <ceed.h>
20c9c2c079SJeremy L Thompson #include <math.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  *
30*175f00a6SJames Wright  * Assumed that prof_wd[i+1] > prof_wd[i] and prof_wd[0] = 0
31*175f00a6SJames Wright  * If wall_dist > prof_wd[-1], then the interpolation takes the values at prof_wd[-1]
32ba6664aeSJames Wright  *
33*175f00a6SJames Wright  * @param[in]  wall_dist Distance to the nearest wall
34*175f00a6SJames Wright  * @param[out] ubar      Mean velocity at wall_dist
35*175f00a6SJames Wright  * @param[out] cij       Cholesky decomposition at wall_dist
36*175f00a6SJames Wright  * @param[out] eps       Turbulent dissipation at wall_dist
37*175f00a6SJames Wright  * @param[out] lt        Turbulent length scale at wall_dist
38ba6664aeSJames Wright  * @param[in]  stg_ctx   STGShur14Context for the problem
39ba6664aeSJames Wright  */
40*175f00a6SJames Wright CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar wall_dist,
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;
45*175f00a6SJames Wright   const CeedScalar *prof_wd   = &stg_ctx->data[stg_ctx->offsets.wall_dist];
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++) {
53*175f00a6SJames Wright     if (wall_dist < prof_wd[i]) {
54ba6664aeSJames Wright       idx = i;
55ba6664aeSJames Wright       break;
56ba6664aeSJames Wright     }
57ba6664aeSJames Wright   }
58ba6664aeSJames Wright 
59*175f00a6SJames Wright   if (idx > 0) { // y within the bounds of prof_wd
60ba6664aeSJames Wright     //*INDENT-OFF*
61*175f00a6SJames Wright     CeedScalar coeff = (wall_dist - prof_wd[idx-1]) / (prof_wd[idx] - prof_wd[idx -1]);
62*175f00a6SJames Wright 
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*
75*175f00a6SJames Wright   } else { // y outside bounds of prof_wd
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 /*
91e159aeacSJames Wright  * @brief Calculate spectrum coefficient, qn
92e159aeacSJames Wright  *
93e159aeacSJames Wright  * Calculates q_n at a given distance to the wall
94e159aeacSJames Wright  *
95e159aeacSJames Wright  * @param[in]  kappa  nth wavenumber
96e159aeacSJames Wright  * @param[in]  dkappa Difference between wavenumbers
97e159aeacSJames Wright  * @param[in]  keta   Dissipation wavenumber
98e159aeacSJames Wright  * @param[in]  kcut   Mesh-induced cutoff wavenumber
99e159aeacSJames Wright  * @param[in]  ke     Energy-containing wavenumber
100e159aeacSJames Wright  * @param[in]  Ektot  Total turbulent kinetic energy of spectrum
101e159aeacSJames Wright  * @returns    qn     Spectrum coefficient
102e159aeacSJames Wright  */
103e159aeacSJames Wright CeedScalar CEED_QFUNCTION_HELPER(Calc_qn)(const CeedScalar kappa,
104e159aeacSJames Wright     const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut,
10562e628f8SJames Wright     const CeedScalar ke, const CeedScalar Ektot_inv) {
106e159aeacSJames Wright   const CeedScalar feta_x_fcut   = exp(-Square(12*kappa/keta)
107e159aeacSJames Wright                                        -Cube(4*Max(kappa - 0.9*kcut, 0)/kcut) );
108e159aeacSJames Wright   return pow(kappa/ke, 4.) * pow(1 + 2.4*Square(kappa/ke),-17./6)
10962e628f8SJames Wright          *feta_x_fcut*dkappa * Ektot_inv;
110e159aeacSJames Wright }
111e159aeacSJames Wright 
112e159aeacSJames Wright // Calculate hmax, ke, keta, and kcut
113*175f00a6SJames Wright void CEED_QFUNCTION_HELPER(SpectrumConstants)(const CeedScalar wall_dist,
114e159aeacSJames Wright     const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3],
115e159aeacSJames Wright     const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke,
116e159aeacSJames Wright     CeedScalar *keta, CeedScalar *kcut) {
117e159aeacSJames Wright   *hmax = Max( Max(h[0], h[1]), h[2]);
118*175f00a6SJames Wright   *ke   = wall_dist==0 ? 1e16 : 2*M_PI/Min(2*wall_dist, 3*lt);
119e159aeacSJames Wright   *keta = 2*M_PI*pow(Cube(nu)/eps, -0.25);
120*175f00a6SJames Wright   *kcut = M_PI/ Min( Max(Max(h[1], h[2]), 0.3*(*hmax)) + 0.1*wall_dist, *hmax );
121e159aeacSJames Wright }
122e159aeacSJames Wright 
123e159aeacSJames 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  *
128*175f00a6SJames Wright  * @param[in]  wall_dist Distance to the nearest wall
129*175f00a6SJames Wright  * @param[in]  eps       Turbulent dissipation w/rt wall_dist
130*175f00a6SJames Wright  * @param[in]  lt        Turbulent length scale w/rt wall_dist
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  */
136*175f00a6SJames Wright void CEED_QFUNCTION_HELPER(CalcSpectrum)(const CeedScalar wall_dist,
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];
142e159aeacSJames Wright   CeedScalar hmax, ke, keta, kcut, Ektot=0.0;
143*175f00a6SJames Wright   SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut);
144ba6664aeSJames Wright 
145ba6664aeSJames Wright   for(CeedInt n=0; n<nmodes; n++) {
146e159aeacSJames Wright     const CeedScalar dkappa = n==0 ? kappa[0] : kappa[n] - kappa[n-1];
147e159aeacSJames 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 
198b277271eSJames Wright /******************************************************
199b277271eSJames Wright  * @brief Calculate u(x,t) for STG inflow condition
200b277271eSJames Wright  *
201b277271eSJames Wright  * @param[in]  X         Location to evaluate u(X,t)
202b277271eSJames Wright  * @param[in]  t         Time to evaluate u(X,t)
203b277271eSJames Wright  * @param[in]  ubar      Mean velocity at X
204b277271eSJames Wright  * @param[in]  cij       Cholesky decomposition at X
205*175f00a6SJames Wright  * @param[in]  Ektot     Total spectrum energy at this location
206*175f00a6SJames Wright  * @param[in]  h         Element size in 3 directions
207*175f00a6SJames Wright  * @param[in]  wall_dist Distance to closest wall
208*175f00a6SJames Wright  * @param[in]  eps       Turbulent dissipation
209*175f00a6SJames Wright  * @param[in]  lt        Turbulent length scale
210b277271eSJames Wright  * @param[out] u         Velocity at X and t
211b277271eSJames Wright  * @param[in]  stg_ctx   STGShur14Context for the problem
212b277271eSJames Wright  */
213b277271eSJames Wright void CEED_QFUNCTION_HELPER(STGShur14_Calc_PrecompEktot)(const CeedScalar X[3],
214b277271eSJames Wright     const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6],
215*175f00a6SJames Wright     const CeedScalar Ektot, const CeedScalar h[3], const CeedScalar wall_dist,
216b277271eSJames Wright     const CeedScalar eps, const CeedScalar lt, const CeedScalar nu, CeedScalar u[3],
217b277271eSJames Wright     const STGShur14Context stg_ctx) {
218b277271eSJames Wright 
219b277271eSJames Wright   //*INDENT-OFF*
220b277271eSJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
221b277271eSJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
222b277271eSJames Wright   const CeedScalar *phi   = &stg_ctx->data[stg_ctx->offsets.phi];
223b277271eSJames Wright   const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma];
224b277271eSJames Wright   const CeedScalar *d     = &stg_ctx->data[stg_ctx->offsets.d];
225b277271eSJames Wright   //*INDENT-ON*
226b277271eSJames Wright   CeedScalar hmax, ke, keta, kcut;
227*175f00a6SJames Wright   SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut);
228b277271eSJames Wright   CeedScalar xdotd, vp[3] = {0.};
229b277271eSJames Wright   CeedScalar xhat[] = {0., X[1], X[2]};
230b277271eSJames Wright 
231b277271eSJames Wright   CeedPragmaSIMD
232b277271eSJames Wright   for(CeedInt n=0; n<nmodes; n++) {
233b277271eSJames Wright     xhat[0] = (X[0] - stg_ctx->u0*t)*Max(2*kappa[0]/kappa[n], 0.1);
234b277271eSJames Wright     xdotd = 0.;
235b277271eSJames Wright     for(CeedInt i=0; i<3; i++) xdotd += d[i*nmodes+n]*xhat[i];
236b277271eSJames Wright     const CeedScalar cos_kxdp = cos(kappa[n]*xdotd + phi[n]);
237b277271eSJames Wright     const CeedScalar dkappa   = n==0 ? kappa[0] : kappa[n] - kappa[n-1];
238b277271eSJames Wright     const CeedScalar qn       = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot);
239b277271eSJames Wright     vp[0] += sqrt(qn)*sigma[0*nmodes+n] * cos_kxdp;
240b277271eSJames Wright     vp[1] += sqrt(qn)*sigma[1*nmodes+n] * cos_kxdp;
241b277271eSJames Wright     vp[2] += sqrt(qn)*sigma[2*nmodes+n] * cos_kxdp;
242b277271eSJames Wright   }
243b277271eSJames Wright   for(CeedInt i=0; i<3; i++) vp[i] *= 2*sqrt(1.5);
244b277271eSJames Wright 
245b277271eSJames Wright   u[0] = ubar[0] + cij[0]*vp[0];
246b277271eSJames Wright   u[1] = ubar[1] + cij[3]*vp[0] + cij[1]*vp[1];
247b277271eSJames Wright   u[2] = ubar[2] + cij[4]*vp[0] + cij[5]*vp[1] + cij[2]*vp[2];
248b277271eSJames Wright }
249b277271eSJames Wright 
25062e628f8SJames Wright // Create preprocessed input for the stg calculation
25162e628f8SJames Wright //
25262e628f8SJames Wright // stg_data[0] = 1 / Ektot (inverse of total spectrum energy)
253b277271eSJames Wright CEED_QFUNCTION(Preprocess_STGShur14)(void *ctx, CeedInt Q,
254b277271eSJames Wright                                      const CeedScalar *const *in, CeedScalar *const *out) {
255b277271eSJames Wright   //*INDENT-OFF*
256b277271eSJames Wright   const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
257b277271eSJames Wright                    (*x)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[1];
258b277271eSJames Wright 
259b277271eSJames Wright   CeedScalar (*stg_data) = (CeedScalar(*)) out[0];
260b277271eSJames Wright 
261b277271eSJames Wright   //*INDENT-ON*
262b277271eSJames Wright   CeedScalar ubar[3], cij[6], eps, lt;
263b277271eSJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
264b277271eSJames Wright   const CeedScalar dx     = stg_ctx->dx;
265b277271eSJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
266b277271eSJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
267b277271eSJames Wright   const CeedScalar P0     = stg_ctx->P0;
268b277271eSJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
269b277271eSJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
270b277271eSJames Wright   const CeedScalar Rd     = cp - cv;
271b277271eSJames Wright   const CeedScalar rho    = P0 / (Rd * theta0);
272b277271eSJames Wright   const CeedScalar nu     = mu / rho;
273b277271eSJames Wright 
274b277271eSJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
275b277271eSJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
2765dc40723SJames Wright   CeedScalar hmax, ke, keta, kcut;
277b277271eSJames Wright 
278b277271eSJames Wright   CeedPragmaSIMD
279b277271eSJames Wright   for(CeedInt i=0; i<Q; i++) {
280*175f00a6SJames Wright     const CeedScalar wall_dist = x[1][i];
281b277271eSJames Wright     const CeedScalar dXdx[2][3] = {
282b277271eSJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
283b277271eSJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
284b277271eSJames Wright     };
285b277271eSJames Wright 
286b277271eSJames Wright     CeedScalar h[3];
287b277271eSJames Wright     h[0] = dx;
288b277271eSJames Wright     for (CeedInt j=1; j<3; j++)
289b277271eSJames Wright       h[j] = 2/sqrt(dXdx[0][j]*dXdx[0][j] + dXdx[1][j]*dXdx[1][j]);
290b277271eSJames Wright 
291*175f00a6SJames Wright     InterpolateProfile(wall_dist, ubar, cij, &eps, &lt, stg_ctx);
292*175f00a6SJames Wright     SpectrumConstants(wall_dist, eps, lt, h, nu, &hmax, &ke, &keta, &kcut);
293b277271eSJames Wright 
294b277271eSJames Wright     // Calculate total TKE per spectrum
295d97dc904SJames Wright     CeedScalar Ek_tot=0;
296b277271eSJames Wright     CeedPragmaSIMD
297b277271eSJames Wright     for(CeedInt n=0; n<nmodes; n++) {
298b277271eSJames Wright       const CeedScalar dkappa = n==0 ? kappa[0] : kappa[n] - kappa[n-1];
299d97dc904SJames Wright       Ek_tot += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0);
300b277271eSJames Wright     }
301d97dc904SJames Wright     // avoid underflowed and poorly defined spectrum coefficients
302d97dc904SJames Wright     stg_data[i] = Ek_tot != 0 ? 1/Ek_tot : 0;
303b277271eSJames Wright   }
304b277271eSJames Wright   return 0;
305b277271eSJames Wright }
306b277271eSJames Wright 
307b77c53c9SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition
308b77c53c9SJames Wright CEED_QFUNCTION(ICsSTG)(void *ctx, CeedInt Q,
309b77c53c9SJames Wright                        const CeedScalar *const *in, CeedScalar *const *out) {
310b77c53c9SJames Wright   // Inputs
311b77c53c9SJames Wright   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
312b77c53c9SJames Wright 
313b77c53c9SJames Wright   // Outputs
314b77c53c9SJames Wright   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
315b77c53c9SJames Wright 
316b77c53c9SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
317b77c53c9SJames Wright   CeedScalar u[3], cij[6], eps, lt;
318b77c53c9SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
319b77c53c9SJames Wright   const CeedScalar P0     = stg_ctx->P0;
320b77c53c9SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
321b77c53c9SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
322b77c53c9SJames Wright   const CeedScalar Rd     = cp - cv;
323b77c53c9SJames Wright   const CeedScalar rho    = P0 / (Rd * theta0);
324b77c53c9SJames Wright 
325b77c53c9SJames Wright   CeedPragmaSIMD
326b77c53c9SJames Wright   for(CeedInt i=0; i<Q; i++) {
327b77c53c9SJames Wright     InterpolateProfile(X[1][i], u, cij, &eps, &lt, stg_ctx);
328b77c53c9SJames Wright 
32997baf651SJames Wright     switch (stg_ctx->newtonian_ctx.state_var) {
33097baf651SJames Wright     case STATEVAR_CONSERVATIVE:
331b77c53c9SJames Wright       q0[0][i] = rho;
332b77c53c9SJames Wright       q0[1][i] = u[0] * rho;
333b77c53c9SJames Wright       q0[2][i] = u[1] * rho;
334b77c53c9SJames Wright       q0[3][i] = u[2] * rho;
335b77c53c9SJames Wright       q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0);
33697baf651SJames Wright       break;
33797baf651SJames Wright 
33897baf651SJames Wright     case STATEVAR_PRIMITIVE:
33997baf651SJames Wright       q0[0][i] = P0;
34097baf651SJames Wright       q0[1][i] = u[0];
34197baf651SJames Wright       q0[2][i] = u[1];
34297baf651SJames Wright       q0[3][i] = u[2];
34397baf651SJames Wright       q0[4][i] = theta0;
34497baf651SJames Wright       break;
3457c4551aaSJames Wright     }
346b77c53c9SJames Wright   } // End of Quadrature Point Loop
347b77c53c9SJames Wright   return 0;
348b77c53c9SJames Wright }
349b77c53c9SJames Wright 
350ba6664aeSJames Wright /********************************************************************
351ba6664aeSJames Wright  * @brief QFunction to calculate the inflow boundary condition
352ba6664aeSJames Wright  *
353ba6664aeSJames Wright  * This will loop through quadrature points, calculate the wavemode amplitudes
354ba6664aeSJames Wright  * at each location, then calculate the actual velocity.
355ba6664aeSJames Wright  */
356ba6664aeSJames Wright CEED_QFUNCTION(STGShur14_Inflow)(void *ctx, CeedInt Q,
35797baf651SJames Wright                                  const CeedScalar *const *in, CeedScalar *const *out) {
358ba6664aeSJames Wright 
359ba6664aeSJames Wright   //*INDENT-OFF*
360ba6664aeSJames Wright   const CeedScalar (*q)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
361e8b03feeSJames Wright                    (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[2],
362e8b03feeSJames Wright                    (*X)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[3];
363ba6664aeSJames Wright 
3644dbab5e5SJames Wright   CeedScalar(*v)[CEED_Q_VLA]            = (CeedScalar(*)[CEED_Q_VLA]) out[0],
3654dbab5e5SJames Wright             (*jac_data_sur)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[1];
366ba6664aeSJames Wright 
367ba6664aeSJames Wright   //*INDENT-ON*
368ba6664aeSJames Wright 
369ba6664aeSJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
370ba6664aeSJames Wright   CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
371ba6664aeSJames Wright   const bool is_implicit  = stg_ctx->is_implicit;
372ba6664aeSJames Wright   const bool mean_only    = stg_ctx->mean_only;
373ba6664aeSJames Wright   const bool prescribe_T  = stg_ctx->prescribe_T;
374ba6664aeSJames Wright   const CeedScalar dx     = stg_ctx->dx;
375ba6664aeSJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
376ba6664aeSJames Wright   const CeedScalar time   = stg_ctx->time;
377ba6664aeSJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
378ba6664aeSJames Wright   const CeedScalar P0     = stg_ctx->P0;
379ba6664aeSJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
380ba6664aeSJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
381ba6664aeSJames Wright   const CeedScalar Rd     = cp - cv;
382ba6664aeSJames Wright   const CeedScalar gamma  = cp/cv;
383ba6664aeSJames Wright 
384ba6664aeSJames Wright   CeedPragmaSIMD
385ba6664aeSJames Wright   for(CeedInt i=0; i<Q; i++) {
386ba6664aeSJames Wright     const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0);
387ba6664aeSJames Wright     const CeedScalar x[] = { X[0][i], X[1][i], X[2][i] };
388ba6664aeSJames Wright     const CeedScalar dXdx[2][3] = {
389ba6664aeSJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
390ba6664aeSJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
391ba6664aeSJames Wright     };
392ba6664aeSJames Wright 
393ba6664aeSJames Wright     CeedScalar h[3];
394ba6664aeSJames Wright     h[0] = dx;
395a939fbabSJames Wright     for (CeedInt j=1; j<3; j++)
396a939fbabSJames Wright       h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]));
397ba6664aeSJames Wright 
398ba6664aeSJames Wright     InterpolateProfile(X[1][i], ubar, cij, &eps, &lt, stg_ctx);
399ba6664aeSJames Wright     if (!mean_only) {
400ba6664aeSJames Wright       CalcSpectrum(X[1][i], eps, lt, h, mu/rho, qn, stg_ctx);
401ba6664aeSJames Wright       STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
402ba6664aeSJames Wright     } else {
403ba6664aeSJames Wright       for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
404ba6664aeSJames Wright     }
405ba6664aeSJames Wright 
4064dbab5e5SJames Wright     const CeedScalar E_kinetic = .5 * rho * Dot3(u, u);
407ba6664aeSJames Wright     CeedScalar E_internal, P;
408ba6664aeSJames Wright     if (prescribe_T) {
409ba6664aeSJames Wright       // Temperature is being set weakly (theta0) and for constant cv this sets E_internal
410ba6664aeSJames Wright       E_internal = rho * cv * theta0;
411ba6664aeSJames Wright       // Find pressure using
412ba6664aeSJames Wright       P = rho * Rd * theta0; // interior rho with exterior T
413ba6664aeSJames Wright     } else {
414ba6664aeSJames Wright       E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution
415ba6664aeSJames Wright       P = E_internal * (gamma - 1.);
416ba6664aeSJames Wright     }
417ba6664aeSJames Wright 
418ba6664aeSJames Wright     const CeedScalar wdetJb  = (is_implicit ? -1. : 1.) * q_data_sur[0][i];
419ba6664aeSJames Wright     // ---- Normal vect
420ba6664aeSJames Wright     const CeedScalar norm[3] = {q_data_sur[1][i],
421ba6664aeSJames Wright                                 q_data_sur[2][i],
422ba6664aeSJames Wright                                 q_data_sur[3][i]
423ba6664aeSJames Wright                                };
424ba6664aeSJames Wright 
425ba6664aeSJames Wright     const CeedScalar E = E_internal + E_kinetic;
426ba6664aeSJames Wright 
427ba6664aeSJames Wright     // Velocity normal to the boundary
4284dbab5e5SJames Wright     const CeedScalar u_normal = Dot3(norm, u);
4294dbab5e5SJames Wright 
430ba6664aeSJames Wright     // The Physics
431ba6664aeSJames Wright     // Zero v so all future terms can safely sum into it
432ba6664aeSJames Wright     for (CeedInt j=0; j<5; j++) v[j][i] = 0.;
433ba6664aeSJames Wright 
434ba6664aeSJames Wright     // The Physics
435ba6664aeSJames Wright     // -- Density
436ba6664aeSJames Wright     v[0][i] -= wdetJb * rho * u_normal;
437ba6664aeSJames Wright 
438ba6664aeSJames Wright     // -- Momentum
439ba6664aeSJames Wright     for (CeedInt j=0; j<3; j++)
440ba6664aeSJames Wright       v[j+1][i] -= wdetJb *(rho * u_normal * u[j] +
441ba6664aeSJames Wright                             norm[j] * P);
442ba6664aeSJames Wright 
443ba6664aeSJames Wright     // -- Total Energy Density
444ba6664aeSJames Wright     v[4][i] -= wdetJb * u_normal * (E + P);
4454dbab5e5SJames Wright 
4464dbab5e5SJames Wright     jac_data_sur[0][i] = rho;
4474dbab5e5SJames Wright     jac_data_sur[1][i] = u[0];
4484dbab5e5SJames Wright     jac_data_sur[2][i] = u[1];
4494dbab5e5SJames Wright     jac_data_sur[3][i] = u[2];
4504dbab5e5SJames Wright     jac_data_sur[4][i] = E;
4514dbab5e5SJames Wright     for (int j=0; j<6; j++) jac_data_sur[5+j][i] = 0.;
452ba6664aeSJames Wright   }
453ba6664aeSJames Wright   return 0;
454ba6664aeSJames Wright }
455ba6664aeSJames Wright 
4564dbab5e5SJames Wright CEED_QFUNCTION(STGShur14_Inflow_Jacobian)(void *ctx, CeedInt Q,
45797baf651SJames Wright     const CeedScalar *const *in, CeedScalar *const *out) {
4584dbab5e5SJames Wright   // *INDENT-OFF*
4594dbab5e5SJames Wright   // Inputs
4604dbab5e5SJames Wright   const CeedScalar (*dq)[CEED_Q_VLA]           = (const CeedScalar(*)[CEED_Q_VLA])in[0],
4614dbab5e5SJames Wright                    (*q_data_sur)[CEED_Q_VLA]   = (const CeedScalar(*)[CEED_Q_VLA])in[2],
4624dbab5e5SJames Wright                    (*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4];
4634dbab5e5SJames Wright   // Outputs
4644dbab5e5SJames Wright   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
4654dbab5e5SJames Wright   // *INDENT-ON*
4664dbab5e5SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
4674dbab5e5SJames Wright   const bool implicit     = stg_ctx->is_implicit;
4684dbab5e5SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
4694dbab5e5SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
4704dbab5e5SJames Wright   const CeedScalar Rd     = cp - cv;
4714dbab5e5SJames Wright   const CeedScalar gamma  = cp/cv;
4724dbab5e5SJames Wright 
4734dbab5e5SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
4744dbab5e5SJames Wright   const bool prescribe_T  = stg_ctx->prescribe_T;
4754dbab5e5SJames Wright 
4764dbab5e5SJames Wright   CeedPragmaSIMD
4774dbab5e5SJames Wright   // Quadrature Point Loop
4784dbab5e5SJames Wright   for (CeedInt i=0; i<Q; i++) {
4794dbab5e5SJames Wright     // Setup
4804dbab5e5SJames Wright     // -- Interp-to-Interp q_data
4814dbab5e5SJames Wright     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
4824dbab5e5SJames Wright     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
4834dbab5e5SJames Wright     // We can effect this by swapping the sign on this weight
4844dbab5e5SJames Wright     const CeedScalar wdetJb  = (implicit ? -1. : 1.) * q_data_sur[0][i];
4854dbab5e5SJames Wright 
4864dbab5e5SJames Wright     // Calculate inflow values
4874dbab5e5SJames Wright     CeedScalar velocity[3];
4884dbab5e5SJames Wright     for (CeedInt j=0; j<3; j++) velocity[j] = jac_data_sur[5+j][i];
4894dbab5e5SJames Wright 
4904dbab5e5SJames Wright     // enabling user to choose between weak T and weak rho inflow
4914dbab5e5SJames Wright     CeedScalar drho, dE, dP;
4924dbab5e5SJames Wright     if (prescribe_T) {
4934dbab5e5SJames Wright       // rho should be from the current solution
4944dbab5e5SJames Wright       drho = dq[0][i];
4954dbab5e5SJames Wright       CeedScalar dE_internal = drho * cv * theta0;
4964dbab5e5SJames Wright       CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity);
4974dbab5e5SJames Wright       dE = dE_internal + dE_kinetic;
4984dbab5e5SJames Wright       dP = drho * Rd * theta0; // interior rho with exterior T
4994dbab5e5SJames Wright     } else { // rho specified, E_internal from solution
5004dbab5e5SJames Wright       drho = 0;
5014dbab5e5SJames Wright       dE = dq[4][i];
5024dbab5e5SJames Wright       dP = dE * (gamma - 1.);
5034dbab5e5SJames Wright     }
5044dbab5e5SJames Wright     const CeedScalar norm[3] = {q_data_sur[1][i],
5054dbab5e5SJames Wright                                 q_data_sur[2][i],
5064dbab5e5SJames Wright                                 q_data_sur[3][i]
5074dbab5e5SJames Wright                                };
5084dbab5e5SJames Wright 
5094dbab5e5SJames Wright     const CeedScalar u_normal = Dot3(norm, velocity);
5104dbab5e5SJames Wright 
5114dbab5e5SJames Wright     v[0][i] = - wdetJb * drho * u_normal;
5124dbab5e5SJames Wright     for (int j=0; j<3; j++)
5134dbab5e5SJames Wright       v[j+1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP);
5144dbab5e5SJames Wright     v[4][i] = - wdetJb * u_normal * (dE + dP);
5154dbab5e5SJames Wright   } // End Quadrature Point Loop
5164dbab5e5SJames Wright   return 0;
5174dbab5e5SJames Wright }
5184dbab5e5SJames Wright 
5190a6353c2SJames Wright /********************************************************************
5200a6353c2SJames Wright  * @brief QFunction to calculate the strongly enforce inflow BC
5210a6353c2SJames Wright  *
5220a6353c2SJames Wright  * This QF is for the strong application of STG via libCEED (rather than
5230a6353c2SJames Wright  * through the native PETSc `DMAddBoundary` -> `bcFunc` method.
5240a6353c2SJames Wright  */
5250a6353c2SJames Wright CEED_QFUNCTION(STGShur14_Inflow_StrongQF)(void *ctx, CeedInt Q,
5260a6353c2SJames Wright     const CeedScalar *const *in, CeedScalar *const *out) {
5270a6353c2SJames Wright 
5280a6353c2SJames Wright   //*INDENT-OFF*
5290a6353c2SJames Wright   const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
5300a6353c2SJames Wright                    (*coords)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA]) in[1],
5315dc40723SJames Wright                    (*scale)                  = (const CeedScalar(*)) in[2],
5325dc40723SJames Wright                    (*stg_data)[CEED_Q_VLA]   = (const CeedScalar(*)[CEED_Q_VLA]) in[3];
5330a6353c2SJames Wright 
5340a6353c2SJames Wright   CeedScalar(*bcval)[CEED_Q_VLA]            = (CeedScalar(*)[CEED_Q_VLA]) out[0];
5350a6353c2SJames Wright   //*INDENT-ON*
5360a6353c2SJames Wright 
5370a6353c2SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
53862e628f8SJames Wright   CeedScalar u[3], ubar[3], cij[6], eps, lt;
5390a6353c2SJames Wright   const bool mean_only    = stg_ctx->mean_only;
5400a6353c2SJames Wright   const CeedScalar dx     = stg_ctx->dx;
5410a6353c2SJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
5420a6353c2SJames Wright   const CeedScalar time   = stg_ctx->time;
5430a6353c2SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
5440a6353c2SJames Wright   const CeedScalar P0     = stg_ctx->P0;
5450a6353c2SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
5460a6353c2SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
5470a6353c2SJames Wright   const CeedScalar Rd     = cp - cv;
5480a6353c2SJames Wright   const CeedScalar rho    = P0 / (Rd * theta0);
5490a6353c2SJames Wright 
5500a6353c2SJames Wright   CeedPragmaSIMD
5510a6353c2SJames Wright   for(CeedInt i=0; i<Q; i++) {
5520a6353c2SJames Wright     const CeedScalar x[] = { coords[0][i], coords[1][i], coords[2][i] };
5530a6353c2SJames Wright     const CeedScalar dXdx[2][3] = {
5540a6353c2SJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
5550a6353c2SJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
5560a6353c2SJames Wright     };
5570a6353c2SJames Wright 
5580a6353c2SJames Wright     CeedScalar h[3];
5590a6353c2SJames Wright     h[0] = dx;
560a939fbabSJames Wright     for (CeedInt j=1; j<3; j++)
561a939fbabSJames Wright       h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]));
5620a6353c2SJames Wright 
5630a6353c2SJames Wright     InterpolateProfile(coords[1][i], ubar, cij, &eps, &lt, stg_ctx);
5640a6353c2SJames Wright     if (!mean_only) {
56562e628f8SJames Wright       if (1) {
5665dc40723SJames Wright         STGShur14_Calc_PrecompEktot(x, time, ubar, cij, stg_data[0][i],
5675dc40723SJames Wright                                     h, x[1], eps, lt, mu/rho, u, stg_ctx);
56862e628f8SJames Wright       } else { // Original way
56962e628f8SJames Wright         CeedScalar qn[STG_NMODES_MAX];
57062e628f8SJames Wright         CalcSpectrum(coords[1][i], eps, lt, h, mu/rho, qn, stg_ctx);
57162e628f8SJames Wright         STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
57262e628f8SJames Wright       }
5730a6353c2SJames Wright     } else {
5740a6353c2SJames Wright       for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
5750a6353c2SJames Wright     }
5760a6353c2SJames Wright 
57797baf651SJames Wright     switch (stg_ctx->newtonian_ctx.state_var) {
57897baf651SJames Wright     case STATEVAR_CONSERVATIVE:
5790a6353c2SJames Wright       bcval[0][i] = scale[i] * rho;
5800a6353c2SJames Wright       bcval[1][i] = scale[i] * rho * u[0];
5810a6353c2SJames Wright       bcval[2][i] = scale[i] * rho * u[1];
5820a6353c2SJames Wright       bcval[3][i] = scale[i] * rho * u[2];
583cf3d54ffSJames Wright       bcval[4][i] = 0.;
58497baf651SJames Wright       break;
58597baf651SJames Wright 
58697baf651SJames Wright     case STATEVAR_PRIMITIVE:
58797baf651SJames Wright       bcval[0][i] = 0;
58897baf651SJames Wright       bcval[1][i] = scale[i] * u[0];
58997baf651SJames Wright       bcval[2][i] = scale[i] * u[1];
59097baf651SJames Wright       bcval[3][i] = scale[i] * u[2];
59197baf651SJames Wright       bcval[4][i] = scale[i] * theta0;
59297baf651SJames Wright       break;
5930a6353c2SJames Wright     }
5947c4551aaSJames Wright   }
5950a6353c2SJames Wright   return 0;
5960a6353c2SJames Wright }
5970a6353c2SJames Wright 
598ba6664aeSJames Wright #endif // stg_shur14_h
599