Function for computing the implicit contribution to the vertical diffusion of momentum, over terrain.
This function (explicitly instantiated below) handles staggering in x, y, or z through the template parameter, stagdir. NOTE: implicit diffusion of w has remains an experimental feature and has not been tested yet with terrain.
239 BL_PROFILE_VAR(
"ImplicitDiffForMom_T()",ImplicitDiffForMom_T);
256 constexpr
int ioff = (stagdir == 0) ? 1 : 0;
257 constexpr
int joff = (stagdir == 1) ? 1 : 0;
260 int ilo = bx.smallEnd(0);
261 int ihi = bx.bigEnd(0);
262 int jlo = bx.smallEnd(1);
263 int jhi = bx.bigEnd(1);
264 int klo = bx.smallEnd(2);
265 int khi = bx.bigEnd(2);
266 amrex::ignore_unused(ilo, ihi, jlo, jhi);
270 amrex::FArrayBox RHS_fab, soln_fab, coeffG_fab;
271 RHS_fab.resize(bx,1, amrex::The_Async_Arena());
272 soln_fab.resize(bx,1, amrex::The_Async_Arena());
273 coeffG_fab.resize(bx,1, amrex::The_Async_Arena());
274 auto const& RHS_a = RHS_fab.array();
275 auto const& soln_a = soln_fab.array();
276 auto const& coeffG_a = coeffG_fab.array();
278 Real dz_inv = cellSizeInv[2];
286 amrex::ignore_unused(foextrap_on_zhi);
289 "Unexpected lower BC used with implicit vertical diffusion");
291 "Unexpected upper BC used with implicit vertical diffusion");
292 if (stagdir < 2 && (ext_dir_on_zlo || ext_dir_on_zhi)) {
293 amrex::Warning(
"No-slip walls have not been fully tested");
296 Real Fact = implicit_fac * dt * dz_inv;
299 ParallelFor(makeSlab(bx,2,0), [=] AMREX_GPU_DEVICE (
int i,
int j,
int)
302 for (
int j(jlo); j<=jhi; ++j) {
303 for (
int i(ilo); i<=ihi; ++i) {
344 Real rhoface, rhoAlpha_lo, rhoAlpha_hi;
345 Real detJface, met_h_zeta_lo, met_h_zeta_hi;
346 Real a_tmp, b_tmp, c_tmp, inv_b2_tmp;
348 detJface =
myhalf * (detJ(i,j,klo) + detJ(i-ioff,j-joff,klo));
351 cell_data, mu_turb, mu_eff,
360 c_tmp = -Fact * gfac * rhoAlpha_hi * dz_inv / met_h_zeta_hi;
362 RHS_a(i,j,klo) = detJface * face_data(i,j,klo);
365 if (ext_dir_on_zlo) {
366 RHS_a(i,j,klo) += Fact * gfac * (tau_corr(i,j,klo+1) - tau_corr(i,j,klo));
369 RHS_a(i,j,klo) =
zero;
371 a_tmp = -
two * Fact * rhoAlpha_lo * dz_inv / met_h_zeta_lo;
372 RHS_a(i,j,klo) +=
two * rhoAlpha_lo * face_data(i,j,klo-1) * dz_inv * dz_inv / met_h_zeta_lo;
374 }
else if (use_SurfLayer) {
376 RHS_a(i,j,klo) += Fact * gfac * (tau_corr(i,j,klo+1) - tau(i,j,klo));
377 RHS_a(i,j,klo) += Fact * tau(i,j,klo);
380 b_tmp = detJface * rhoface - a_tmp - c_tmp;
383 RHS_a(i,j,klo) /= b_tmp;
384 coeffG_a(i,j,klo) = c_tmp / b_tmp;
389 for (
int k(klo+1); k <
khi; k++) {
390 detJface =
myhalf * (detJ(i,j,k) + detJ(i-ioff,j-joff,k));
393 cell_data, mu_turb, mu_eff,
401 a_tmp = -Fact * rhoAlpha_lo * dz_inv / met_h_zeta_lo;
402 c_tmp = -Fact * rhoAlpha_hi * dz_inv / met_h_zeta_hi;
403 b_tmp = detJface * rhoface - a_tmp - c_tmp;
404 inv_b2_tmp =
one / (b_tmp - a_tmp * coeffG_a(i,j,k-1));
406 RHS_a(i,j,k) = detJface * face_data(i,j,k);
407 RHS_a(i,j,k) += Fact * gfac * (tau_corr(i,j,k+1) - tau_corr(i,j,k));
409 RHS_a(i,j,k) = (RHS_a(i,j,k) - a_tmp * RHS_a(i,j,k-1)) * inv_b2_tmp;
410 coeffG_a(i,j,k) = c_tmp * inv_b2_tmp;
416 detJface =
myhalf * (detJ(i,j,
khi) + detJ(i-ioff,j-joff,
khi));
419 cell_data, mu_turb, mu_eff,
427 a_tmp = -Fact * gfac * rhoAlpha_lo * dz_inv / met_h_zeta_lo;
430 RHS_a(i,j,
khi) = detJface * face_data(i,j,
khi);
431 RHS_a(i,j,
khi) += Fact * gfac * (tau_corr(i,j,
khi+1) - tau_corr(i,j,
khi));
434 if (ext_dir_on_zhi) {
439 c_tmp = -
two * Fact * rhoAlpha_hi * dz_inv / met_h_zeta_hi;
440 RHS_a(i,j,
khi) +=
two * rhoAlpha_hi * face_data(i,j,
khi+1) * dz_inv * dz_inv / met_h_zeta_hi;
444 b_tmp = detJface * rhoface - a_tmp - c_tmp;
445 inv_b2_tmp =
one / (b_tmp - a_tmp * coeffG_a(i,j,
khi-1));
448 soln_a(i,j,
khi) = (RHS_a(i,j,
khi) - a_tmp * RHS_a(i,j,
khi-1)) * inv_b2_tmp;
453 for (
int k(
khi-1); k>=klo; --k) {
454 soln_a(i,j,k) = RHS_a(i,j,k) - coeffG_a(i,j,k) * soln_a(i,j,k+1);
459 for (
int k(klo); k<=
khi; ++k) {
461 face_data(i,j,k) = rhoface * soln_a(i,j,k);
constexpr amrex::Real three
Definition: ERF_Constants.H:9
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void getRhoAlphaForFaces(int i, int j, int k, int ioff, int joff, amrex::Real &rhoAlpha_lo, amrex::Real &rhoAlpha_hi, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &mu_turb, const amrex::Real mu_eff, bool l_consA, bool l_turb)
Definition: ERF_GetRhoAlphaForFaces.H:5
#define Rho_comp
Definition: ERF_IndexDefines.H:36
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+myhalf) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=one) { Real dT=T_pert *(std::cos(PI *L)+one)/two;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })
bool l_turb
Definition: ERF_SetupVertDiff.H:9
bool l_consA
Definition: ERF_SetupVertDiff.H:8
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtKface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:184
AMREX_ASSERT_WITH_MESSAGE(wbar_cutoff_min > wbar_cutoff_max, "ERROR: wbar_cutoff_min < wbar_cutoff_max")
@ xvel_bc
Definition: ERF_IndexDefines.H:87
@ foextrap
Definition: ERF_IndexDefines.H:226
@ ext_dir
Definition: ERF_IndexDefines.H:227
@ ext_dir_prim
Definition: ERF_IndexDefines.H:229
Definition: ERF_DiffStruct.H:19
amrex::Real rho0_trans
Definition: ERF_DiffStruct.H:91
MolecDiffType molec_diff_type
Definition: ERF_DiffStruct.H:84
amrex::Real dynamic_viscosity
Definition: ERF_DiffStruct.H:96
DiffChoice diffChoice
Definition: ERF_DataStruct.H:1088
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:1091
Definition: ERF_TurbStruct.H:42
bool use_kturb
Definition: ERF_TurbStruct.H:424