xref: /libCEED/examples/fluids/qfunctions/stg_shur14.h (revision b277271efb1cc32fd1d2796fb2c801bd089feba2)
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 /*
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,
105e159aeacSJames Wright     const CeedScalar ke, const CeedScalar Ektot) {
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)
109e159aeacSJames Wright          *feta_x_fcut*dkappa/Ektot;
110e159aeacSJames Wright }
111e159aeacSJames Wright 
112e159aeacSJames Wright // Calculate hmax, ke, keta, and kcut
113e159aeacSJames Wright void CEED_QFUNCTION_HELPER(SpectrumConstants)(const CeedScalar dw,
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]);
118e159aeacSJames Wright   *ke   = dw==0 ? 1e16 : 2*M_PI/Min(2*dw, 3*lt);
119e159aeacSJames Wright   *keta = 2*M_PI*pow(Cube(nu)/eps, -0.25);
120e159aeacSJames Wright   *kcut = M_PI/ Min( Max(Max(h[1], h[2]), 0.3*(*hmax)) + 0.1*dw, *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  *
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];
142e159aeacSJames Wright   CeedScalar hmax, ke, keta, kcut, Ektot=0.0;
143e159aeacSJames Wright   SpectrumConstants(dw, 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 
198*b277271eSJames Wright /******************************************************
199*b277271eSJames Wright  * @brief Calculate u(x,t) for STG inflow condition
200*b277271eSJames Wright  *
201*b277271eSJames Wright  * @param[in]  X       Location to evaluate u(X,t)
202*b277271eSJames Wright  * @param[in]  t       Time to evaluate u(X,t)
203*b277271eSJames Wright  * @param[in]  ubar    Mean velocity at X
204*b277271eSJames Wright  * @param[in]  cij     Cholesky decomposition at X
205*b277271eSJames Wright  * @param[in]  qn      Wavemode amplitudes at X, [nmodes]
206*b277271eSJames Wright  * @param[out] u       Velocity at X and t
207*b277271eSJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
208*b277271eSJames Wright  */
209*b277271eSJames Wright void CEED_QFUNCTION_HELPER(STGShur14_Calc_PrecompEktot)(const CeedScalar X[3],
210*b277271eSJames Wright     const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6],
211*b277271eSJames Wright     const CeedScalar Ektot, const CeedScalar h[3], const CeedScalar dw,
212*b277271eSJames Wright     const CeedScalar eps, const CeedScalar lt, const CeedScalar nu, CeedScalar u[3],
213*b277271eSJames Wright     const STGShur14Context stg_ctx) {
214*b277271eSJames Wright 
215*b277271eSJames Wright   //*INDENT-OFF*
216*b277271eSJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
217*b277271eSJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
218*b277271eSJames Wright   const CeedScalar *phi   = &stg_ctx->data[stg_ctx->offsets.phi];
219*b277271eSJames Wright   const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma];
220*b277271eSJames Wright   const CeedScalar *d     = &stg_ctx->data[stg_ctx->offsets.d];
221*b277271eSJames Wright   //*INDENT-ON*
222*b277271eSJames Wright   CeedScalar hmax, ke, keta, kcut;
223*b277271eSJames Wright   SpectrumConstants(dw, eps, lt, h, nu, &hmax, &ke, &keta, &kcut);
224*b277271eSJames Wright   CeedScalar xdotd, vp[3] = {0.};
225*b277271eSJames Wright   CeedScalar xhat[] = {0., X[1], X[2]};
226*b277271eSJames Wright 
227*b277271eSJames Wright   CeedPragmaSIMD
228*b277271eSJames Wright   for(CeedInt n=0; n<nmodes; n++) {
229*b277271eSJames Wright     xhat[0] = (X[0] - stg_ctx->u0*t)*Max(2*kappa[0]/kappa[n], 0.1);
230*b277271eSJames Wright     xdotd = 0.;
231*b277271eSJames Wright     for(CeedInt i=0; i<3; i++) xdotd += d[i*nmodes+n]*xhat[i];
232*b277271eSJames Wright     const CeedScalar cos_kxdp = cos(kappa[n]*xdotd + phi[n]);
233*b277271eSJames Wright     const CeedScalar dkappa   = n==0 ? kappa[0] : kappa[n] - kappa[n-1];
234*b277271eSJames Wright     const CeedScalar qn       = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot);
235*b277271eSJames Wright     vp[0] += sqrt(qn)*sigma[0*nmodes+n] * cos_kxdp;
236*b277271eSJames Wright     vp[1] += sqrt(qn)*sigma[1*nmodes+n] * cos_kxdp;
237*b277271eSJames Wright     vp[2] += sqrt(qn)*sigma[2*nmodes+n] * cos_kxdp;
238*b277271eSJames Wright   }
239*b277271eSJames Wright   for(CeedInt i=0; i<3; i++) vp[i] *= 2*sqrt(1.5);
240*b277271eSJames Wright 
241*b277271eSJames Wright   u[0] = ubar[0] + cij[0]*vp[0];
242*b277271eSJames Wright   u[1] = ubar[1] + cij[3]*vp[0] + cij[1]*vp[1];
243*b277271eSJames Wright   u[2] = ubar[2] + cij[4]*vp[0] + cij[5]*vp[1] + cij[2]*vp[2];
244*b277271eSJames Wright }
245*b277271eSJames Wright 
246*b277271eSJames Wright CEED_QFUNCTION(Preprocess_STGShur14)(void *ctx, CeedInt Q,
247*b277271eSJames Wright                        const CeedScalar *const *in, CeedScalar *const *out) {
248*b277271eSJames Wright   //*INDENT-OFF*
249*b277271eSJames Wright   const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
250*b277271eSJames Wright                    (*x)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[1];
251*b277271eSJames Wright 
252*b277271eSJames Wright   CeedScalar (*stg_data) = (CeedScalar(*)) out[0];
253*b277271eSJames Wright 
254*b277271eSJames Wright   //*INDENT-ON*
255*b277271eSJames Wright   CeedScalar ubar[3], cij[6], eps, lt;
256*b277271eSJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
257*b277271eSJames Wright   const CeedScalar dx     = stg_ctx->dx;
258*b277271eSJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
259*b277271eSJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
260*b277271eSJames Wright   const CeedScalar P0     = stg_ctx->P0;
261*b277271eSJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
262*b277271eSJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
263*b277271eSJames Wright   const CeedScalar Rd     = cp - cv;
264*b277271eSJames Wright   const CeedScalar rho    = P0 / (Rd * theta0);
265*b277271eSJames Wright   const CeedScalar nu     = mu / rho;
266*b277271eSJames Wright 
267*b277271eSJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
268*b277271eSJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
269*b277271eSJames Wright   CeedScalar hmax, ke, keta, kcut, Ektot=0.0;
270*b277271eSJames Wright 
271*b277271eSJames Wright   CeedPragmaSIMD
272*b277271eSJames Wright   for(CeedInt i=0; i<Q; i++) {
273*b277271eSJames Wright     const CeedScalar dw = x[1][i];
274*b277271eSJames Wright     const CeedScalar dXdx[2][3] = {
275*b277271eSJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
276*b277271eSJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
277*b277271eSJames Wright     };
278*b277271eSJames Wright 
279*b277271eSJames Wright     CeedScalar h[3];
280*b277271eSJames Wright     h[0] = dx;
281*b277271eSJames Wright     for (CeedInt j=1; j<3; j++)
282*b277271eSJames Wright       h[j] = 2/sqrt(dXdx[0][j]*dXdx[0][j] + dXdx[1][j]*dXdx[1][j]);
283*b277271eSJames Wright 
284*b277271eSJames Wright     InterpolateProfile(dw, ubar, cij, &eps, &lt, stg_ctx);
285*b277271eSJames Wright     SpectrumConstants(dw, eps, lt, h, nu, &hmax, &ke, &keta, &kcut);
286*b277271eSJames Wright 
287*b277271eSJames Wright     // Calculate total TKE per spectrum
288*b277271eSJames Wright     stg_data[i] = 0.;
289*b277271eSJames Wright     CeedPragmaSIMD
290*b277271eSJames Wright     for(CeedInt n=0; n<nmodes; n++) {
291*b277271eSJames Wright       const CeedScalar dkappa = n==0 ? kappa[0] : kappa[n] - kappa[n-1];
292*b277271eSJames Wright       stg_data[i] += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0);
293*b277271eSJames Wright     }
294*b277271eSJames Wright   }
295*b277271eSJames Wright   return 0;
296*b277271eSJames Wright }
297*b277271eSJames Wright 
298b77c53c9SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition
299b77c53c9SJames Wright CEED_QFUNCTION(ICsSTG)(void *ctx, CeedInt Q,
300b77c53c9SJames Wright                        const CeedScalar *const *in, CeedScalar *const *out) {
301b77c53c9SJames Wright   // Inputs
302b77c53c9SJames Wright   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
303b77c53c9SJames Wright 
304b77c53c9SJames Wright   // Outputs
305b77c53c9SJames Wright   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
306b77c53c9SJames Wright 
307b77c53c9SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
308b77c53c9SJames Wright   CeedScalar u[3], cij[6], eps, lt;
309b77c53c9SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
310b77c53c9SJames Wright   const CeedScalar P0     = stg_ctx->P0;
311b77c53c9SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
312b77c53c9SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
313b77c53c9SJames Wright   const CeedScalar Rd     = cp - cv;
314b77c53c9SJames Wright   const CeedScalar rho    = P0 / (Rd * theta0);
315b77c53c9SJames Wright 
316b77c53c9SJames Wright   CeedPragmaSIMD
317b77c53c9SJames Wright   for(CeedInt i=0; i<Q; i++) {
318b77c53c9SJames Wright     InterpolateProfile(X[1][i], u, cij, &eps, &lt, stg_ctx);
319b77c53c9SJames Wright 
320b77c53c9SJames Wright     q0[0][i] = rho;
321b77c53c9SJames Wright     q0[1][i] = u[0] * rho;
322b77c53c9SJames Wright     q0[2][i] = u[1] * rho;
323b77c53c9SJames Wright     q0[3][i] = u[2] * rho;
324b77c53c9SJames Wright     q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0);
325b77c53c9SJames Wright   } // End of Quadrature Point Loop
326b77c53c9SJames Wright   return 0;
327b77c53c9SJames Wright }
328b77c53c9SJames Wright 
329ba6664aeSJames Wright /********************************************************************
330ba6664aeSJames Wright  * @brief QFunction to calculate the inflow boundary condition
331ba6664aeSJames Wright  *
332ba6664aeSJames Wright  * This will loop through quadrature points, calculate the wavemode amplitudes
333ba6664aeSJames Wright  * at each location, then calculate the actual velocity.
334ba6664aeSJames Wright  */
335ba6664aeSJames Wright CEED_QFUNCTION(STGShur14_Inflow)(void *ctx, CeedInt Q,
336ba6664aeSJames Wright                                  const CeedScalar *const *in,
337ba6664aeSJames Wright                                  CeedScalar *const *out) {
338ba6664aeSJames Wright 
339ba6664aeSJames Wright   //*INDENT-OFF*
340ba6664aeSJames Wright   const CeedScalar (*q)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
341e8b03feeSJames Wright                    (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[2],
342e8b03feeSJames Wright                    (*X)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[3];
343ba6664aeSJames Wright 
3444dbab5e5SJames Wright   CeedScalar(*v)[CEED_Q_VLA]            = (CeedScalar(*)[CEED_Q_VLA]) out[0],
3454dbab5e5SJames Wright             (*jac_data_sur)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[1];
346ba6664aeSJames Wright 
347ba6664aeSJames Wright   //*INDENT-ON*
348ba6664aeSJames Wright 
349ba6664aeSJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
350ba6664aeSJames Wright   CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
351ba6664aeSJames Wright   const bool is_implicit  = stg_ctx->is_implicit;
352ba6664aeSJames Wright   const bool mean_only    = stg_ctx->mean_only;
353ba6664aeSJames Wright   const bool prescribe_T  = stg_ctx->prescribe_T;
354ba6664aeSJames Wright   const CeedScalar dx     = stg_ctx->dx;
355ba6664aeSJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
356ba6664aeSJames Wright   const CeedScalar time   = stg_ctx->time;
357ba6664aeSJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
358ba6664aeSJames Wright   const CeedScalar P0     = stg_ctx->P0;
359ba6664aeSJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
360ba6664aeSJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
361ba6664aeSJames Wright   const CeedScalar Rd     = cp - cv;
362ba6664aeSJames Wright   const CeedScalar gamma  = cp/cv;
363ba6664aeSJames Wright 
364ba6664aeSJames Wright   CeedPragmaSIMD
365ba6664aeSJames Wright   for(CeedInt i=0; i<Q; i++) {
366ba6664aeSJames Wright     const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0);
367ba6664aeSJames Wright     const CeedScalar x[] = { X[0][i], X[1][i], X[2][i] };
368ba6664aeSJames Wright     const CeedScalar dXdx[2][3] = {
369ba6664aeSJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
370ba6664aeSJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
371ba6664aeSJames Wright     };
372ba6664aeSJames Wright 
373ba6664aeSJames Wright     CeedScalar h[3];
374ba6664aeSJames Wright     h[0] = dx;
375a939fbabSJames Wright     for (CeedInt j=1; j<3; j++)
376a939fbabSJames Wright       h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]));
377ba6664aeSJames Wright 
378ba6664aeSJames Wright     InterpolateProfile(X[1][i], ubar, cij, &eps, &lt, stg_ctx);
379ba6664aeSJames Wright     if (!mean_only) {
380ba6664aeSJames Wright       CalcSpectrum(X[1][i], eps, lt, h, mu/rho, qn, stg_ctx);
381ba6664aeSJames Wright       STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
382ba6664aeSJames Wright     } else {
383ba6664aeSJames Wright       for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
384ba6664aeSJames Wright     }
385ba6664aeSJames Wright 
3864dbab5e5SJames Wright     const CeedScalar E_kinetic = .5 * rho * Dot3(u, u);
387ba6664aeSJames Wright     CeedScalar E_internal, P;
388ba6664aeSJames Wright     if (prescribe_T) {
389ba6664aeSJames Wright       // Temperature is being set weakly (theta0) and for constant cv this sets E_internal
390ba6664aeSJames Wright       E_internal = rho * cv * theta0;
391ba6664aeSJames Wright       // Find pressure using
392ba6664aeSJames Wright       P = rho * Rd * theta0; // interior rho with exterior T
393ba6664aeSJames Wright     } else {
394ba6664aeSJames Wright       E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution
395ba6664aeSJames Wright       P = E_internal * (gamma - 1.);
396ba6664aeSJames Wright     }
397ba6664aeSJames Wright 
398ba6664aeSJames Wright     const CeedScalar wdetJb  = (is_implicit ? -1. : 1.) * q_data_sur[0][i];
399ba6664aeSJames Wright     // ---- Normal vect
400ba6664aeSJames Wright     const CeedScalar norm[3] = {q_data_sur[1][i],
401ba6664aeSJames Wright                                 q_data_sur[2][i],
402ba6664aeSJames Wright                                 q_data_sur[3][i]
403ba6664aeSJames Wright                                };
404ba6664aeSJames Wright 
405ba6664aeSJames Wright     const CeedScalar E = E_internal + E_kinetic;
406ba6664aeSJames Wright 
407ba6664aeSJames Wright     // Velocity normal to the boundary
4084dbab5e5SJames Wright     const CeedScalar u_normal = Dot3(norm, u);
4094dbab5e5SJames Wright 
410ba6664aeSJames Wright     // The Physics
411ba6664aeSJames Wright     // Zero v so all future terms can safely sum into it
412ba6664aeSJames Wright     for (CeedInt j=0; j<5; j++) v[j][i] = 0.;
413ba6664aeSJames Wright 
414ba6664aeSJames Wright     // The Physics
415ba6664aeSJames Wright     // -- Density
416ba6664aeSJames Wright     v[0][i] -= wdetJb * rho * u_normal;
417ba6664aeSJames Wright 
418ba6664aeSJames Wright     // -- Momentum
419ba6664aeSJames Wright     for (CeedInt j=0; j<3; j++)
420ba6664aeSJames Wright       v[j+1][i] -= wdetJb *(rho * u_normal * u[j] +
421ba6664aeSJames Wright                             norm[j] * P);
422ba6664aeSJames Wright 
423ba6664aeSJames Wright     // -- Total Energy Density
424ba6664aeSJames Wright     v[4][i] -= wdetJb * u_normal * (E + P);
4254dbab5e5SJames Wright 
4264dbab5e5SJames Wright     jac_data_sur[0][i] = rho;
4274dbab5e5SJames Wright     jac_data_sur[1][i] = u[0];
4284dbab5e5SJames Wright     jac_data_sur[2][i] = u[1];
4294dbab5e5SJames Wright     jac_data_sur[3][i] = u[2];
4304dbab5e5SJames Wright     jac_data_sur[4][i] = E;
4314dbab5e5SJames Wright     for (int j=0; j<6; j++) jac_data_sur[5+j][i] = 0.;
432ba6664aeSJames Wright   }
433ba6664aeSJames Wright   return 0;
434ba6664aeSJames Wright }
435ba6664aeSJames Wright 
4364dbab5e5SJames Wright CEED_QFUNCTION(STGShur14_Inflow_Jacobian)(void *ctx, CeedInt Q,
4374dbab5e5SJames Wright     const CeedScalar *const *in,
4384dbab5e5SJames Wright     CeedScalar *const *out) {
4394dbab5e5SJames Wright   // *INDENT-OFF*
4404dbab5e5SJames Wright   // Inputs
4414dbab5e5SJames Wright   const CeedScalar (*dq)[CEED_Q_VLA]           = (const CeedScalar(*)[CEED_Q_VLA])in[0],
4424dbab5e5SJames Wright                    (*q_data_sur)[CEED_Q_VLA]   = (const CeedScalar(*)[CEED_Q_VLA])in[2],
4434dbab5e5SJames Wright                    (*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4];
4444dbab5e5SJames Wright   // Outputs
4454dbab5e5SJames Wright   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
4464dbab5e5SJames Wright   // *INDENT-ON*
4474dbab5e5SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
4484dbab5e5SJames Wright   const bool implicit     = stg_ctx->is_implicit;
4494dbab5e5SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
4504dbab5e5SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
4514dbab5e5SJames Wright   const CeedScalar Rd     = cp - cv;
4524dbab5e5SJames Wright   const CeedScalar gamma  = cp/cv;
4534dbab5e5SJames Wright 
4544dbab5e5SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
4554dbab5e5SJames Wright   const bool prescribe_T  = stg_ctx->prescribe_T;
4564dbab5e5SJames Wright 
4574dbab5e5SJames Wright   CeedPragmaSIMD
4584dbab5e5SJames Wright   // Quadrature Point Loop
4594dbab5e5SJames Wright   for (CeedInt i=0; i<Q; i++) {
4604dbab5e5SJames Wright     // Setup
4614dbab5e5SJames Wright     // Setup
4624dbab5e5SJames Wright     // -- Interp-to-Interp q_data
4634dbab5e5SJames Wright     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
4644dbab5e5SJames Wright     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
4654dbab5e5SJames Wright     // We can effect this by swapping the sign on this weight
4664dbab5e5SJames Wright     const CeedScalar wdetJb  = (implicit ? -1. : 1.) * q_data_sur[0][i];
4674dbab5e5SJames Wright 
4684dbab5e5SJames Wright     // Calculate inflow values
4694dbab5e5SJames Wright     CeedScalar velocity[3];
4704dbab5e5SJames Wright     for (CeedInt j=0; j<3; j++) velocity[j] = jac_data_sur[5+j][i];
4714dbab5e5SJames Wright 
4724dbab5e5SJames Wright     // enabling user to choose between weak T and weak rho inflow
4734dbab5e5SJames Wright     CeedScalar drho, dE, dP;
4744dbab5e5SJames Wright     if (prescribe_T) {
4754dbab5e5SJames Wright       // rho should be from the current solution
4764dbab5e5SJames Wright       drho = dq[0][i];
4774dbab5e5SJames Wright       CeedScalar dE_internal = drho * cv * theta0;
4784dbab5e5SJames Wright       CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity);
4794dbab5e5SJames Wright       dE = dE_internal + dE_kinetic;
4804dbab5e5SJames Wright       dP = drho * Rd * theta0; // interior rho with exterior T
4814dbab5e5SJames Wright     } else { // rho specified, E_internal from solution
4824dbab5e5SJames Wright       drho = 0;
4834dbab5e5SJames Wright       dE = dq[4][i];
4844dbab5e5SJames Wright       dP = dE * (gamma - 1.);
4854dbab5e5SJames Wright     }
4864dbab5e5SJames Wright     const CeedScalar norm[3] = {q_data_sur[1][i],
4874dbab5e5SJames Wright                                 q_data_sur[2][i],
4884dbab5e5SJames Wright                                 q_data_sur[3][i]
4894dbab5e5SJames Wright                                };
4904dbab5e5SJames Wright 
4914dbab5e5SJames Wright     const CeedScalar u_normal = Dot3(norm, velocity);
4924dbab5e5SJames Wright 
4934dbab5e5SJames Wright     v[0][i] = - wdetJb * drho * u_normal;
4944dbab5e5SJames Wright     for (int j=0; j<3; j++)
4954dbab5e5SJames Wright       v[j+1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP);
4964dbab5e5SJames Wright     v[4][i] = - wdetJb * u_normal * (dE + dP);
4974dbab5e5SJames Wright   } // End Quadrature Point Loop
4984dbab5e5SJames Wright   return 0;
4994dbab5e5SJames Wright }
5004dbab5e5SJames Wright 
5010a6353c2SJames Wright /********************************************************************
5020a6353c2SJames Wright  * @brief QFunction to calculate the strongly enforce inflow BC
5030a6353c2SJames Wright  *
5040a6353c2SJames Wright  * This QF is for the strong application of STG via libCEED (rather than
5050a6353c2SJames Wright  * through the native PETSc `DMAddBoundary` -> `bcFunc` method.
5060a6353c2SJames Wright  */
5070a6353c2SJames Wright CEED_QFUNCTION(STGShur14_Inflow_StrongQF)(void *ctx, CeedInt Q,
5080a6353c2SJames Wright     const CeedScalar *const *in, CeedScalar *const *out) {
5090a6353c2SJames Wright 
5100a6353c2SJames Wright   //*INDENT-OFF*
5110a6353c2SJames Wright   const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
5120a6353c2SJames Wright                    (*coords)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA]) in[1],
5130a6353c2SJames Wright                    (*scale)                  = (const CeedScalar(*)) in[2];
5140a6353c2SJames Wright 
5150a6353c2SJames Wright   CeedScalar(*bcval)[CEED_Q_VLA]            = (CeedScalar(*)[CEED_Q_VLA]) out[0];
5160a6353c2SJames Wright   //*INDENT-ON*
5170a6353c2SJames Wright 
5180a6353c2SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
5190a6353c2SJames Wright   CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
5200a6353c2SJames Wright   const bool mean_only    = stg_ctx->mean_only;
5210a6353c2SJames Wright   const CeedScalar dx     = stg_ctx->dx;
5220a6353c2SJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
5230a6353c2SJames Wright   const CeedScalar time   = stg_ctx->time;
5240a6353c2SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
5250a6353c2SJames Wright   const CeedScalar P0     = stg_ctx->P0;
5260a6353c2SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
5270a6353c2SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
5280a6353c2SJames Wright   const CeedScalar Rd     = cp - cv;
5290a6353c2SJames Wright   const CeedScalar rho    = P0 / (Rd * theta0);
5300a6353c2SJames Wright 
5310a6353c2SJames Wright   CeedPragmaSIMD
5320a6353c2SJames Wright   for(CeedInt i=0; i<Q; i++) {
5330a6353c2SJames Wright     const CeedScalar x[] = { coords[0][i], coords[1][i], coords[2][i] };
5340a6353c2SJames Wright     const CeedScalar dXdx[2][3] = {
5350a6353c2SJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
5360a6353c2SJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
5370a6353c2SJames Wright     };
5380a6353c2SJames Wright 
5390a6353c2SJames Wright     CeedScalar h[3];
5400a6353c2SJames Wright     h[0] = dx;
541a939fbabSJames Wright     for (CeedInt j=1; j<3; j++)
542a939fbabSJames Wright       h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]));
5430a6353c2SJames Wright 
5440a6353c2SJames Wright     InterpolateProfile(coords[1][i], ubar, cij, &eps, &lt, stg_ctx);
5450a6353c2SJames Wright     if (!mean_only) {
5460a6353c2SJames Wright       CalcSpectrum(coords[1][i], eps, lt, h, mu/rho, qn, stg_ctx);
5470a6353c2SJames Wright       STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
5480a6353c2SJames Wright     } else {
5490a6353c2SJames Wright       for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
5500a6353c2SJames Wright     }
5510a6353c2SJames Wright 
5520a6353c2SJames Wright     bcval[0][i] = scale[i] * rho;
5530a6353c2SJames Wright     bcval[1][i] = scale[i] * rho * u[0];
5540a6353c2SJames Wright     bcval[2][i] = scale[i] * rho * u[1];
5550a6353c2SJames Wright     bcval[3][i] = scale[i] * rho * u[2];
556cf3d54ffSJames Wright     bcval[4][i] = 0.;
5570a6353c2SJames Wright   }
5580a6353c2SJames Wright   return 0;
5590a6353c2SJames Wright }
5600a6353c2SJames Wright 
561ba6664aeSJames Wright #endif // stg_shur14_h
562