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 /* 91ba6664aeSJames Wright * @brief Calculate spectrum coefficients for STG 92ba6664aeSJames Wright * 93ba6664aeSJames Wright * Calculates q_n at a given distance to the wall 94ba6664aeSJames Wright * 95ba6664aeSJames Wright * @param[in] dw Distance to the nearest wall 96ba6664aeSJames Wright * @param[in] eps Turbulent dissipation w/rt dw 97ba6664aeSJames Wright * @param[in] lt Turbulent length scale w/rt dw 98ba6664aeSJames Wright * @param[in] h Element lengths in coordinate directions 99ba6664aeSJames Wright * @param[in] nu Dynamic Viscosity; 100ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 101ba6664aeSJames Wright * @param[out] qn Spectrum coefficients, [nmodes] 102ba6664aeSJames Wright */ 103ba6664aeSJames Wright void CEED_QFUNCTION_HELPER(CalcSpectrum)(const CeedScalar dw, 104ba6664aeSJames Wright const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 105ba6664aeSJames Wright const CeedScalar nu, CeedScalar qn[], const STGShur14Context stg_ctx) { 106ba6664aeSJames Wright 107ba6664aeSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 108ba6664aeSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 109ba6664aeSJames Wright 110ba6664aeSJames Wright const CeedScalar hmax = Max( Max(h[0], h[1]), h[2]); 111cfcf1481SJames Wright const CeedScalar ke = dw==0 ? 1e16 : 2*M_PI/Min(2*dw, 3*lt); 11213fa47b2SJames Wright const CeedScalar keta = 2*M_PI*pow(Cube(nu)/eps, -0.25); 113ba6664aeSJames Wright const CeedScalar kcut = 114ba6664aeSJames Wright M_PI/ Min( Max(Max(h[1], h[2]), 0.3*hmax) + 0.1*dw, hmax ); 115ba6664aeSJames Wright CeedScalar fcut, feta, Ektot=0.0; 116ba6664aeSJames Wright 117ba6664aeSJames Wright for(CeedInt n=0; n<nmodes; n++) { 118ba6664aeSJames Wright feta = exp(-Square(12*kappa[n]/keta)); 119ba6664aeSJames Wright fcut = exp( -pow(4*Max(kappa[n] - 0.9*kcut, 0)/kcut, 3.) ); 120ba6664aeSJames Wright qn[n] = pow(kappa[n]/ke, 4.) 121ba6664aeSJames Wright * pow(1 + 2.4*Square(kappa[n]/ke),-17./6)*feta*fcut; 122ba6664aeSJames Wright qn[n] *= n==0 ? kappa[0] : kappa[n] - kappa[n-1]; 123ba6664aeSJames Wright Ektot += qn[n]; 124ba6664aeSJames Wright } 125ba6664aeSJames Wright 126961c9c98SJames Wright if (Ektot == 0) return; 127ba6664aeSJames Wright for(CeedInt n=0; n<nmodes; n++) qn[n] /= Ektot; 128ba6664aeSJames Wright } 129ba6664aeSJames Wright 130ba6664aeSJames Wright /****************************************************** 131ba6664aeSJames Wright * @brief Calculate u(x,t) for STG inflow condition 132ba6664aeSJames Wright * 133ba6664aeSJames Wright * @param[in] X Location to evaluate u(X,t) 134ba6664aeSJames Wright * @param[in] t Time to evaluate u(X,t) 135ba6664aeSJames Wright * @param[in] ubar Mean velocity at X 136ba6664aeSJames Wright * @param[in] cij Cholesky decomposition at X 137ba6664aeSJames Wright * @param[in] qn Wavemode amplitudes at X, [nmodes] 138ba6664aeSJames Wright * @param[out] u Velocity at X and t 139ba6664aeSJames Wright * @param[in] stg_ctx STGShur14Context for the problem 140ba6664aeSJames Wright */ 141ba6664aeSJames Wright void CEED_QFUNCTION_HELPER(STGShur14_Calc)(const CeedScalar X[3], 142ba6664aeSJames Wright const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 143ba6664aeSJames Wright const CeedScalar qn[], CeedScalar u[3], 144ba6664aeSJames Wright const STGShur14Context stg_ctx) { 145ba6664aeSJames Wright 146ba6664aeSJames Wright //*INDENT-OFF* 147ba6664aeSJames Wright const CeedInt nmodes = stg_ctx->nmodes; 148ba6664aeSJames Wright const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 149ba6664aeSJames Wright const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 150ba6664aeSJames Wright const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 151ba6664aeSJames Wright const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 152ba6664aeSJames Wright //*INDENT-ON* 153ba6664aeSJames Wright CeedScalar xdotd, vp[3] = {0.}; 154ba6664aeSJames Wright CeedScalar xhat[] = {0., X[1], X[2]}; 155ba6664aeSJames Wright 156ba6664aeSJames Wright CeedPragmaSIMD 157ba6664aeSJames Wright for(CeedInt n=0; n<nmodes; n++) { 158ba6664aeSJames Wright xhat[0] = (X[0] - stg_ctx->u0*t)*Max(2*kappa[0]/kappa[n], 0.1); 159ba6664aeSJames Wright xdotd = 0.; 160ba6664aeSJames Wright for(CeedInt i=0; i<3; i++) xdotd += d[i*nmodes+n]*xhat[i]; 161ba6664aeSJames Wright const CeedScalar cos_kxdp = cos(kappa[n]*xdotd + phi[n]); 162961c9c98SJames Wright vp[0] += sqrt(qn[n])*sigma[0*nmodes+n] * cos_kxdp; 163961c9c98SJames Wright vp[1] += sqrt(qn[n])*sigma[1*nmodes+n] * cos_kxdp; 164961c9c98SJames Wright vp[2] += sqrt(qn[n])*sigma[2*nmodes+n] * cos_kxdp; 165ba6664aeSJames Wright } 166961c9c98SJames Wright for(CeedInt i=0; i<3; i++) vp[i] *= 2*sqrt(1.5); 167ba6664aeSJames Wright 168ba6664aeSJames Wright u[0] = ubar[0] + cij[0]*vp[0]; 169ba6664aeSJames Wright u[1] = ubar[1] + cij[3]*vp[0] + cij[1]*vp[1]; 170ba6664aeSJames Wright u[2] = ubar[2] + cij[4]*vp[0] + cij[5]*vp[1] + cij[2]*vp[2]; 171ba6664aeSJames Wright } 172ba6664aeSJames Wright 173b77c53c9SJames Wright // Extrude the STGInflow profile through out the domain for an initial condition 174b77c53c9SJames Wright CEED_QFUNCTION(ICsSTG)(void *ctx, CeedInt Q, 175b77c53c9SJames Wright const CeedScalar *const *in, CeedScalar *const *out) { 176b77c53c9SJames Wright // Inputs 177b77c53c9SJames Wright const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 178b77c53c9SJames Wright 179b77c53c9SJames Wright // Outputs 180b77c53c9SJames Wright CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 181b77c53c9SJames Wright 182b77c53c9SJames Wright const STGShur14Context stg_ctx = (STGShur14Context) ctx; 183b77c53c9SJames Wright CeedScalar u[3], cij[6], eps, lt; 184b77c53c9SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 185b77c53c9SJames Wright const CeedScalar P0 = stg_ctx->P0; 186b77c53c9SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 187b77c53c9SJames Wright const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 188b77c53c9SJames Wright const CeedScalar Rd = cp - cv; 189b77c53c9SJames Wright const CeedScalar rho = P0 / (Rd * theta0); 190b77c53c9SJames Wright 191b77c53c9SJames Wright CeedPragmaSIMD 192b77c53c9SJames Wright for(CeedInt i=0; i<Q; i++) { 193b77c53c9SJames Wright InterpolateProfile(X[1][i], u, cij, &eps, <, stg_ctx); 194b77c53c9SJames Wright 195b77c53c9SJames Wright q0[0][i] = rho; 196b77c53c9SJames Wright q0[1][i] = u[0] * rho; 197b77c53c9SJames Wright q0[2][i] = u[1] * rho; 198b77c53c9SJames Wright q0[3][i] = u[2] * rho; 199b77c53c9SJames Wright q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0); 200b77c53c9SJames Wright } // End of Quadrature Point Loop 201b77c53c9SJames Wright return 0; 202b77c53c9SJames Wright } 203b77c53c9SJames Wright 204ba6664aeSJames Wright /******************************************************************** 205ba6664aeSJames Wright * @brief QFunction to calculate the inflow boundary condition 206ba6664aeSJames Wright * 207ba6664aeSJames Wright * This will loop through quadrature points, calculate the wavemode amplitudes 208ba6664aeSJames Wright * at each location, then calculate the actual velocity. 209ba6664aeSJames Wright */ 210ba6664aeSJames Wright CEED_QFUNCTION(STGShur14_Inflow)(void *ctx, CeedInt Q, 211ba6664aeSJames Wright const CeedScalar *const *in, 212ba6664aeSJames Wright CeedScalar *const *out) { 213ba6664aeSJames Wright 214ba6664aeSJames Wright //*INDENT-OFF* 215ba6664aeSJames Wright const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0], 216e8b03feeSJames Wright (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[2], 217e8b03feeSJames Wright (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[3]; 218ba6664aeSJames Wright 2194dbab5e5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[0], 2204dbab5e5SJames Wright (*jac_data_sur)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[1]; 221ba6664aeSJames Wright 222ba6664aeSJames Wright //*INDENT-ON* 223ba6664aeSJames Wright 224ba6664aeSJames Wright const STGShur14Context stg_ctx = (STGShur14Context) ctx; 225ba6664aeSJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 226ba6664aeSJames Wright const bool is_implicit = stg_ctx->is_implicit; 227ba6664aeSJames Wright const bool mean_only = stg_ctx->mean_only; 228ba6664aeSJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 229ba6664aeSJames Wright const CeedScalar dx = stg_ctx->dx; 230ba6664aeSJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 231ba6664aeSJames Wright const CeedScalar time = stg_ctx->time; 232ba6664aeSJames Wright const CeedScalar theta0 = stg_ctx->theta0; 233ba6664aeSJames Wright const CeedScalar P0 = stg_ctx->P0; 234ba6664aeSJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 235ba6664aeSJames Wright const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 236ba6664aeSJames Wright const CeedScalar Rd = cp - cv; 237ba6664aeSJames Wright const CeedScalar gamma = cp/cv; 238ba6664aeSJames Wright 239ba6664aeSJames Wright CeedPragmaSIMD 240ba6664aeSJames Wright for(CeedInt i=0; i<Q; i++) { 241ba6664aeSJames Wright const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0); 242ba6664aeSJames Wright const CeedScalar x[] = { X[0][i], X[1][i], X[2][i] }; 243ba6664aeSJames Wright const CeedScalar dXdx[2][3] = { 244ba6664aeSJames Wright {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 245ba6664aeSJames Wright {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 246ba6664aeSJames Wright }; 247ba6664aeSJames Wright 248ba6664aeSJames Wright CeedScalar h[3]; 249ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 25013fa47b2SJames Wright h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 251ba6664aeSJames Wright h[0] = dx; 252ba6664aeSJames Wright 253ba6664aeSJames Wright InterpolateProfile(X[1][i], ubar, cij, &eps, <, stg_ctx); 254ba6664aeSJames Wright if (!mean_only) { 255ba6664aeSJames Wright CalcSpectrum(X[1][i], eps, lt, h, mu/rho, qn, stg_ctx); 256ba6664aeSJames Wright STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 257ba6664aeSJames Wright } else { 258ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) u[j] = ubar[j]; 259ba6664aeSJames Wright } 260ba6664aeSJames Wright 2614dbab5e5SJames Wright const CeedScalar E_kinetic = .5 * rho * Dot3(u, u); 262ba6664aeSJames Wright CeedScalar E_internal, P; 263ba6664aeSJames Wright if (prescribe_T) { 264ba6664aeSJames Wright // Temperature is being set weakly (theta0) and for constant cv this sets E_internal 265ba6664aeSJames Wright E_internal = rho * cv * theta0; 266ba6664aeSJames Wright // Find pressure using 267ba6664aeSJames Wright P = rho * Rd * theta0; // interior rho with exterior T 268ba6664aeSJames Wright } else { 269ba6664aeSJames Wright E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution 270ba6664aeSJames Wright P = E_internal * (gamma - 1.); 271ba6664aeSJames Wright } 272ba6664aeSJames Wright 273ba6664aeSJames Wright const CeedScalar wdetJb = (is_implicit ? -1. : 1.) * q_data_sur[0][i]; 274ba6664aeSJames Wright // ---- Normal vect 275ba6664aeSJames Wright const CeedScalar norm[3] = {q_data_sur[1][i], 276ba6664aeSJames Wright q_data_sur[2][i], 277ba6664aeSJames Wright q_data_sur[3][i] 278ba6664aeSJames Wright }; 279ba6664aeSJames Wright 280ba6664aeSJames Wright const CeedScalar E = E_internal + E_kinetic; 281ba6664aeSJames Wright 282ba6664aeSJames Wright // Velocity normal to the boundary 2834dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, u); 2844dbab5e5SJames Wright 285ba6664aeSJames Wright // The Physics 286ba6664aeSJames Wright // Zero v so all future terms can safely sum into it 287ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) v[j][i] = 0.; 288ba6664aeSJames Wright 289ba6664aeSJames Wright // The Physics 290ba6664aeSJames Wright // -- Density 291ba6664aeSJames Wright v[0][i] -= wdetJb * rho * u_normal; 292ba6664aeSJames Wright 293ba6664aeSJames Wright // -- Momentum 294ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 295ba6664aeSJames Wright v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + 296ba6664aeSJames Wright norm[j] * P); 297ba6664aeSJames Wright 298ba6664aeSJames Wright // -- Total Energy Density 299ba6664aeSJames Wright v[4][i] -= wdetJb * u_normal * (E + P); 3004dbab5e5SJames Wright 3014dbab5e5SJames Wright jac_data_sur[0][i] = rho; 3024dbab5e5SJames Wright jac_data_sur[1][i] = u[0]; 3034dbab5e5SJames Wright jac_data_sur[2][i] = u[1]; 3044dbab5e5SJames Wright jac_data_sur[3][i] = u[2]; 3054dbab5e5SJames Wright jac_data_sur[4][i] = E; 3064dbab5e5SJames Wright for (int j=0; j<6; j++) jac_data_sur[5+j][i] = 0.; 307ba6664aeSJames Wright } 308ba6664aeSJames Wright return 0; 309ba6664aeSJames Wright } 310ba6664aeSJames Wright 3114dbab5e5SJames Wright CEED_QFUNCTION(STGShur14_Inflow_Jacobian)(void *ctx, CeedInt Q, 3124dbab5e5SJames Wright const CeedScalar *const *in, 3134dbab5e5SJames Wright CeedScalar *const *out) { 3144dbab5e5SJames Wright // *INDENT-OFF* 3154dbab5e5SJames Wright // Inputs 3164dbab5e5SJames Wright const CeedScalar (*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 3174dbab5e5SJames Wright (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 3184dbab5e5SJames Wright (*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 3194dbab5e5SJames Wright // Outputs 3204dbab5e5SJames Wright CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 3214dbab5e5SJames Wright // *INDENT-ON* 3224dbab5e5SJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 3234dbab5e5SJames Wright const bool implicit = stg_ctx->is_implicit; 3244dbab5e5SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 3254dbab5e5SJames Wright const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 3264dbab5e5SJames Wright const CeedScalar Rd = cp - cv; 3274dbab5e5SJames Wright const CeedScalar gamma = cp/cv; 3284dbab5e5SJames Wright 3294dbab5e5SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 3304dbab5e5SJames Wright const bool prescribe_T = stg_ctx->prescribe_T; 3314dbab5e5SJames Wright 3324dbab5e5SJames Wright CeedPragmaSIMD 3334dbab5e5SJames Wright // Quadrature Point Loop 3344dbab5e5SJames Wright for (CeedInt i=0; i<Q; i++) { 3354dbab5e5SJames Wright // Setup 3364dbab5e5SJames Wright // Setup 3374dbab5e5SJames Wright // -- Interp-to-Interp q_data 3384dbab5e5SJames Wright // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 3394dbab5e5SJames Wright // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 3404dbab5e5SJames Wright // We can effect this by swapping the sign on this weight 3414dbab5e5SJames Wright const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 3424dbab5e5SJames Wright 3434dbab5e5SJames Wright // Calculate inflow values 3444dbab5e5SJames Wright CeedScalar velocity[3]; 3454dbab5e5SJames Wright for (CeedInt j=0; j<3; j++) velocity[j] = jac_data_sur[5+j][i]; 3464dbab5e5SJames Wright 3474dbab5e5SJames Wright // enabling user to choose between weak T and weak rho inflow 3484dbab5e5SJames Wright CeedScalar drho, dE, dP; 3494dbab5e5SJames Wright if (prescribe_T) { 3504dbab5e5SJames Wright // rho should be from the current solution 3514dbab5e5SJames Wright drho = dq[0][i]; 3524dbab5e5SJames Wright CeedScalar dE_internal = drho * cv * theta0; 3534dbab5e5SJames Wright CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity); 3544dbab5e5SJames Wright dE = dE_internal + dE_kinetic; 3554dbab5e5SJames Wright dP = drho * Rd * theta0; // interior rho with exterior T 3564dbab5e5SJames Wright } else { // rho specified, E_internal from solution 3574dbab5e5SJames Wright drho = 0; 3584dbab5e5SJames Wright dE = dq[4][i]; 3594dbab5e5SJames Wright dP = dE * (gamma - 1.); 3604dbab5e5SJames Wright } 3614dbab5e5SJames Wright const CeedScalar norm[3] = {q_data_sur[1][i], 3624dbab5e5SJames Wright q_data_sur[2][i], 3634dbab5e5SJames Wright q_data_sur[3][i] 3644dbab5e5SJames Wright }; 3654dbab5e5SJames Wright 3664dbab5e5SJames Wright const CeedScalar u_normal = Dot3(norm, velocity); 3674dbab5e5SJames Wright 3684dbab5e5SJames Wright v[0][i] = - wdetJb * drho * u_normal; 3694dbab5e5SJames Wright for (int j=0; j<3; j++) 3704dbab5e5SJames Wright v[j+1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP); 3714dbab5e5SJames Wright v[4][i] = - wdetJb * u_normal * (dE + dP); 3724dbab5e5SJames Wright } // End Quadrature Point Loop 3734dbab5e5SJames Wright return 0; 3744dbab5e5SJames Wright } 3754dbab5e5SJames Wright 376*0a6353c2SJames Wright /******************************************************************** 377*0a6353c2SJames Wright * @brief QFunction to calculate the strongly enforce inflow BC 378*0a6353c2SJames Wright * 379*0a6353c2SJames Wright * This QF is for the strong application of STG via libCEED (rather than 380*0a6353c2SJames Wright * through the native PETSc `DMAddBoundary` -> `bcFunc` method. 381*0a6353c2SJames Wright */ 382*0a6353c2SJames Wright CEED_QFUNCTION(STGShur14_Inflow_StrongQF)(void *ctx, CeedInt Q, 383*0a6353c2SJames Wright const CeedScalar *const *in, CeedScalar *const *out) { 384*0a6353c2SJames Wright 385*0a6353c2SJames Wright //*INDENT-OFF* 386*0a6353c2SJames Wright const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0], 387*0a6353c2SJames Wright (*coords)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[1], 388*0a6353c2SJames Wright (*scale) = (const CeedScalar(*)) in[2]; 389*0a6353c2SJames Wright 390*0a6353c2SJames Wright CeedScalar(*bcval)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[0]; 391*0a6353c2SJames Wright //*INDENT-ON* 392*0a6353c2SJames Wright 393*0a6353c2SJames Wright const STGShur14Context stg_ctx = (STGShur14Context) ctx; 394*0a6353c2SJames Wright CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 395*0a6353c2SJames Wright const bool mean_only = stg_ctx->mean_only; 396*0a6353c2SJames Wright const CeedScalar dx = stg_ctx->dx; 397*0a6353c2SJames Wright const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 398*0a6353c2SJames Wright const CeedScalar time = stg_ctx->time; 399*0a6353c2SJames Wright const CeedScalar theta0 = stg_ctx->theta0; 400*0a6353c2SJames Wright const CeedScalar P0 = stg_ctx->P0; 401*0a6353c2SJames Wright const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 402*0a6353c2SJames Wright const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 403*0a6353c2SJames Wright const CeedScalar Rd = cp - cv; 404*0a6353c2SJames Wright const CeedScalar rho = P0 / (Rd * theta0); 405*0a6353c2SJames Wright 406*0a6353c2SJames Wright CeedPragmaSIMD 407*0a6353c2SJames Wright for(CeedInt i=0; i<Q; i++) { 408*0a6353c2SJames Wright const CeedScalar x[] = { coords[0][i], coords[1][i], coords[2][i] }; 409*0a6353c2SJames Wright const CeedScalar dXdx[2][3] = { 410*0a6353c2SJames Wright {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 411*0a6353c2SJames Wright {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 412*0a6353c2SJames Wright }; 413*0a6353c2SJames Wright 414*0a6353c2SJames Wright CeedScalar h[3]; 415*0a6353c2SJames Wright for (CeedInt j=0; j<3; j++) 416*0a6353c2SJames Wright h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 417*0a6353c2SJames Wright h[0] = dx; 418*0a6353c2SJames Wright 419*0a6353c2SJames Wright InterpolateProfile(coords[1][i], ubar, cij, &eps, <, stg_ctx); 420*0a6353c2SJames Wright if (!mean_only) { 421*0a6353c2SJames Wright CalcSpectrum(coords[1][i], eps, lt, h, mu/rho, qn, stg_ctx); 422*0a6353c2SJames Wright STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 423*0a6353c2SJames Wright } else { 424*0a6353c2SJames Wright for (CeedInt j=0; j<3; j++) u[j] = ubar[j]; 425*0a6353c2SJames Wright } 426*0a6353c2SJames Wright 427*0a6353c2SJames Wright bcval[0][i] = scale[i] * rho; 428*0a6353c2SJames Wright bcval[1][i] = scale[i] * rho * u[0]; 429*0a6353c2SJames Wright bcval[2][i] = scale[i] * rho * u[1]; 430*0a6353c2SJames Wright bcval[3][i] = scale[i] * rho * u[2]; 431*0a6353c2SJames Wright } 432*0a6353c2SJames Wright return 0; 433*0a6353c2SJames Wright } 434*0a6353c2SJames Wright 435ba6664aeSJames Wright #endif // stg_shur14_h 436