1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Implementation of the Synthetic Turbulence Generation (STG) algorithm 10 /// presented in Shur et al. 2014 11 // 12 /// SetupSTG_Rand reads in the input files and fills in STGShur14Context. Then 13 /// STGShur14_CalcQF is run over quadrature points. Before the program exits, 14 /// TearDownSTG is run to free the memory of the allocated arrays. 15 16 #ifndef stg_shur14_h 17 #define stg_shur14_h 18 19 #include <math.h> 20 #include <ceed.h> 21 #include <stdlib.h> 22 #include "stg_shur14_type.h" 23 #include "utils.h" 24 25 #define STG_NMODES_MAX 1024 26 27 /* 28 * @brief Interpolate quantities from input profile to given location 29 * 30 * Assumed that prof_dw[i+1] > prof_dw[i] and prof_dw[0] = 0 31 * If dw > prof_dw[-1], then the interpolation takes the values at prof_dw[-1] 32 * 33 * @param[in] dw Distance to the nearest wall 34 * @param[out] ubar Mean velocity at dw 35 * @param[out] cij Cholesky decomposition at dw 36 * @param[out] eps Turbulent dissipation at dw 37 * @param[out] lt Turbulent length scale at dw 38 * @param[in] stg_ctx STGShur14Context for the problem 39 */ 40 CEED_QFUNCTION_HELPER void InterpolateProfile(const CeedScalar dw, 41 CeedScalar ubar[3], CeedScalar cij[6], CeedScalar *eps, CeedScalar *lt, 42 const STGShur14Context stg_ctx) { 43 44 const CeedInt nprofs = stg_ctx->nprofs; 45 const CeedScalar *prof_dw = &stg_ctx->data[stg_ctx->offsets.prof_dw]; 46 const CeedScalar *prof_eps = &stg_ctx->data[stg_ctx->offsets.eps]; 47 const CeedScalar *prof_lt = &stg_ctx->data[stg_ctx->offsets.lt]; 48 const CeedScalar *prof_ubar = &stg_ctx->data[stg_ctx->offsets.ubar]; 49 const CeedScalar *prof_cij = &stg_ctx->data[stg_ctx->offsets.cij]; 50 CeedInt idx=-1; 51 52 for(CeedInt i=0; i<nprofs; i++) { 53 if (dw < prof_dw[i]) { 54 idx = i; 55 break; 56 } 57 } 58 59 if (idx > 0) { // y within the bounds of prof_dw 60 CeedScalar coeff = (dw - prof_dw[idx-1]) / (prof_dw[idx] - prof_dw[idx-1]); 61 62 //*INDENT-OFF* 63 ubar[0] = prof_ubar[0*nprofs+idx-1] + coeff*( prof_ubar[0*nprofs+idx] - prof_ubar[0*nprofs+idx-1] ); 64 ubar[1] = prof_ubar[1*nprofs+idx-1] + coeff*( prof_ubar[1*nprofs+idx] - prof_ubar[1*nprofs+idx-1] ); 65 ubar[2] = prof_ubar[2*nprofs+idx-1] + coeff*( prof_ubar[2*nprofs+idx] - prof_ubar[2*nprofs+idx-1] ); 66 cij[0] = prof_cij[0*nprofs+idx-1] + coeff*( prof_cij[0*nprofs+idx] - prof_cij[0*nprofs+idx-1] ); 67 cij[1] = prof_cij[1*nprofs+idx-1] + coeff*( prof_cij[1*nprofs+idx] - prof_cij[1*nprofs+idx-1] ); 68 cij[2] = prof_cij[2*nprofs+idx-1] + coeff*( prof_cij[2*nprofs+idx] - prof_cij[2*nprofs+idx-1] ); 69 cij[3] = prof_cij[3*nprofs+idx-1] + coeff*( prof_cij[3*nprofs+idx] - prof_cij[3*nprofs+idx-1] ); 70 cij[4] = prof_cij[4*nprofs+idx-1] + coeff*( prof_cij[4*nprofs+idx] - prof_cij[4*nprofs+idx-1] ); 71 cij[5] = prof_cij[5*nprofs+idx-1] + coeff*( prof_cij[5*nprofs+idx] - prof_cij[5*nprofs+idx-1] ); 72 *eps = prof_eps[idx-1] + coeff*( prof_eps[idx] - prof_eps[idx-1] ); 73 *lt = prof_lt[idx-1] + coeff*( prof_lt[idx] - prof_lt[idx-1] ); 74 //*INDENT-ON* 75 } else { // y outside bounds of prof_dw 76 ubar[0] = prof_ubar[1*nprofs-1]; 77 ubar[1] = prof_ubar[2*nprofs-1]; 78 ubar[2] = prof_ubar[3*nprofs-1]; 79 cij[0] = prof_cij[1*nprofs-1]; 80 cij[1] = prof_cij[2*nprofs-1]; 81 cij[2] = prof_cij[3*nprofs-1]; 82 cij[3] = prof_cij[4*nprofs-1]; 83 cij[4] = prof_cij[5*nprofs-1]; 84 cij[5] = prof_cij[6*nprofs-1]; 85 *eps = prof_eps[nprofs-1]; 86 *lt = prof_lt[nprofs-1]; 87 } 88 } 89 90 /* 91 * @brief Calculate spectrum coefficient, qn 92 * 93 * Calculates q_n at a given distance to the wall 94 * 95 * @param[in] kappa nth wavenumber 96 * @param[in] dkappa Difference between wavenumbers 97 * @param[in] keta Dissipation wavenumber 98 * @param[in] kcut Mesh-induced cutoff wavenumber 99 * @param[in] ke Energy-containing wavenumber 100 * @param[in] Ektot Total turbulent kinetic energy of spectrum 101 * @returns qn Spectrum coefficient 102 */ 103 CeedScalar CEED_QFUNCTION_HELPER(Calc_qn)(const CeedScalar kappa, 104 const CeedScalar dkappa, const CeedScalar keta, const CeedScalar kcut, 105 const CeedScalar ke, const CeedScalar Ektot_inv) { 106 const CeedScalar feta_x_fcut = exp(-Square(12*kappa/keta) 107 -Cube(4*Max(kappa - 0.9*kcut, 0)/kcut) ); 108 return pow(kappa/ke, 4.) * pow(1 + 2.4*Square(kappa/ke),-17./6) 109 *feta_x_fcut*dkappa * Ektot_inv; 110 } 111 112 // Calculate hmax, ke, keta, and kcut 113 void CEED_QFUNCTION_HELPER(SpectrumConstants)(const CeedScalar dw, 114 const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 115 const CeedScalar nu, CeedScalar *hmax, CeedScalar *ke, 116 CeedScalar *keta, CeedScalar *kcut) { 117 *hmax = Max( Max(h[0], h[1]), h[2]); 118 *ke = dw==0 ? 1e16 : 2*M_PI/Min(2*dw, 3*lt); 119 *keta = 2*M_PI*pow(Cube(nu)/eps, -0.25); 120 *kcut = M_PI/ Min( Max(Max(h[1], h[2]), 0.3*(*hmax)) + 0.1*dw, *hmax ); 121 } 122 123 /* 124 * @brief Calculate spectrum coefficients for STG 125 * 126 * Calculates q_n at a given distance to the wall 127 * 128 * @param[in] dw Distance to the nearest wall 129 * @param[in] eps Turbulent dissipation w/rt dw 130 * @param[in] lt Turbulent length scale w/rt dw 131 * @param[in] h Element lengths in coordinate directions 132 * @param[in] nu Dynamic Viscosity; 133 * @param[in] stg_ctx STGShur14Context for the problem 134 * @param[out] qn Spectrum coefficients, [nmodes] 135 */ 136 void CEED_QFUNCTION_HELPER(CalcSpectrum)(const CeedScalar dw, 137 const CeedScalar eps, const CeedScalar lt, const CeedScalar h[3], 138 const CeedScalar nu, CeedScalar qn[], const STGShur14Context stg_ctx) { 139 140 const CeedInt nmodes = stg_ctx->nmodes; 141 const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 142 CeedScalar hmax, ke, keta, kcut, Ektot=0.0; 143 SpectrumConstants(dw, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 144 145 for(CeedInt n=0; n<nmodes; n++) { 146 const CeedScalar dkappa = n==0 ? kappa[0] : kappa[n] - kappa[n-1]; 147 qn[n] = Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 148 Ektot += qn[n]; 149 } 150 151 if (Ektot == 0) return; 152 for(CeedInt n=0; n<nmodes; n++) qn[n] /= Ektot; 153 } 154 155 /****************************************************** 156 * @brief Calculate u(x,t) for STG inflow condition 157 * 158 * @param[in] X Location to evaluate u(X,t) 159 * @param[in] t Time to evaluate u(X,t) 160 * @param[in] ubar Mean velocity at X 161 * @param[in] cij Cholesky decomposition at X 162 * @param[in] qn Wavemode amplitudes at X, [nmodes] 163 * @param[out] u Velocity at X and t 164 * @param[in] stg_ctx STGShur14Context for the problem 165 */ 166 void CEED_QFUNCTION_HELPER(STGShur14_Calc)(const CeedScalar X[3], 167 const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 168 const CeedScalar qn[], CeedScalar u[3], 169 const STGShur14Context stg_ctx) { 170 171 //*INDENT-OFF* 172 const CeedInt nmodes = stg_ctx->nmodes; 173 const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 174 const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 175 const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 176 const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 177 //*INDENT-ON* 178 CeedScalar xdotd, vp[3] = {0.}; 179 CeedScalar xhat[] = {0., X[1], X[2]}; 180 181 CeedPragmaSIMD 182 for(CeedInt n=0; n<nmodes; n++) { 183 xhat[0] = (X[0] - stg_ctx->u0*t)*Max(2*kappa[0]/kappa[n], 0.1); 184 xdotd = 0.; 185 for(CeedInt i=0; i<3; i++) xdotd += d[i*nmodes+n]*xhat[i]; 186 const CeedScalar cos_kxdp = cos(kappa[n]*xdotd + phi[n]); 187 vp[0] += sqrt(qn[n])*sigma[0*nmodes+n] * cos_kxdp; 188 vp[1] += sqrt(qn[n])*sigma[1*nmodes+n] * cos_kxdp; 189 vp[2] += sqrt(qn[n])*sigma[2*nmodes+n] * cos_kxdp; 190 } 191 for(CeedInt i=0; i<3; i++) vp[i] *= 2*sqrt(1.5); 192 193 u[0] = ubar[0] + cij[0]*vp[0]; 194 u[1] = ubar[1] + cij[3]*vp[0] + cij[1]*vp[1]; 195 u[2] = ubar[2] + cij[4]*vp[0] + cij[5]*vp[1] + cij[2]*vp[2]; 196 } 197 198 /****************************************************** 199 * @brief Calculate u(x,t) for STG inflow condition 200 * 201 * @param[in] X Location to evaluate u(X,t) 202 * @param[in] t Time to evaluate u(X,t) 203 * @param[in] ubar Mean velocity at X 204 * @param[in] cij Cholesky decomposition at X 205 * @param[in] qn Wavemode amplitudes at X, [nmodes] 206 * @param[out] u Velocity at X and t 207 * @param[in] stg_ctx STGShur14Context for the problem 208 */ 209 void CEED_QFUNCTION_HELPER(STGShur14_Calc_PrecompEktot)(const CeedScalar X[3], 210 const CeedScalar t, const CeedScalar ubar[3], const CeedScalar cij[6], 211 const CeedScalar Ektot, const CeedScalar h[3], const CeedScalar dw, 212 const CeedScalar eps, const CeedScalar lt, const CeedScalar nu, CeedScalar u[3], 213 const STGShur14Context stg_ctx) { 214 215 //*INDENT-OFF* 216 const CeedInt nmodes = stg_ctx->nmodes; 217 const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 218 const CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 219 const CeedScalar *sigma = &stg_ctx->data[stg_ctx->offsets.sigma]; 220 const CeedScalar *d = &stg_ctx->data[stg_ctx->offsets.d]; 221 //*INDENT-ON* 222 CeedScalar hmax, ke, keta, kcut; 223 SpectrumConstants(dw, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 224 CeedScalar xdotd, vp[3] = {0.}; 225 CeedScalar xhat[] = {0., X[1], X[2]}; 226 227 CeedPragmaSIMD 228 for(CeedInt n=0; n<nmodes; n++) { 229 xhat[0] = (X[0] - stg_ctx->u0*t)*Max(2*kappa[0]/kappa[n], 0.1); 230 xdotd = 0.; 231 for(CeedInt i=0; i<3; i++) xdotd += d[i*nmodes+n]*xhat[i]; 232 const CeedScalar cos_kxdp = cos(kappa[n]*xdotd + phi[n]); 233 const CeedScalar dkappa = n==0 ? kappa[0] : kappa[n] - kappa[n-1]; 234 const CeedScalar qn = Calc_qn(kappa[n], dkappa, keta, kcut, ke, Ektot); 235 vp[0] += sqrt(qn)*sigma[0*nmodes+n] * cos_kxdp; 236 vp[1] += sqrt(qn)*sigma[1*nmodes+n] * cos_kxdp; 237 vp[2] += sqrt(qn)*sigma[2*nmodes+n] * cos_kxdp; 238 } 239 for(CeedInt i=0; i<3; i++) vp[i] *= 2*sqrt(1.5); 240 241 u[0] = ubar[0] + cij[0]*vp[0]; 242 u[1] = ubar[1] + cij[3]*vp[0] + cij[1]*vp[1]; 243 u[2] = ubar[2] + cij[4]*vp[0] + cij[5]*vp[1] + cij[2]*vp[2]; 244 } 245 246 // Create preprocessed input for the stg calculation 247 // 248 // stg_data[0] = 1 / Ektot (inverse of total spectrum energy) 249 CEED_QFUNCTION(Preprocess_STGShur14)(void *ctx, CeedInt Q, 250 const CeedScalar *const *in, CeedScalar *const *out) { 251 //*INDENT-OFF* 252 const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0], 253 (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[1]; 254 255 CeedScalar (*stg_data) = (CeedScalar(*)) out[0]; 256 257 //*INDENT-ON* 258 CeedScalar ubar[3], cij[6], eps, lt; 259 const STGShur14Context stg_ctx = (STGShur14Context) ctx; 260 const CeedScalar dx = stg_ctx->dx; 261 const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 262 const CeedScalar theta0 = stg_ctx->theta0; 263 const CeedScalar P0 = stg_ctx->P0; 264 const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 265 const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 266 const CeedScalar Rd = cp - cv; 267 const CeedScalar rho = P0 / (Rd * theta0); 268 const CeedScalar nu = mu / rho; 269 270 const CeedInt nmodes = stg_ctx->nmodes; 271 const CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 272 CeedScalar hmax, ke, keta, kcut; 273 274 CeedPragmaSIMD 275 for(CeedInt i=0; i<Q; i++) { 276 const CeedScalar dw = x[1][i]; 277 const CeedScalar dXdx[2][3] = { 278 {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 279 {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 280 }; 281 282 CeedScalar h[3]; 283 h[0] = dx; 284 for (CeedInt j=1; j<3; j++) 285 h[j] = 2/sqrt(dXdx[0][j]*dXdx[0][j] + dXdx[1][j]*dXdx[1][j]); 286 287 InterpolateProfile(dw, ubar, cij, &eps, <, stg_ctx); 288 SpectrumConstants(dw, eps, lt, h, nu, &hmax, &ke, &keta, &kcut); 289 290 // Calculate total TKE per spectrum 291 stg_data[i] = 0.; 292 CeedPragmaSIMD 293 for(CeedInt n=0; n<nmodes; n++) { 294 const CeedScalar dkappa = n==0 ? kappa[0] : kappa[n] - kappa[n-1]; 295 stg_data[i] += Calc_qn(kappa[n], dkappa, keta, kcut, ke, 1.0); 296 } 297 stg_data[i] = 1/stg_data[i]; 298 } 299 return 0; 300 } 301 302 // Extrude the STGInflow profile through out the domain for an initial condition 303 CEED_QFUNCTION(ICsSTG)(void *ctx, CeedInt Q, 304 const CeedScalar *const *in, CeedScalar *const *out) { 305 // Inputs 306 const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 307 308 // Outputs 309 CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 310 311 const STGShur14Context stg_ctx = (STGShur14Context) ctx; 312 CeedScalar u[3], cij[6], eps, lt; 313 const CeedScalar theta0 = stg_ctx->theta0; 314 const CeedScalar P0 = stg_ctx->P0; 315 const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 316 const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 317 const CeedScalar Rd = cp - cv; 318 const CeedScalar rho = P0 / (Rd * theta0); 319 320 CeedPragmaSIMD 321 for(CeedInt i=0; i<Q; i++) { 322 InterpolateProfile(X[1][i], u, cij, &eps, <, stg_ctx); 323 324 if (stg_ctx->newtonian_ctx.is_primitive) { 325 q0[0][i] = P0; 326 q0[1][i] = u[0]; 327 q0[2][i] = u[1]; 328 q0[3][i] = u[2]; 329 q0[4][i] = theta0; 330 } else { 331 q0[0][i] = rho; 332 q0[1][i] = u[0] * rho; 333 q0[2][i] = u[1] * rho; 334 q0[3][i] = u[2] * rho; 335 q0[4][i] = rho * (0.5 * Dot3(u, u) + cv * theta0); 336 } 337 } // End of Quadrature Point Loop 338 return 0; 339 } 340 341 /******************************************************************** 342 * @brief QFunction to calculate the inflow boundary condition 343 * 344 * This will loop through quadrature points, calculate the wavemode amplitudes 345 * at each location, then calculate the actual velocity. 346 */ 347 CEED_QFUNCTION(STGShur14_Inflow)(void *ctx, CeedInt Q, 348 const CeedScalar *const *in, 349 CeedScalar *const *out) { 350 351 //*INDENT-OFF* 352 const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0], 353 (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[2], 354 (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[3]; 355 356 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[0], 357 (*jac_data_sur)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[1]; 358 359 //*INDENT-ON* 360 361 const STGShur14Context stg_ctx = (STGShur14Context) ctx; 362 CeedScalar qn[STG_NMODES_MAX], u[3], ubar[3], cij[6], eps, lt; 363 const bool is_implicit = stg_ctx->is_implicit; 364 const bool mean_only = stg_ctx->mean_only; 365 const bool prescribe_T = stg_ctx->prescribe_T; 366 const CeedScalar dx = stg_ctx->dx; 367 const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 368 const CeedScalar time = stg_ctx->time; 369 const CeedScalar theta0 = stg_ctx->theta0; 370 const CeedScalar P0 = stg_ctx->P0; 371 const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 372 const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 373 const CeedScalar Rd = cp - cv; 374 const CeedScalar gamma = cp/cv; 375 376 CeedPragmaSIMD 377 for(CeedInt i=0; i<Q; i++) { 378 const CeedScalar rho = prescribe_T ? q[0][i] : P0 / (Rd * theta0); 379 const CeedScalar x[] = { X[0][i], X[1][i], X[2][i] }; 380 const CeedScalar dXdx[2][3] = { 381 {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 382 {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 383 }; 384 385 CeedScalar h[3]; 386 h[0] = dx; 387 for (CeedInt j=1; j<3; j++) 388 h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 389 390 InterpolateProfile(X[1][i], ubar, cij, &eps, <, stg_ctx); 391 if (!mean_only) { 392 CalcSpectrum(X[1][i], eps, lt, h, mu/rho, qn, stg_ctx); 393 STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 394 } else { 395 for (CeedInt j=0; j<3; j++) u[j] = ubar[j]; 396 } 397 398 const CeedScalar E_kinetic = .5 * rho * Dot3(u, u); 399 CeedScalar E_internal, P; 400 if (prescribe_T) { 401 // Temperature is being set weakly (theta0) and for constant cv this sets E_internal 402 E_internal = rho * cv * theta0; 403 // Find pressure using 404 P = rho * Rd * theta0; // interior rho with exterior T 405 } else { 406 E_internal = q[4][i] - E_kinetic; // uses prescribed rho and u, E from solution 407 P = E_internal * (gamma - 1.); 408 } 409 410 const CeedScalar wdetJb = (is_implicit ? -1. : 1.) * q_data_sur[0][i]; 411 // ---- Normal vect 412 const CeedScalar norm[3] = {q_data_sur[1][i], 413 q_data_sur[2][i], 414 q_data_sur[3][i] 415 }; 416 417 const CeedScalar E = E_internal + E_kinetic; 418 419 // Velocity normal to the boundary 420 const CeedScalar u_normal = Dot3(norm, u); 421 422 // The Physics 423 // Zero v so all future terms can safely sum into it 424 for (CeedInt j=0; j<5; j++) v[j][i] = 0.; 425 426 // The Physics 427 // -- Density 428 v[0][i] -= wdetJb * rho * u_normal; 429 430 // -- Momentum 431 for (CeedInt j=0; j<3; j++) 432 v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + 433 norm[j] * P); 434 435 // -- Total Energy Density 436 v[4][i] -= wdetJb * u_normal * (E + P); 437 438 jac_data_sur[0][i] = rho; 439 jac_data_sur[1][i] = u[0]; 440 jac_data_sur[2][i] = u[1]; 441 jac_data_sur[3][i] = u[2]; 442 jac_data_sur[4][i] = E; 443 for (int j=0; j<6; j++) jac_data_sur[5+j][i] = 0.; 444 } 445 return 0; 446 } 447 448 CEED_QFUNCTION(STGShur14_Inflow_Jacobian)(void *ctx, CeedInt Q, 449 const CeedScalar *const *in, 450 CeedScalar *const *out) { 451 // *INDENT-OFF* 452 // Inputs 453 const CeedScalar (*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 454 (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 455 (*jac_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[4]; 456 // Outputs 457 CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 458 // *INDENT-ON* 459 const STGShur14Context stg_ctx = (STGShur14Context)ctx; 460 const bool implicit = stg_ctx->is_implicit; 461 const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 462 const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 463 const CeedScalar Rd = cp - cv; 464 const CeedScalar gamma = cp/cv; 465 466 const CeedScalar theta0 = stg_ctx->theta0; 467 const bool prescribe_T = stg_ctx->prescribe_T; 468 469 CeedPragmaSIMD 470 // Quadrature Point Loop 471 for (CeedInt i=0; i<Q; i++) { 472 // Setup 473 // Setup 474 // -- Interp-to-Interp q_data 475 // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 476 // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 477 // We can effect this by swapping the sign on this weight 478 const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 479 480 // Calculate inflow values 481 CeedScalar velocity[3]; 482 for (CeedInt j=0; j<3; j++) velocity[j] = jac_data_sur[5+j][i]; 483 484 // enabling user to choose between weak T and weak rho inflow 485 CeedScalar drho, dE, dP; 486 if (prescribe_T) { 487 // rho should be from the current solution 488 drho = dq[0][i]; 489 CeedScalar dE_internal = drho * cv * theta0; 490 CeedScalar dE_kinetic = .5 * drho * Dot3(velocity, velocity); 491 dE = dE_internal + dE_kinetic; 492 dP = drho * Rd * theta0; // interior rho with exterior T 493 } else { // rho specified, E_internal from solution 494 drho = 0; 495 dE = dq[4][i]; 496 dP = dE * (gamma - 1.); 497 } 498 const CeedScalar norm[3] = {q_data_sur[1][i], 499 q_data_sur[2][i], 500 q_data_sur[3][i] 501 }; 502 503 const CeedScalar u_normal = Dot3(norm, velocity); 504 505 v[0][i] = - wdetJb * drho * u_normal; 506 for (int j=0; j<3; j++) 507 v[j+1][i] = -wdetJb * (drho * u_normal * velocity[j] + norm[j] * dP); 508 v[4][i] = - wdetJb * u_normal * (dE + dP); 509 } // End Quadrature Point Loop 510 return 0; 511 } 512 513 /******************************************************************** 514 * @brief QFunction to calculate the strongly enforce inflow BC 515 * 516 * This QF is for the strong application of STG via libCEED (rather than 517 * through the native PETSc `DMAddBoundary` -> `bcFunc` method. 518 */ 519 CEED_QFUNCTION(STGShur14_Inflow_StrongQF)(void *ctx, CeedInt Q, 520 const CeedScalar *const *in, CeedScalar *const *out) { 521 522 //*INDENT-OFF* 523 const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[0], 524 (*coords)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[1], 525 (*scale) = (const CeedScalar(*)) in[2], 526 (*stg_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA]) in[3]; 527 528 CeedScalar(*bcval)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[0]; 529 //*INDENT-ON* 530 531 const STGShur14Context stg_ctx = (STGShur14Context) ctx; 532 CeedScalar u[3], ubar[3], cij[6], eps, lt; 533 const bool mean_only = stg_ctx->mean_only; 534 const CeedScalar dx = stg_ctx->dx; 535 const CeedScalar mu = stg_ctx->newtonian_ctx.mu; 536 const CeedScalar time = stg_ctx->time; 537 const CeedScalar theta0 = stg_ctx->theta0; 538 const CeedScalar P0 = stg_ctx->P0; 539 const CeedScalar cv = stg_ctx->newtonian_ctx.cv; 540 const CeedScalar cp = stg_ctx->newtonian_ctx.cp; 541 const CeedScalar Rd = cp - cv; 542 const CeedScalar rho = P0 / (Rd * theta0); 543 544 CeedPragmaSIMD 545 for(CeedInt i=0; i<Q; i++) { 546 const CeedScalar x[] = { coords[0][i], coords[1][i], coords[2][i] }; 547 const CeedScalar dXdx[2][3] = { 548 {q_data_sur[4][i], q_data_sur[5][i], q_data_sur[6][i]}, 549 {q_data_sur[7][i], q_data_sur[8][i], q_data_sur[9][i]} 550 }; 551 552 CeedScalar h[3]; 553 h[0] = dx; 554 for (CeedInt j=1; j<3; j++) 555 h[j] = 2/sqrt(Square(dXdx[0][j]) + Square(dXdx[1][j])); 556 557 InterpolateProfile(coords[1][i], ubar, cij, &eps, <, stg_ctx); 558 if (!mean_only) { 559 if (1) { 560 STGShur14_Calc_PrecompEktot(x, time, ubar, cij, stg_data[0][i], 561 h, x[1], eps, lt, mu/rho, u, stg_ctx); 562 } else { // Original way 563 CeedScalar qn[STG_NMODES_MAX]; 564 CalcSpectrum(coords[1][i], eps, lt, h, mu/rho, qn, stg_ctx); 565 STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 566 } 567 } else { 568 for (CeedInt j=0; j<3; j++) u[j] = ubar[j]; 569 } 570 571 if (stg_ctx->newtonian_ctx.is_primitive) { 572 bcval[0][i] = 0; 573 bcval[1][i] = scale[i] * u[0]; 574 bcval[2][i] = scale[i] * u[1]; 575 bcval[3][i] = scale[i] * u[2]; 576 bcval[4][i] = scale[i] * theta0; 577 578 } else { 579 bcval[0][i] = scale[i] * rho; 580 bcval[1][i] = scale[i] * rho * u[0]; 581 bcval[2][i] = scale[i] * rho * u[1]; 582 bcval[3][i] = scale[i] * rho * u[2]; 583 bcval[4][i] = 0.; 584 } 585 } 586 return 0; 587 } 588 589 #endif // stg_shur14_h 590