xref: /honee/qfunctions/stg_shur14.h (revision 9eeef72b195e5d6f3ec9a875fa44c8a602a1f5fd)
1493642f1SJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2493642f1SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3493642f1SJames Wright //
4493642f1SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5493642f1SJames Wright //
6493642f1SJames Wright // This file is part of CEED:  http://github.com/ceed
7493642f1SJames Wright 
8493642f1SJames Wright /// @file
9493642f1SJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm
10493642f1SJames Wright /// presented in Shur et al. 2014
11493642f1SJames Wright //
12493642f1SJames Wright /// SetupSTG_Rand reads in the input files and fills in STGShur14Context. Then
13493642f1SJames Wright /// STGShur14_CalcQF is run over quadrature points. Before the program exits,
14493642f1SJames Wright /// TearDownSTG is run to free the memory of the allocated arrays.
15493642f1SJames Wright 
16493642f1SJames Wright #ifndef stg_shur14_h
17493642f1SJames Wright #define stg_shur14_h
18493642f1SJames Wright 
19493642f1SJames Wright #include <math.h>
20493642f1SJames Wright #include <ceed.h>
21493642f1SJames Wright #include <stdlib.h>
22493642f1SJames Wright #include "stg_shur14_type.h"
23704b8bbeSJames Wright #include "utils.h"
24493642f1SJames Wright 
25493642f1SJames Wright #define STG_NMODES_MAX 1024
26493642f1SJames Wright 
27493642f1SJames Wright /*
28493642f1SJames Wright  * @brief Interpolate quantities from input profile to given location
29493642f1SJames Wright  *
30493642f1SJames Wright  * Assumed that prof_dw[i+1] > prof_dw[i] and prof_dw[0] = 0
31493642f1SJames Wright  * If dw > prof_dw[-1], then the interpolation takes the values at prof_dw[-1]
32493642f1SJames Wright  *
33493642f1SJames Wright  * @param[in]  dw      Distance to the nearest wall
34493642f1SJames Wright  * @param[out] ubar    Mean velocity at dw
35493642f1SJames Wright  * @param[out] cij     Cholesky decomposition at dw
36493642f1SJames Wright  * @param[out] eps     Turbulent dissipation at dw
37493642f1SJames Wright  * @param[out] lt      Turbulent length scale at dw
38493642f1SJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
39493642f1SJames Wright  */
40493642f1SJames Wright CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar dw,
41493642f1SJames Wright     CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt,
42493642f1SJames Wright     const STGShur14Context stg_ctx) {
43493642f1SJames Wright 
44493642f1SJames Wright   const CeedInt    nprofs    = stg_ctx->nprofs;
45493642f1SJames Wright   const CeedScalar *prof_dw  = &stg_ctx->data[stg_ctx->offsets.prof_dw];
46493642f1SJames Wright   const CeedScalar *prof_eps = &stg_ctx->data[stg_ctx->offsets.eps];
47493642f1SJames Wright   const CeedScalar *prof_lt  = &stg_ctx->data[stg_ctx->offsets.lt];
48493642f1SJames Wright   const CeedScalar *prof_ubar = &stg_ctx->data[stg_ctx->offsets.ubar];
49493642f1SJames Wright   const CeedScalar *prof_cij  = &stg_ctx->data[stg_ctx->offsets.cij];
50493642f1SJames Wright   CeedInt idx=-1;
51493642f1SJames Wright 
52493642f1SJames Wright   for(CeedInt i=0; i<nprofs; i++) {
53493642f1SJames Wright     if (dw < prof_dw[i]) {
54493642f1SJames Wright       idx = i;
55493642f1SJames Wright       break;
56493642f1SJames Wright     }
57493642f1SJames Wright   }
58493642f1SJames Wright 
59493642f1SJames Wright   if (idx > 0) { // y within the bounds of prof_dw
60493642f1SJames Wright     CeedScalar coeff = (dw - prof_dw[idx-1]) / (prof_dw[idx] - prof_dw[idx-1]);
61493642f1SJames Wright 
62493642f1SJames Wright     //*INDENT-OFF*
63493642f1SJames Wright     ubar[0] = prof_ubar[0*nprofs+idx-1] + coeff*( prof_ubar[0*nprofs+idx] - prof_ubar[0*nprofs+idx-1] );
64493642f1SJames Wright     ubar[1] = prof_ubar[1*nprofs+idx-1] + coeff*( prof_ubar[1*nprofs+idx] - prof_ubar[1*nprofs+idx-1] );
65493642f1SJames Wright     ubar[2] = prof_ubar[2*nprofs+idx-1] + coeff*( prof_ubar[2*nprofs+idx] - prof_ubar[2*nprofs+idx-1] );
66493642f1SJames Wright     cij[0]  = prof_cij[0*nprofs+idx-1]  + coeff*( prof_cij[0*nprofs+idx]  - prof_cij[0*nprofs+idx-1] );
67493642f1SJames Wright     cij[1]  = prof_cij[1*nprofs+idx-1]  + coeff*( prof_cij[1*nprofs+idx]  - prof_cij[1*nprofs+idx-1] );
68493642f1SJames Wright     cij[2]  = prof_cij[2*nprofs+idx-1]  + coeff*( prof_cij[2*nprofs+idx]  - prof_cij[2*nprofs+idx-1] );
69493642f1SJames Wright     cij[3]  = prof_cij[3*nprofs+idx-1]  + coeff*( prof_cij[3*nprofs+idx]  - prof_cij[3*nprofs+idx-1] );
70493642f1SJames Wright     cij[4]  = prof_cij[4*nprofs+idx-1]  + coeff*( prof_cij[4*nprofs+idx]  - prof_cij[4*nprofs+idx-1] );
71493642f1SJames Wright     cij[5]  = prof_cij[5*nprofs+idx-1]  + coeff*( prof_cij[5*nprofs+idx]  - prof_cij[5*nprofs+idx-1] );
72493642f1SJames Wright     *eps    = prof_eps[idx-1]           + coeff*( prof_eps[idx]           - prof_eps[idx-1] );
73493642f1SJames Wright     *lt     = prof_lt[idx-1]            + coeff*( prof_lt[idx]            - prof_lt[idx-1] );
74493642f1SJames Wright     //*INDENT-ON*
75493642f1SJames Wright   } else { // y outside bounds of prof_dw
76493642f1SJames Wright     ubar[0] = prof_ubar[1*nprofs-1];
77493642f1SJames Wright     ubar[1] = prof_ubar[2*nprofs-1];
78493642f1SJames Wright     ubar[2] = prof_ubar[3*nprofs-1];
79493642f1SJames Wright     cij[0]  = prof_cij[1*nprofs-1];
80493642f1SJames Wright     cij[1]  = prof_cij[2*nprofs-1];
81493642f1SJames Wright     cij[2]  = prof_cij[3*nprofs-1];
82493642f1SJames Wright     cij[3]  = prof_cij[4*nprofs-1];
83493642f1SJames Wright     cij[4]  = prof_cij[5*nprofs-1];
84493642f1SJames Wright     cij[5]  = prof_cij[6*nprofs-1];
85493642f1SJames Wright     *eps    = prof_eps[nprofs-1];
86493642f1SJames Wright     *lt     = prof_lt[nprofs-1];
87493642f1SJames Wright   }
88493642f1SJames Wright }
89493642f1SJames Wright 
90493642f1SJames Wright /*
9171cd6200SJames Wright  * @brief Calculate spectrum coefficient, qn
9271cd6200SJames Wright  *
9371cd6200SJames Wright  * Calculates q_n at a given distance to the wall
9471cd6200SJames Wright  *
9571cd6200SJames Wright  * @param[in]  kappa  nth wavenumber
9671cd6200SJames Wright  * @param[in]  dkappa Difference between wavenumbers
9771cd6200SJames Wright  * @param[in]  keta   Dissipation wavenumber
9871cd6200SJames Wright  * @param[in]  kcut   Mesh-induced cutoff wavenumber
9971cd6200SJames Wright  * @param[in]  ke     Energy-containing wavenumber
10071cd6200SJames Wright  * @param[in]  Ektot  Total turbulent kinetic energy of spectrum
10171cd6200SJames Wright  * @returns    qn     Spectrum coefficient
10271cd6200SJames Wright  */
10371cd6200SJames Wright CeedScalar CEED_QFUNCTION_HELPER(Calc_qn)(const CeedScalar kappa,
10471cd6200SJames Wright     const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut,
10571cd6200SJames Wright     const CeedScalar ke, const CeedScalar Ektot) {
10671cd6200SJames Wright   const CeedScalar feta_x_fcut   = exp(-Square(12*kappa/keta)
10771cd6200SJames Wright                                        -Cube(4*Max(kappa - 0.9*kcut, 0)/kcut) );
10871cd6200SJames Wright   return pow(kappa/ke, 4.) * pow(1 + 2.4*Square(kappa/ke),-17./6)
10971cd6200SJames Wright          *feta_x_fcut*dkappa/Ektot;
11071cd6200SJames Wright }
11171cd6200SJames Wright 
11271cd6200SJames Wright // Calculate hmax, ke, keta, and kcut
11371cd6200SJames Wright void CEED_QFUNCTION_HELPER(SpectrumConstants)(const CeedScalar dw,
11471cd6200SJames Wright     const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3],
11571cd6200SJames Wright     const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke,
11671cd6200SJames Wright     CeedScalar *keta, CeedScalar *kcut) {
11771cd6200SJames Wright   *hmax = Max( Max(h[0], h[1]), h[2]);
11871cd6200SJames Wright   *ke   = dw==0 ? 1e16 : 2*M_PI/Min(2*dw, 3*lt);
11971cd6200SJames Wright   *keta = 2*M_PI*pow(Cube(nu)/eps, -0.25);
12071cd6200SJames Wright   *kcut = M_PI/ Min( Max(Max(h[1], h[2]), 0.3*(*hmax)) + 0.1*dw, *hmax );
12171cd6200SJames Wright }
12271cd6200SJames Wright 
12371cd6200SJames Wright /*
124493642f1SJames Wright  * @brief Calculate spectrum coefficients for STG
125493642f1SJames Wright  *
126493642f1SJames Wright  * Calculates q_n at a given distance to the wall
127493642f1SJames Wright  *
128493642f1SJames Wright  * @param[in]  dw      Distance to the nearest wall
129493642f1SJames Wright  * @param[in]  eps     Turbulent dissipation w/rt dw
130493642f1SJames Wright  * @param[in]  lt      Turbulent length scale w/rt dw
131493642f1SJames Wright  * @param[in]  h       Element lengths in coordinate directions
132493642f1SJames Wright  * @param[in]  nu      Dynamic Viscosity;
133493642f1SJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
134493642f1SJames Wright  * @param[out] qn      Spectrum coefficients, [nmodes]
135493642f1SJames Wright  */
136493642f1SJames Wright void CEED_QFUNCTION_HELPER(CalcSpectrum)(const CeedScalar dw,
137493642f1SJames Wright     const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3],
138493642f1SJames Wright     const CeedScalar nu, CeedScalar qn[], const STGShur14Context stg_ctx) {
139493642f1SJames Wright 
140493642f1SJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
141493642f1SJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
14271cd6200SJames Wright   CeedScalar hmax, ke, keta, kcut, Ektot=0.0;
14371cd6200SJames Wright   SpectrumConstants(dw, eps, lt, h, nu, &hmax, &ke, &keta, &kcut);
144493642f1SJames Wright 
145493642f1SJames Wright   for(CeedInt n=0; n<nmodes; n++) {
14671cd6200SJames Wright     const CeedScalar dkappa = n==0 ? kappa[0] : kappa[n] - kappa[n-1];
14771cd6200SJames Wright     qn[n] = Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0);
148493642f1SJames Wright     Ektot += qn[n];
149493642f1SJames Wright   }
150493642f1SJames Wright 
1510a8dc919SJames Wright   if (Ektot == 0) return;
152493642f1SJames Wright   for(CeedInt n=0; n<nmodes; n++) qn[n] /= Ektot;
153493642f1SJames Wright }
154493642f1SJames Wright 
155493642f1SJames Wright /******************************************************
156493642f1SJames Wright  * @brief Calculate u(x,t) for STG inflow condition
157493642f1SJames Wright  *
158493642f1SJames Wright  * @param[in]  X       Location to evaluate u(X,t)
159493642f1SJames Wright  * @param[in]  t       Time to evaluate u(X,t)
160493642f1SJames Wright  * @param[in]  ubar    Mean velocity at X
161493642f1SJames Wright  * @param[in]  cij     Cholesky decomposition at X
162493642f1SJames Wright  * @param[in]  qn      Wavemode amplitudes at X, [nmodes]
163493642f1SJames Wright  * @param[out] u       Velocity at X and t
164493642f1SJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
165493642f1SJames Wright  */
166493642f1SJames Wright void CEED_QFUNCTION_HELPER(STGShur14_Calc)(const CeedScalar X[3],
167493642f1SJames Wright     const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6],
168493642f1SJames Wright     const CeedScalar qn[], CeedScalar u[3],
169493642f1SJames Wright     const STGShur14Context stg_ctx) {
170493642f1SJames Wright 
171493642f1SJames Wright   //*INDENT-OFF*
172493642f1SJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
173493642f1SJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
174493642f1SJames Wright   const CeedScalar *phi   = &stg_ctx->data[stg_ctx->offsets.phi];
175493642f1SJames Wright   const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma];
176493642f1SJames Wright   const CeedScalar *d     = &stg_ctx->data[stg_ctx->offsets.d];
177493642f1SJames Wright   //*INDENT-ON*
178493642f1SJames Wright   CeedScalar xdotd, vp[3] = {0.};
179493642f1SJames Wright   CeedScalar xhat[] = {0., X[1], X[2]};
180493642f1SJames Wright 
181493642f1SJames Wright   CeedPragmaSIMD
182493642f1SJames Wright   for(CeedInt n=0; n<nmodes; n++) {
183493642f1SJames Wright     xhat[0] = (X[0] - stg_ctx->u0*t)*Max(2*kappa[0]/kappa[n], 0.1);
184493642f1SJames Wright     xdotd = 0.;
185493642f1SJames Wright     for(CeedInt i=0; i<3; i++) xdotd += d[i*nmodes+n]*xhat[i];
186493642f1SJames Wright     const CeedScalar cos_kxdp = cos(kappa[n]*xdotd + phi[n]);
1870a8dc919SJames Wright     vp[0] += sqrt(qn[n])*sigma[0*nmodes+n] * cos_kxdp;
1880a8dc919SJames Wright     vp[1] += sqrt(qn[n])*sigma[1*nmodes+n] * cos_kxdp;
1890a8dc919SJames Wright     vp[2] += sqrt(qn[n])*sigma[2*nmodes+n] * cos_kxdp;
190493642f1SJames Wright   }
1910a8dc919SJames Wright   for(CeedInt i=0; i<3; i++) vp[i] *= 2*sqrt(1.5);
192493642f1SJames Wright 
193493642f1SJames Wright   u[0] = ubar[0] + cij[0]*vp[0];
194493642f1SJames Wright   u[1] = ubar[1] + cij[3]*vp[0] + cij[1]*vp[1];
195493642f1SJames Wright   u[2] = ubar[2] + cij[4]*vp[0] + cij[5]*vp[1] + cij[2]*vp[2];
196493642f1SJames Wright }
197493642f1SJames Wright 
1988eea80fcSJames Wright /******************************************************
1998eea80fcSJames Wright  * @brief Calculate u(x,t) for STG inflow condition
2008eea80fcSJames Wright  *
2018eea80fcSJames Wright  * @param[in]  X       Location to evaluate u(X,t)
2028eea80fcSJames Wright  * @param[in]  t       Time to evaluate u(X,t)
2038eea80fcSJames Wright  * @param[in]  ubar    Mean velocity at X
2048eea80fcSJames Wright  * @param[in]  cij     Cholesky decomposition at X
2058eea80fcSJames Wright  * @param[in]  qn      Wavemode amplitudes at X, [nmodes]
2068eea80fcSJames Wright  * @param[out] u       Velocity at X and t
2078eea80fcSJames Wright  * @param[in]  stg_ctx STGShur14Context for the problem
2088eea80fcSJames Wright  */
2098eea80fcSJames Wright void CEED_QFUNCTION_HELPER(STGShur14_Calc_PrecompEktot)(const CeedScalar X[3],
2108eea80fcSJames Wright     const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6],
2118eea80fcSJames Wright     const CeedScalar Ektot, const CeedScalar h[3], const CeedScalar dw,
2128eea80fcSJames Wright     const CeedScalar eps, const CeedScalar lt, const CeedScalar nu, CeedScalar u[3],
2138eea80fcSJames Wright     const STGShur14Context stg_ctx) {
2148eea80fcSJames Wright 
2158eea80fcSJames Wright   //*INDENT-OFF*
2168eea80fcSJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
2178eea80fcSJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
2188eea80fcSJames Wright   const CeedScalar *phi   = &stg_ctx->data[stg_ctx->offsets.phi];
2198eea80fcSJames Wright   const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma];
2208eea80fcSJames Wright   const CeedScalar *d     = &stg_ctx->data[stg_ctx->offsets.d];
2218eea80fcSJames Wright   //*INDENT-ON*
2228eea80fcSJames Wright   CeedScalar hmax, ke, keta, kcut;
2238eea80fcSJames Wright   SpectrumConstants(dw, eps, lt, h, nu, &hmax, &ke, &keta, &kcut);
2248eea80fcSJames Wright   CeedScalar xdotd, vp[3] = {0.};
2258eea80fcSJames Wright   CeedScalar xhat[] = {0., X[1], X[2]};
2268eea80fcSJames Wright 
2278eea80fcSJames Wright   CeedPragmaSIMD
2288eea80fcSJames Wright   for(CeedInt n=0; n<nmodes; n++) {
2298eea80fcSJames Wright     xhat[0] = (X[0] - stg_ctx->u0*t)*Max(2*kappa[0]/kappa[n], 0.1);
2308eea80fcSJames Wright     xdotd = 0.;
2318eea80fcSJames Wright     for(CeedInt i=0; i<3; i++) xdotd += d[i*nmodes+n]*xhat[i];
2328eea80fcSJames Wright     const CeedScalar cos_kxdp = cos(kappa[n]*xdotd + phi[n]);
2338eea80fcSJames Wright     const CeedScalar dkappa   = n==0 ? kappa[0] : kappa[n] - kappa[n-1];
2348eea80fcSJames Wright     const CeedScalar qn       = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot);
2358eea80fcSJames Wright     vp[0] += sqrt(qn)*sigma[0*nmodes+n] * cos_kxdp;
2368eea80fcSJames Wright     vp[1] += sqrt(qn)*sigma[1*nmodes+n] * cos_kxdp;
2378eea80fcSJames Wright     vp[2] += sqrt(qn)*sigma[2*nmodes+n] * cos_kxdp;
2388eea80fcSJames Wright   }
2398eea80fcSJames Wright   for(CeedInt i=0; i<3; i++) vp[i] *= 2*sqrt(1.5);
2408eea80fcSJames Wright 
2418eea80fcSJames Wright   u[0] = ubar[0] + cij[0]*vp[0];
2428eea80fcSJames Wright   u[1] = ubar[1] + cij[3]*vp[0] + cij[1]*vp[1];
2438eea80fcSJames Wright   u[2] = ubar[2] + cij[4]*vp[0] + cij[5]*vp[1] + cij[2]*vp[2];
2448eea80fcSJames Wright }
2458eea80fcSJames Wright 
2468eea80fcSJames Wright CEED_QFUNCTION(Preprocess_STGShur14)(void *ctx, CeedInt Q,
2478eea80fcSJames Wright                                      const CeedScalar *const *in, CeedScalar *const *out) {
2488eea80fcSJames Wright   //*INDENT-OFF*
2498eea80fcSJames Wright   const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
2508eea80fcSJames Wright                    (*x)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[1];
2518eea80fcSJames Wright 
2528eea80fcSJames Wright   CeedScalar (*stg_data) = (CeedScalar(*)) out[0];
2538eea80fcSJames Wright 
2548eea80fcSJames Wright   //*INDENT-ON*
2558eea80fcSJames Wright   CeedScalar ubar[3], cij[6], eps, lt;
2568eea80fcSJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
2578eea80fcSJames Wright   const CeedScalar dx     = stg_ctx->dx;
2588eea80fcSJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
2598eea80fcSJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
2608eea80fcSJames Wright   const CeedScalar P0     = stg_ctx->P0;
2618eea80fcSJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
2628eea80fcSJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
2638eea80fcSJames Wright   const CeedScalar Rd     = cp - cv;
2648eea80fcSJames Wright   const CeedScalar rho    = P0 / (Rd * theta0);
2658eea80fcSJames Wright   const CeedScalar nu     = mu / rho;
2668eea80fcSJames Wright 
2678eea80fcSJames Wright   const CeedInt    nmodes = stg_ctx->nmodes;
2688eea80fcSJames Wright   const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
269*9eeef72bSJames Wright   CeedScalar hmax, ke, keta, kcut;
2708eea80fcSJames Wright 
2718eea80fcSJames Wright   CeedPragmaSIMD
2728eea80fcSJames Wright   for(CeedInt i=0; i<Q; i++) {
2738eea80fcSJames Wright     const CeedScalar dw = x[1][i];
2748eea80fcSJames Wright     const CeedScalar dXdx[2][3] = {
2758eea80fcSJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
2768eea80fcSJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
2778eea80fcSJames Wright     };
2788eea80fcSJames Wright 
2798eea80fcSJames Wright     CeedScalar h[3];
2808eea80fcSJames Wright     h[0] = dx;
2818eea80fcSJames Wright     for (CeedInt j=1; j<3; j++)
2828eea80fcSJames Wright       h[j] = 2/sqrt(dXdx[0][j]*dXdx[0][j] + dXdx[1][j]*dXdx[1][j]);
2838eea80fcSJames Wright 
2848eea80fcSJames Wright     InterpolateProfile(dw, ubar, cij, &eps, &lt, stg_ctx);
2858eea80fcSJames Wright     SpectrumConstants(dw, eps, lt, h, nu, &hmax, &ke, &keta, &kcut);
2868eea80fcSJames Wright 
2878eea80fcSJames Wright     // Calculate total TKE per spectrum
2888eea80fcSJames Wright     stg_data[i] = 0.;
2898eea80fcSJames Wright     CeedPragmaSIMD
2908eea80fcSJames Wright     for(CeedInt n=0; n<nmodes; n++) {
2918eea80fcSJames Wright       const CeedScalar dkappa = n==0 ? kappa[0] : kappa[n] - kappa[n-1];
2928eea80fcSJames Wright       stg_data[i] += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0);
2938eea80fcSJames Wright     }
2948eea80fcSJames Wright   }
2958eea80fcSJames Wright   return 0;
2968eea80fcSJames Wright }
2978eea80fcSJames Wright 
29843bff553SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition
29943bff553SJames Wright CEED_QFUNCTION(ICsSTG)(void *ctx, CeedInt Q,
30043bff553SJames Wright                        const CeedScalar *const *in, CeedScalar *const *out) {
30143bff553SJames Wright   // Inputs
30243bff553SJames Wright   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
30343bff553SJames Wright 
30443bff553SJames Wright   // Outputs
30543bff553SJames Wright   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
30643bff553SJames Wright 
30743bff553SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
30843bff553SJames Wright   CeedScalar u[3], cij[6], eps, lt;
30943bff553SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
31043bff553SJames Wright   const CeedScalar P0     = stg_ctx->P0;
31143bff553SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
31243bff553SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
31343bff553SJames Wright   const CeedScalar Rd     = cp - cv;
31443bff553SJames Wright   const CeedScalar rho    = P0 / (Rd * theta0);
31543bff553SJames Wright 
31643bff553SJames Wright   CeedPragmaSIMD
31743bff553SJames Wright   for(CeedInt i=0; i<Q; i++) {
31843bff553SJames Wright     InterpolateProfile(X[1][i], u, cij, &eps, &lt, stg_ctx);
31943bff553SJames Wright 
32043bff553SJames Wright     q0[0][i] = rho;
32143bff553SJames Wright     q0[1][i] = u[0] * rho;
32243bff553SJames Wright     q0[2][i] = u[1] * rho;
32343bff553SJames Wright     q0[3][i] = u[2] * rho;
32443bff553SJames Wright     q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0);
32543bff553SJames Wright   } // End of Quadrature Point Loop
32643bff553SJames Wright   return 0;
32743bff553SJames Wright }
32843bff553SJames Wright 
329493642f1SJames Wright /********************************************************************
330493642f1SJames Wright  * @brief QFunction to calculate the inflow boundary condition
331493642f1SJames Wright  *
332493642f1SJames Wright  * This will loop through quadrature points, calculate the wavemode amplitudes
333493642f1SJames Wright  * at each location, then calculate the actual velocity.
334493642f1SJames Wright  */
335493642f1SJames Wright CEED_QFUNCTION(STGShur14_Inflow)(void *ctx, CeedInt Q,
336493642f1SJames Wright                                  const CeedScalar *const *in,
337493642f1SJames Wright                                  CeedScalar *const *out) {
338493642f1SJames Wright 
339493642f1SJames Wright   //*INDENT-OFF*
340493642f1SJames Wright   const CeedScalar (*q)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
341dd64951cSJames Wright                    (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[2],
342dd64951cSJames Wright                    (*X)[CEED_Q_VLA]          = (const CeedScalar(*)[CEED_Q_VLA]) in[3];
343493642f1SJames Wright 
344a6e8f989SJames Wright   CeedScalar(*v)[CEED_Q_VLA]            = (CeedScalar(*)[CEED_Q_VLA]) out[0],
345a6e8f989SJames Wright             (*jac_data_sur)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[1];
346493642f1SJames Wright 
347493642f1SJames Wright   //*INDENT-ON*
348493642f1SJames Wright 
349493642f1SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
350493642f1SJames Wright   CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
351493642f1SJames Wright   const bool is_implicit  = stg_ctx->is_implicit;
352493642f1SJames Wright   const bool mean_only    = stg_ctx->mean_only;
353493642f1SJames Wright   const bool prescribe_T  = stg_ctx->prescribe_T;
354493642f1SJames Wright   const CeedScalar dx     = stg_ctx->dx;
355493642f1SJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
356493642f1SJames Wright   const CeedScalar time   = stg_ctx->time;
357493642f1SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
358493642f1SJames Wright   const CeedScalar P0     = stg_ctx->P0;
359493642f1SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
360493642f1SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
361493642f1SJames Wright   const CeedScalar Rd     = cp - cv;
362493642f1SJames Wright   const CeedScalar gamma  = cp/cv;
363493642f1SJames Wright 
364493642f1SJames Wright   CeedPragmaSIMD
365493642f1SJames Wright   for(CeedInt i=0; i<Q; i++) {
366493642f1SJames Wright     const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0);
367493642f1SJames Wright     const CeedScalar x[] = { X[0][i], X[1][i], X[2][i] };
368493642f1SJames Wright     const CeedScalar dXdx[2][3] = {
369493642f1SJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
370493642f1SJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
371493642f1SJames Wright     };
372493642f1SJames Wright 
373493642f1SJames Wright     CeedScalar h[3];
374493642f1SJames Wright     h[0] = dx;
375f6438f20SJames Wright     for (CeedInt j=1; j<3; j++)
376f6438f20SJames Wright       h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]));
377493642f1SJames Wright 
378493642f1SJames Wright     InterpolateProfile(X[1][i], ubar, cij, &eps, &lt, stg_ctx);
379493642f1SJames Wright     if (!mean_only) {
380493642f1SJames Wright       CalcSpectrum(X[1][i], eps, lt, h, mu/rho, qn, stg_ctx);
381493642f1SJames Wright       STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
382493642f1SJames Wright     } else {
383493642f1SJames Wright       for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
384493642f1SJames Wright     }
385493642f1SJames Wright 
386a6e8f989SJames Wright     const CeedScalar E_kinetic = .5 * rho * Dot3(u, u);
387493642f1SJames Wright     CeedScalar E_internal, P;
388493642f1SJames Wright     if (prescribe_T) {
389493642f1SJames Wright       // Temperature is being set weakly (theta0) and for constant cv this sets E_internal
390493642f1SJames Wright       E_internal = rho * cv * theta0;
391493642f1SJames Wright       // Find pressure using
392493642f1SJames Wright       P = rho * Rd * theta0; // interior rho with exterior T
393493642f1SJames Wright     } else {
394493642f1SJames Wright       E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution
395493642f1SJames Wright       P = E_internal * (gamma - 1.);
396493642f1SJames Wright     }
397493642f1SJames Wright 
398493642f1SJames Wright     const CeedScalar wdetJb  = (is_implicit ? -1. : 1.) * q_data_sur[0][i];
399493642f1SJames Wright     // ---- Normal vect
400493642f1SJames Wright     const CeedScalar norm[3] = {q_data_sur[1][i],
401493642f1SJames Wright                                 q_data_sur[2][i],
402493642f1SJames Wright                                 q_data_sur[3][i]
403493642f1SJames Wright                                };
404493642f1SJames Wright 
405493642f1SJames Wright     const CeedScalar E = E_internal + E_kinetic;
406493642f1SJames Wright 
407493642f1SJames Wright     // Velocity normal to the boundary
408a6e8f989SJames Wright     const CeedScalar u_normal = Dot3(norm, u);
409a6e8f989SJames Wright 
410493642f1SJames Wright     // The Physics
411493642f1SJames Wright     // Zero v so all future terms can safely sum into it
412493642f1SJames Wright     for (CeedInt j=0; j<5; j++) v[j][i] = 0.;
413493642f1SJames Wright 
414493642f1SJames Wright     // The Physics
415493642f1SJames Wright     // -- Density
416493642f1SJames Wright     v[0][i] -= wdetJb * rho * u_normal;
417493642f1SJames Wright 
418493642f1SJames Wright     // -- Momentum
419493642f1SJames Wright     for (CeedInt j=0; j<3; j++)
420493642f1SJames Wright       v[j+1][i] -= wdetJb *(rho * u_normal * u[j] +
421493642f1SJames Wright                             norm[j] * P);
422493642f1SJames Wright 
423493642f1SJames Wright     // -- Total Energy Density
424493642f1SJames Wright     v[4][i] -= wdetJb * u_normal * (E + P);
425a6e8f989SJames Wright 
426a6e8f989SJames Wright     jac_data_sur[0][i] = rho;
427a6e8f989SJames Wright     jac_data_sur[1][i] = u[0];
428a6e8f989SJames Wright     jac_data_sur[2][i] = u[1];
429a6e8f989SJames Wright     jac_data_sur[3][i] = u[2];
430a6e8f989SJames Wright     jac_data_sur[4][i] = E;
431a6e8f989SJames Wright     for (int j=0; j<6; j++) jac_data_sur[5+j][i] = 0.;
432493642f1SJames Wright   }
433493642f1SJames Wright   return 0;
434493642f1SJames Wright }
435493642f1SJames Wright 
436a6e8f989SJames Wright CEED_QFUNCTION(STGShur14_Inflow_Jacobian)(void *ctx, CeedInt Q,
437a6e8f989SJames Wright     const CeedScalar *const *in,
438a6e8f989SJames Wright     CeedScalar *const *out) {
439a6e8f989SJames Wright   // *INDENT-OFF*
440a6e8f989SJames Wright   // Inputs
441a6e8f989SJames Wright   const CeedScalar (*dq)[CEED_Q_VLA]           = (const CeedScalar(*)[CEED_Q_VLA])in[0],
442a6e8f989SJames Wright                    (*q_data_sur)[CEED_Q_VLA]   = (const CeedScalar(*)[CEED_Q_VLA])in[2],
443a6e8f989SJames Wright                    (*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4];
444a6e8f989SJames Wright   // Outputs
445a6e8f989SJames Wright   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
446a6e8f989SJames Wright   // *INDENT-ON*
447a6e8f989SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
448a6e8f989SJames Wright   const bool implicit     = stg_ctx->is_implicit;
449a6e8f989SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
450a6e8f989SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
451a6e8f989SJames Wright   const CeedScalar Rd     = cp - cv;
452a6e8f989SJames Wright   const CeedScalar gamma  = cp/cv;
453a6e8f989SJames Wright 
454a6e8f989SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
455a6e8f989SJames Wright   const bool prescribe_T  = stg_ctx->prescribe_T;
456a6e8f989SJames Wright 
457a6e8f989SJames Wright   CeedPragmaSIMD
458a6e8f989SJames Wright   // Quadrature Point Loop
459a6e8f989SJames Wright   for (CeedInt i=0; i<Q; i++) {
460a6e8f989SJames Wright     // Setup
461a6e8f989SJames Wright     // Setup
462a6e8f989SJames Wright     // -- Interp-to-Interp q_data
463a6e8f989SJames Wright     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
464a6e8f989SJames Wright     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
465a6e8f989SJames Wright     // We can effect this by swapping the sign on this weight
466a6e8f989SJames Wright     const CeedScalar wdetJb  = (implicit ? -1. : 1.) * q_data_sur[0][i];
467a6e8f989SJames Wright 
468a6e8f989SJames Wright     // Calculate inflow values
469a6e8f989SJames Wright     CeedScalar velocity[3];
470a6e8f989SJames Wright     for (CeedInt j=0; j<3; j++) velocity[j] = jac_data_sur[5+j][i];
471a6e8f989SJames Wright 
472a6e8f989SJames Wright     // enabling user to choose between weak T and weak rho inflow
473a6e8f989SJames Wright     CeedScalar drho, dE, dP;
474a6e8f989SJames Wright     if (prescribe_T) {
475a6e8f989SJames Wright       // rho should be from the current solution
476a6e8f989SJames Wright       drho = dq[0][i];
477a6e8f989SJames Wright       CeedScalar dE_internal = drho * cv * theta0;
478a6e8f989SJames Wright       CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity);
479a6e8f989SJames Wright       dE = dE_internal + dE_kinetic;
480a6e8f989SJames Wright       dP = drho * Rd * theta0; // interior rho with exterior T
481a6e8f989SJames Wright     } else { // rho specified, E_internal from solution
482a6e8f989SJames Wright       drho = 0;
483a6e8f989SJames Wright       dE = dq[4][i];
484a6e8f989SJames Wright       dP = dE * (gamma - 1.);
485a6e8f989SJames Wright     }
486a6e8f989SJames Wright     const CeedScalar norm[3] = {q_data_sur[1][i],
487a6e8f989SJames Wright                                 q_data_sur[2][i],
488a6e8f989SJames Wright                                 q_data_sur[3][i]
489a6e8f989SJames Wright                                };
490a6e8f989SJames Wright 
491a6e8f989SJames Wright     const CeedScalar u_normal = Dot3(norm, velocity);
492a6e8f989SJames Wright 
493a6e8f989SJames Wright     v[0][i] = - wdetJb * drho * u_normal;
494a6e8f989SJames Wright     for (int j=0; j<3; j++)
495a6e8f989SJames Wright       v[j+1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP);
496a6e8f989SJames Wright     v[4][i] = - wdetJb * u_normal * (dE + dP);
497a6e8f989SJames Wright   } // End Quadrature Point Loop
498a6e8f989SJames Wright   return 0;
499a6e8f989SJames Wright }
500a6e8f989SJames Wright 
501b7190ff7SJames Wright /********************************************************************
502b7190ff7SJames Wright  * @brief QFunction to calculate the strongly enforce inflow BC
503b7190ff7SJames Wright  *
504b7190ff7SJames Wright  * This QF is for the strong application of STG via libCEED (rather than
505b7190ff7SJames Wright  * through the native PETSc `DMAddBoundary` -> `bcFunc` method.
506b7190ff7SJames Wright  */
507b7190ff7SJames Wright CEED_QFUNCTION(STGShur14_Inflow_StrongQF)(void *ctx, CeedInt Q,
508b7190ff7SJames Wright     const CeedScalar *const *in, CeedScalar *const *out) {
509b7190ff7SJames Wright 
510b7190ff7SJames Wright   //*INDENT-OFF*
511b7190ff7SJames Wright   const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0],
512b7190ff7SJames Wright                    (*coords)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA]) in[1],
513*9eeef72bSJames Wright                    (*scale)                  = (const CeedScalar(*)) in[2],
514*9eeef72bSJames Wright                    (*stg_data)[CEED_Q_VLA]   = (const CeedScalar(*)[CEED_Q_VLA]) in[3];
515b7190ff7SJames Wright 
516b7190ff7SJames Wright   CeedScalar(*bcval)[CEED_Q_VLA]            = (CeedScalar(*)[CEED_Q_VLA]) out[0];
517b7190ff7SJames Wright   //*INDENT-ON*
518b7190ff7SJames Wright 
519b7190ff7SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
520b7190ff7SJames Wright   CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt;
521b7190ff7SJames Wright   const bool mean_only    = stg_ctx->mean_only;
522b7190ff7SJames Wright   const CeedScalar dx     = stg_ctx->dx;
523b7190ff7SJames Wright   const CeedScalar mu     = stg_ctx->newtonian_ctx.mu;
524b7190ff7SJames Wright   const CeedScalar time   = stg_ctx->time;
525b7190ff7SJames Wright   const CeedScalar theta0 = stg_ctx->theta0;
526b7190ff7SJames Wright   const CeedScalar P0     = stg_ctx->P0;
527b7190ff7SJames Wright   const CeedScalar cv     = stg_ctx->newtonian_ctx.cv;
528b7190ff7SJames Wright   const CeedScalar cp     = stg_ctx->newtonian_ctx.cp;
529b7190ff7SJames Wright   const CeedScalar Rd     = cp - cv;
530b7190ff7SJames Wright   const CeedScalar rho    = P0 / (Rd * theta0);
531b7190ff7SJames Wright 
532b7190ff7SJames Wright   CeedPragmaSIMD
533b7190ff7SJames Wright   for(CeedInt i=0; i<Q; i++) {
534b7190ff7SJames Wright     const CeedScalar x[] = { coords[0][i], coords[1][i], coords[2][i] };
535b7190ff7SJames Wright     const CeedScalar dXdx[2][3] = {
536b7190ff7SJames Wright       {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]},
537b7190ff7SJames Wright       {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]}
538b7190ff7SJames Wright     };
539b7190ff7SJames Wright 
540b7190ff7SJames Wright     CeedScalar h[3];
541b7190ff7SJames Wright     h[0] = dx;
542f6438f20SJames Wright     for (CeedInt j=1; j<3; j++)
543f6438f20SJames Wright       h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j]));
544b7190ff7SJames Wright 
545b7190ff7SJames Wright     InterpolateProfile(coords[1][i], ubar, cij, &eps, &lt, stg_ctx);
546b7190ff7SJames Wright     if (!mean_only) {
547*9eeef72bSJames Wright       STGShur14_Calc_PrecompEktot(x, time, ubar, cij, stg_data[0][i],
548*9eeef72bSJames Wright                                   h, x[1], eps, lt, mu/rho, u, stg_ctx);
549*9eeef72bSJames Wright       // CalcSpectrum(coords[1][i], eps, lt, h, mu/rho, qn, stg_ctx);
550*9eeef72bSJames Wright       // STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
551b7190ff7SJames Wright     } else {
552b7190ff7SJames Wright       for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
553b7190ff7SJames Wright     }
554b7190ff7SJames Wright 
555b7190ff7SJames Wright     bcval[0][i] = scale[i] * rho;
556b7190ff7SJames Wright     bcval[1][i] = scale[i] * rho * u[0];
557b7190ff7SJames Wright     bcval[2][i] = scale[i] * rho * u[1];
558b7190ff7SJames Wright     bcval[3][i] = scale[i] * rho * u[2];
55966531c8bSJames Wright     bcval[4][i] = 0.;
560b7190ff7SJames Wright   }
561b7190ff7SJames Wright   return 0;
562b7190ff7SJames Wright }
563b7190ff7SJames Wright 
564493642f1SJames Wright #endif // stg_shur14_h
565