ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_ComputeTurbulentViscosity.cpp File Reference
#include "ERF_SurfaceLayer.H"
#include "ERF_EddyViscosity.H"
#include "ERF_Diffusion.H"
#include "ERF_PBLModels.H"
#include "ERF_TileNoZ.H"
#include "ERF_TerrainMetrics.H"
#include "ERF_MoistUtils.H"
#include "ERF_RichardsonNumber.H"
Include dependency graph for ERF_ComputeTurbulentViscosity.cpp:

Functions

void ComputeTurbulentViscosityLES (Vector< std::unique_ptr< MultiFab >> &Tau_lev, const MultiFab &cons_in, MultiFab &eddyViscosity, MultiFab &Hfx1, MultiFab &Hfx2, MultiFab &Hfx3, MultiFab &Diss, const Geometry &geom, bool use_terrain_fitted_coords, Vector< std::unique_ptr< MultiFab >> &mapfac, const std::unique_ptr< MultiFab > &z_phys_nd, const TurbChoice &turbChoice, const Real const_grav, std::unique_ptr< SurfaceLayer > &, const MoistureComponentIndices &moisture_indices, const MultiFab *xvel, const MultiFab *yvel)
 
void ComputeTurbulentViscosityLES_EB (Vector< std::unique_ptr< MultiFab >> &Tau_lev, const MultiFab &cons_in, MultiFab &eddyViscosity, MultiFab &Hfx1, MultiFab &Hfx2, MultiFab &Hfx3, const Geometry &geom, Vector< std::unique_ptr< MultiFab >> &mapfac, const TurbChoice &turbChoice, const Real const_grav, [[maybe_unused]] const SolverChoice &solverChoice, std::unique_ptr< SurfaceLayer > &, const MoistureComponentIndices &moisture_indices, const eb_ &ebfact, const MultiFab *xvel, const MultiFab *yvel)
 
void ComputeTurbulentViscosityRANS (Vector< std::unique_ptr< MultiFab >> &, const MultiFab &cons_in, const MultiFab &wdist, MultiFab &eddyViscosity, MultiFab &Hfx1, MultiFab &Hfx2, MultiFab &Hfx3, MultiFab &Diss, const Geometry &geom, bool use_terrain_fitted_coords, Vector< std::unique_ptr< MultiFab >> &, const std::unique_ptr< MultiFab > &z_phys_nd, const TurbChoice &turbChoice, const Real const_grav, std::unique_ptr< SurfaceLayer > &SurfLayer, const MultiFab *z_0)
 
void ComputeTurbulentViscosity (Real dt, const MultiFab &xvel, const MultiFab &yvel, Vector< std::unique_ptr< MultiFab >> &Tau_lev, MultiFab &cons_in, const MultiFab &wdist, MultiFab &eddyViscosity, MultiFab &Hfx1, MultiFab &Hfx2, MultiFab &Hfx3, MultiFab &Diss, const Geometry &geom, Vector< std::unique_ptr< MultiFab >> &mapfac, const std::unique_ptr< MultiFab > &z_phys_nd, const SolverChoice &solverChoice, std::unique_ptr< SurfaceLayer > &SurfLayer, const MultiFab *z_0, const bool &use_terrain_fitted_coords, const bool &use_moisture, int level, const BCRec *bc_ptr, const eb_ &ebfact, bool vert_only)
 

Function Documentation

◆ ComputeTurbulentViscosity()

void ComputeTurbulentViscosity ( Real  dt,
const MultiFab &  xvel,
const MultiFab &  yvel,
Vector< std::unique_ptr< MultiFab >> &  Tau_lev,
MultiFab &  cons_in,
const MultiFab &  wdist,
MultiFab &  eddyViscosity,
MultiFab &  Hfx1,
MultiFab &  Hfx2,
MultiFab &  Hfx3,
MultiFab &  Diss,
const Geometry &  geom,
Vector< std::unique_ptr< MultiFab >> &  mapfac,
const std::unique_ptr< MultiFab > &  z_phys_nd,
const SolverChoice solverChoice,
std::unique_ptr< SurfaceLayer > &  SurfLayer,
const MultiFab *  z_0,
const bool &  use_terrain_fitted_coords,
const bool &  use_moisture,
int  level,
const BCRec *  bc_ptr,
const eb_ ebfact,
bool  vert_only 
)

Wrapper to compute turbulent viscosity with LES or PBL.

Parameters
[in]xvelvelocity in x-dir
[in]yvelvelocity in y-dir
[in]Tau_levstrain at this level
[in]cons_incell center conserved quantities
[out]eddyViscosityturbulent viscosity
[in]Hfx1heat flux in x-dir
[in]Hfx2heat flux in y-dir
[in]Hfx3heat flux in z-dir
[in]Dissdissipation of turbulent kinetic energy
[in]geomproblem geometry
[in]mapfacmap factors
[in]turbChoicecontainer with turbulence parameters
[in]mostpointer to Monin-Obukhov class if instantiated
[in]vert_onlyflag for vertical components of eddyViscosity
807 {
808  BL_PROFILE_VAR("ComputeTurbulentViscosity()",ComputeTurbulentViscosity);
809  //
810  // In LES mode, the turbulent viscosity is isotropic (unless mix_isotropic is set to false), so
811  // the LES model sets both horizontal and vertical viscosities
812  //
813  // In PBL mode, the primary purpose of the PBL model is to control vertical transport, so the PBL model sets the vertical viscosity.
814  // Optionally, the PBL model can be run in conjunction with an LES model that sets the horizontal viscosity
815  // (this isn’t truly LES, but the model form is the same as Smagorinsky).
816  //
817  // ComputeTurbulentViscosityLES populates the LES viscosity for both horizontal and vertical components.
818  // ComputeTurbulentViscosityPBL computes the PBL viscosity just for the vertical component.
819  //
820 
821  TurbChoice turbChoice = solverChoice.turbChoice[level];
822  const Real const_grav = solverChoice.gravity;
823 
824  if (!SurfLayer) {
825  AMREX_ALWAYS_ASSERT(!vert_only);
826  }
827 
828  bool impose_phys_bcs = false;
829 
830  if (turbChoice.les_type != LESType::None) {
831  impose_phys_bcs = true;
832  if (solverChoice.terrain_type == TerrainType::EB) {
834  cons_in, eddyViscosity,
835  Hfx1, Hfx2, Hfx3,
836  geom,
837  mapfac, turbChoice, const_grav,
838  solverChoice,
839  SurfLayer, solverChoice.moisture_indices,
840  ebfact, &xvel, &yvel);
841  } else {
843  cons_in, eddyViscosity,
844  Hfx1, Hfx2, Hfx3, Diss,
845  geom, use_terrain_fitted_coords,
846  mapfac, z_phys_nd, turbChoice, const_grav,
847  SurfLayer, solverChoice.moisture_indices,
848  &xvel, &yvel);
849  }
850  }
851 
852  if (turbChoice.rans_type != RANSType::None) {
853  impose_phys_bcs = true;
855  cons_in, wdist,
856  eddyViscosity,
857  Hfx1, Hfx2, Hfx3, Diss,
858  geom, use_terrain_fitted_coords,
859  mapfac, z_phys_nd, turbChoice, const_grav,
860  SurfLayer, z_0);
861  }
862 
863  if (turbChoice.pbl_type == PBLType::MYJ) {
864  ComputeDiffusivityMYJ(dt, xvel, yvel, cons_in, eddyViscosity,
865  geom, turbChoice, SurfLayer,
866  use_terrain_fitted_coords, use_moisture,
867  level, bc_ptr, vert_only, z_phys_nd,
868  solverChoice.moisture_indices);
869  } else if (turbChoice.pbl_type == PBLType::MYNN25) {
870  ComputeDiffusivityMYNN25(xvel, yvel, cons_in, eddyViscosity,
871  geom, turbChoice, SurfLayer,
872  use_terrain_fitted_coords, use_moisture,
873  level, bc_ptr, vert_only, z_phys_nd,
874  solverChoice.moisture_indices);
875  } else if (turbChoice.pbl_type == PBLType::MYNNEDMF) {
876  ComputeDiffusivityMYNNEDMF(xvel, yvel, cons_in, eddyViscosity,
877  geom, turbChoice, SurfLayer,
878  use_terrain_fitted_coords, use_moisture,
879  level, bc_ptr, vert_only, z_phys_nd,
880  solverChoice.moisture_indices);
881  } else if (turbChoice.pbl_type == PBLType::YSU) {
882  ComputeDiffusivityYSU(xvel, yvel, cons_in, eddyViscosity,
883  geom, turbChoice, SurfLayer,
884  use_terrain_fitted_coords, use_moisture,
885  level, bc_ptr, vert_only, z_phys_nd,
886  solverChoice.moisture_indices);
887  } else if (turbChoice.pbl_type == PBLType::MRF) {
888  ComputeDiffusivityMRF(xvel, yvel, cons_in, eddyViscosity,
889  geom, turbChoice, SurfLayer,
890  use_terrain_fitted_coords, use_moisture,
891  level, bc_ptr, vert_only, z_phys_nd,
892  solverChoice.moisture_indices);
893  } else if (turbChoice.pbl_type == PBLType::SHOC) {
894  // NOTE: Nothing to do here. The SHOC class handles setting the vertical
895  // components of eddyDiffs in slow RHS pre.
896  }
897 
898  //
899  // At all levels we need to fill values outside the physical boundary for the LES coeffs.
900  // In addition, for all cases, if at level > 0, we want to fill fine ghost cell values that
901  // overlie coarse grid cells (and that are not in another fine valid region) with
902  // extrapolated values from the interior, rather than interpolating from the coarser level,
903  // since we may be using a different turbulence model there.
904  //
905  // Note: here "covered" refers to "covered by valid region of another grid at this level"
906  // Note: here "physbnd" refers to "cells outside the domain if not periodic"
907  // Note: here "interior" refers to "valid cells, i.e. inside 'my' grid"
908  //
909  {
910  int is_covered = 0;
911  int is_notcovered = 1;
912  int is_physbnd = 2;
913  int is_interior = 3;
914  iMultiFab cc_mask(eddyViscosity.boxArray(),eddyViscosity.DistributionMap(),1,1);
915  cc_mask.BuildMask(geom.Domain(), geom.periodicity(), is_covered, is_notcovered, is_physbnd, is_interior);
916 
917  Box domain = geom.Domain();
918  for (int i = 0; i < AMREX_SPACEDIM; ++i) {
919  if (geom.isPeriodic(i)) {
920  domain.grow(i,1);
921  }
922  }
923 
924  eddyViscosity.FillBoundary(geom.periodicity());
925 
926  int ncomp = eddyViscosity.nComp();
927 
928 #ifdef _OPENMP
929 #pragma omp parallel if (Gpu::notInLaunchRegion())
930 #endif
931  for (MFIter mfi(eddyViscosity); mfi.isValid(); ++mfi)
932  {
933  Box vbx = mfi.validbox();
934 
935  Box planex_lo = mfi.growntilebox(1); planex_lo.setBig(0, vbx.smallEnd(0)-1);
936  Box planey_lo = mfi.growntilebox(1); planey_lo.setBig(1, vbx.smallEnd(1)-1);
937  Box planez_lo = mfi.growntilebox(1); planez_lo.setBig(2, vbx.smallEnd(2)-1);
938 
939  Box planex_hi = mfi.growntilebox(1); planex_hi.setSmall(0, vbx.bigEnd(0)+1);
940  Box planey_hi = mfi.growntilebox(1); planey_hi.setSmall(1, vbx.bigEnd(1)+1);
941  Box planez_hi = mfi.growntilebox(1); planez_hi.setSmall(2, vbx.bigEnd(2)+1);
942 
943  int i_lo = vbx.smallEnd(0); int i_hi = vbx.bigEnd(0);
944  int j_lo = vbx.smallEnd(1); int j_hi = vbx.bigEnd(1);
945  int k_lo = vbx.smallEnd(2); int k_hi = vbx.bigEnd(2);
946 
947  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
948  const Array4<int>& mask_arr = cc_mask.array(mfi);
949 
950  auto domlo = lbound(domain);
951  auto domhi = ubound(domain);
952 
953  ParallelFor(planex_lo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
954  {
955  int lj = amrex::min(amrex::max(j, domlo.y), domhi.y);
956  int lk = amrex::min(amrex::max(k, domlo.z), domhi.z);
957  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i_lo,lj,lk) != is_notcovered) ||
958  (mask_arr(i,j,k) == is_physbnd && i < domlo.x && impose_phys_bcs)) {
959  for (int n = 0; n < ncomp; n++) {
960  mu_turb(i,j,k,n) = mu_turb(i_lo,lj,lk,n);
961  }
962  }
963  });
964  ParallelFor(planex_hi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
965  {
966  int lj = amrex::min(amrex::max(j, domlo.y), domhi.y);
967  int lk = amrex::min(amrex::max(k, domlo.z), domhi.z);
968  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i_hi,lj,lk) != is_notcovered) ||
969  (mask_arr(i,j,k) == is_physbnd && i > domhi.x && impose_phys_bcs)) {
970  for (int n = 0; n < ncomp; n++) {
971  mu_turb(i,j,k,n) = mu_turb(i_hi,lj,lk,n);
972  }
973  }
974  });
975  ParallelFor(planey_lo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
976  {
977  int lk = amrex::min(amrex::max(k, domlo.z), domhi.z);
978  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i,j_lo,lk) != is_notcovered) ||
979  (mask_arr(i,j,k) == is_physbnd && j < domlo.y && impose_phys_bcs)) {
980  for (int n = 0; n < ncomp; n++) {
981  mu_turb(i,j,k,n) = mu_turb(i,j_lo,lk,n);
982  }
983  }
984  });
985  ParallelFor(planey_hi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
986  {
987  int lk = amrex::min(amrex::max(k, domlo.z), domhi.z);
988  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i,j_hi,lk) != is_notcovered)||
989  (mask_arr(i,j,k) == is_physbnd && j > domhi.y && impose_phys_bcs)) {
990  for (int n = 0; n < ncomp; n++) {
991  mu_turb(i,j,k,n) = mu_turb(i,j_hi,lk,n);
992  }
993  }
994  });
995  ParallelFor(planez_lo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
996  {
997  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i,j,k_lo) != is_notcovered) ||
998  (mask_arr(i,j,k) == is_physbnd && k < domlo.z && impose_phys_bcs)) {
999  if (mask_arr(i,j,k) == is_physbnd) {
1000  }
1001  for (int n = 0; n < ncomp; n++) {
1002  mu_turb(i,j,k,n) = mu_turb(i,j,k_lo,n);
1003  }
1004  }
1005  });
1006  ParallelFor(planez_hi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
1007  {
1008  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i,j,k_hi) != is_notcovered) ||
1009  (mask_arr(i,j,k) == is_physbnd && k > domhi.z && impose_phys_bcs)) {
1010  for (int n = 0; n < ncomp; n++) {
1011  mu_turb(i,j,k,n) = mu_turb(i,j,k_hi,n);
1012  }
1013  }
1014  });
1015  } // mfi
1016  eddyViscosity.FillBoundary(geom.periodicity());
1017  }
1018 }
void ComputeDiffusivityMRF(const MultiFab &xvel, const MultiFab &yvel, const MultiFab &cons_in, MultiFab &eddyViscosity, const Geometry &geom, const TurbChoice &turbChoice, std::unique_ptr< SurfaceLayer > &SurfLayer, bool use_terrain_fitted_coords, bool, int level, const BCRec *bc_ptr, bool, const std::unique_ptr< MultiFab > &z_phys_nd, const MoistureComponentIndices &moisture_indices)
Definition: ERF_ComputeDiffusivityMRF.cpp:12
void ComputeDiffusivityMYJ(Real dt, const MultiFab &xvel, const MultiFab &yvel, MultiFab &cons_in, MultiFab &eddyViscosity, const Geometry &geom, const TurbChoice &, std::unique_ptr< SurfaceLayer > &, bool use_terrain_fitted_coords, bool, int, const BCRec *bc_ptr, bool, const std::unique_ptr< MultiFab > &z_phys_nd, const MoistureComponentIndices &moisture_indices)
Definition: ERF_ComputeDiffusivityMYJ.cpp:13
void ComputeDiffusivityMYNN25(const MultiFab &xvel, const MultiFab &yvel, const MultiFab &cons_in, MultiFab &eddyViscosity, const Geometry &geom, const TurbChoice &turbChoice, std::unique_ptr< SurfaceLayer > &SurfLayer, bool use_terrain_fitted_coords, bool use_moisture, int level, const BCRec *bc_ptr, bool, const std::unique_ptr< MultiFab > &z_phys_nd, const MoistureComponentIndices &moisture_indices)
Definition: ERF_ComputeDiffusivityMYNN25.cpp:13
void ComputeDiffusivityMYNNEDMF(const MultiFab &xvel, const MultiFab &yvel, const MultiFab &cons_in, MultiFab &eddyViscosity, const Geometry &geom, const TurbChoice &turbChoice, std::unique_ptr< SurfaceLayer > &SurfLayer, bool use_terrain_fitted_coords, bool use_moisture, int level, const BCRec *bc_ptr, bool, const std::unique_ptr< MultiFab > &z_phys_nd, const MoistureComponentIndices &moisture_indices)
Definition: ERF_ComputeDiffusivityMYNNEDMF.cpp:4177
void ComputeDiffusivityYSU(const MultiFab &xvel, const MultiFab &yvel, const MultiFab &cons_in, MultiFab &eddyViscosity, const Geometry &geom, const TurbChoice &turbChoice, std::unique_ptr< SurfaceLayer > &SurfLayer, bool use_terrain_fitted_coords, bool, int level, const BCRec *bc_ptr, bool, const std::unique_ptr< MultiFab > &z_phys_nd, const MoistureComponentIndices &moisture_indices)
Definition: ERF_ComputeDiffusivityYSU.cpp:11
void ComputeTurbulentViscosityRANS(Vector< std::unique_ptr< MultiFab >> &, const MultiFab &cons_in, const MultiFab &wdist, MultiFab &eddyViscosity, MultiFab &Hfx1, MultiFab &Hfx2, MultiFab &Hfx3, MultiFab &Diss, const Geometry &geom, bool use_terrain_fitted_coords, Vector< std::unique_ptr< MultiFab >> &, const std::unique_ptr< MultiFab > &z_phys_nd, const TurbChoice &turbChoice, const Real const_grav, std::unique_ptr< SurfaceLayer > &SurfLayer, const MultiFab *z_0)
Definition: ERF_ComputeTurbulentViscosity.cpp:563
void ComputeTurbulentViscosity(Real dt, const MultiFab &xvel, const MultiFab &yvel, Vector< std::unique_ptr< MultiFab >> &Tau_lev, MultiFab &cons_in, const MultiFab &wdist, MultiFab &eddyViscosity, MultiFab &Hfx1, MultiFab &Hfx2, MultiFab &Hfx3, MultiFab &Diss, const Geometry &geom, Vector< std::unique_ptr< MultiFab >> &mapfac, const std::unique_ptr< MultiFab > &z_phys_nd, const SolverChoice &solverChoice, std::unique_ptr< SurfaceLayer > &SurfLayer, const MultiFab *z_0, const bool &use_terrain_fitted_coords, const bool &use_moisture, int level, const BCRec *bc_ptr, const eb_ &ebfact, bool vert_only)
Definition: ERF_ComputeTurbulentViscosity.cpp:788
void ComputeTurbulentViscosityLES(Vector< std::unique_ptr< MultiFab >> &Tau_lev, const MultiFab &cons_in, MultiFab &eddyViscosity, MultiFab &Hfx1, MultiFab &Hfx2, MultiFab &Hfx3, MultiFab &Diss, const Geometry &geom, bool use_terrain_fitted_coords, Vector< std::unique_ptr< MultiFab >> &mapfac, const std::unique_ptr< MultiFab > &z_phys_nd, const TurbChoice &turbChoice, const Real const_grav, std::unique_ptr< SurfaceLayer > &, const MoistureComponentIndices &moisture_indices, const MultiFab *xvel, const MultiFab *yvel)
Definition: ERF_ComputeTurbulentViscosity.cpp:30
void ComputeTurbulentViscosityLES_EB(Vector< std::unique_ptr< MultiFab >> &Tau_lev, const MultiFab &cons_in, MultiFab &eddyViscosity, MultiFab &Hfx1, MultiFab &Hfx2, MultiFab &Hfx3, const Geometry &geom, Vector< std::unique_ptr< MultiFab >> &mapfac, const TurbChoice &turbChoice, const Real const_grav, [[maybe_unused]] const SolverChoice &solverChoice, std::unique_ptr< SurfaceLayer > &, const MoistureComponentIndices &moisture_indices, const eb_ &ebfact, const MultiFab *xvel, const MultiFab *yvel)
Definition: ERF_ComputeTurbulentViscosity.cpp:336
const bool use_moisture
Definition: ERF_InitCustomPert_Bomex.H:14
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
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+0.5) *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<=1.0) { Real dT=T_pert *(std::cos(PI *L)+1.0)/2.0;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);} } })
for(int i=0;i< m_num_species;i++)
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:48
TurbChoice turbChoice
Definition: ERF_SetupVertDiff.H:5
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Real z_0
Definition: ERF_UpdateWSubsidence_Bomex.H:10
@ xvel
Definition: ERF_IndexDefines.H:157
@ yvel
Definition: ERF_IndexDefines.H:158
amrex::Real gravity
Definition: ERF_DataStruct.H:1142
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:1080
static TerrainType terrain_type
Definition: ERF_DataStruct.H:1056
MoistureComponentIndices moisture_indices
Definition: ERF_DataStruct.H:1211
Definition: ERF_TurbStruct.H:42
PBLType pbl_type
Definition: ERF_TurbStruct.H:418
RANSType rans_type
Definition: ERF_TurbStruct.H:413
LESType les_type
Definition: ERF_TurbStruct.H:371

Referenced by ERF::advance_dycore().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeTurbulentViscosityLES()

void ComputeTurbulentViscosityLES ( Vector< std::unique_ptr< MultiFab >> &  Tau_lev,
const MultiFab &  cons_in,
MultiFab &  eddyViscosity,
MultiFab &  Hfx1,
MultiFab &  Hfx2,
MultiFab &  Hfx3,
MultiFab &  Diss,
const Geometry &  geom,
bool  use_terrain_fitted_coords,
Vector< std::unique_ptr< MultiFab >> &  mapfac,
const std::unique_ptr< MultiFab > &  z_phys_nd,
const TurbChoice turbChoice,
const Real  const_grav,
std::unique_ptr< SurfaceLayer > &  ,
const MoistureComponentIndices moisture_indices,
const MultiFab *  xvel,
const MultiFab *  yvel 
)

Function for computing the turbulent viscosity with LES.

Parameters
[in]Tau_levstrain at this level
[in]cons_incell center conserved quantities
[out]eddyViscosityturbulent viscosity
[in]Hfx1heat flux in x-dir
[in]Hfx2heat flux in y-dir
[in]Hfx3heat flux in z-dir
[in]Dissdissipation of turbulent kinetic energy
[in]geomproblem geometry
[in]mapfacmap factors
[in]turbChoicecontainer with turbulence parameters
[in]xvelx-direction velocity (for moist Ri correction)
[in]yvely-direction velocity (for moist Ri correction)
41 {
42  const GpuArray<Real, AMREX_SPACEDIM> cellSizeInv = geom.InvCellSizeArray();
43  const Box& domain = geom.Domain();
44 
45  Real inv_Pr_t = turbChoice.Pr_t_inv;
46  Real inv_Sc_t = turbChoice.Sc_t_inv;
47  Real inv_sigma_k = 1.0 / turbChoice.sigma_k;
48 
49  bool use_thetav_grad = (turbChoice.strat_type == StratType::thetav);
50  bool use_thetal_grad = (turbChoice.strat_type == StratType::thetal);
51 
52  bool isotropic = turbChoice.mix_isotropic;
53 
54  // SMAGORINSKY: Fill Kturb for momentum in horizontal and vertical
55  //***********************************************************************************
56  if (turbChoice.les_type == LESType::Smagorinsky)
57  {
58  Real Cs = turbChoice.Cs;
59  bool smag2d = turbChoice.smag2d;
60 
61  // Define variables required inside device lambdas (scalars only)
62  Real l_abs_g = const_grav;
63  Real l_Ri_crit = turbChoice.Ri_crit;
64  bool l_use_Ri_corr = turbChoice.use_Ri_correction;
65  bool l_has_xvel = (xvel != nullptr);
66  bool l_has_yvel = (yvel != nullptr);
67 
68 #ifdef _OPENMP
69 #pragma omp parallel if (Gpu::notInLaunchRegion())
70 #endif
71  for (MFIter mfi(eddyViscosity,TilingIfNotGPU()); mfi.isValid(); ++mfi)
72  {
73  Box bxcc = mfi.growntilebox(1) & domain;
74  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
75  const Array4<Real>& hfx_x = Hfx1.array(mfi);
76  const Array4<Real>& hfx_y = Hfx2.array(mfi);
77  const Array4<Real>& hfx_z = Hfx3.array(mfi);
78  const Array4<Real const >& cell_data = cons_in.array(mfi);
79  Array4<Real const> tau11 = Tau_lev[TauType::tau11]->array(mfi);
80  Array4<Real const> tau22 = Tau_lev[TauType::tau22]->array(mfi);
81  Array4<Real const> tau33 = Tau_lev[TauType::tau33]->array(mfi);
82  Array4<Real const> tau12 = Tau_lev[TauType::tau12]->array(mfi);
83  Array4<Real const> tau13 = Tau_lev[TauType::tau13]->array(mfi);
84  Array4<Real const> tau23 = Tau_lev[TauType::tau23]->array(mfi);
85  Array4<Real const> mf_u = mapfac[MapFacType::u_x]->const_array(mfi);
86  Array4<Real const> mf_v = mapfac[MapFacType::v_y]->const_array(mfi);
87  Array4<Real const> z_nd_arr = z_phys_nd->const_array(mfi);
88 
89  Array4<Real const> u_arr = (l_has_xvel) ? xvel->const_array(mfi) : Array4<Real const>{};
90  Array4<Real const> v_arr = (l_has_yvel) ? yvel->const_array(mfi) : Array4<Real const>{};
91 
92  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
93  {
94  // =====================================================================
95  // 1. STRAIN RATE MAGNITUDE CALCULATION
96  // =====================================================================
97  Real SmnSmn;
98  if (smag2d) {
99  SmnSmn = ComputeSmnSmn2D(i,j,k,tau11,tau22,tau12);
100  } else {
101  SmnSmn = ComputeSmnSmn(i,j,k,tau11,tau22,tau33,tau12,tau13,tau23);
102  }
103  Real strain_rate_magnitude = std::sqrt(2.0 * SmnSmn);
104 
105  // =====================================================================
106  // 2. GRID SCALE CALCULATION (filter width Δ)
107  // =====================================================================
108  Real dxInv = cellSizeInv[0];
109  Real dyInv = cellSizeInv[1];
110  Real dzInv = cellSizeInv[2];
111  if (use_terrain_fitted_coords) {
112  dzInv /= Compute_h_zeta_AtCellCenter(i,j,k, cellSizeInv, z_nd_arr);
113  }
114 
115  Real Delta;
116  Real DeltaH;
117  if (isotropic) {
118  Real cellVolMsf = 1.0 / (dxInv * mf_u(i,j,0) * dyInv * mf_v(i,j,0) * dzInv);
119  Delta = std::cbrt(cellVolMsf);
120  DeltaH = Delta;
121  } else {
122  Delta = 1.0 / dzInv;
123  DeltaH = std::sqrt(1.0 / (dxInv * mf_u(i,j,0) * dyInv * mf_v(i,j,0)));
124  }
125 
126  Real rho = cell_data(i, j, k, Rho_comp);
127  Real CsDeltaSqr_h = Cs * Cs * DeltaH * DeltaH;
128  Real CsDeltaSqr_v = Cs * Cs * Delta * Delta;
129 
130  Real nu_turb_base_h = CsDeltaSqr_h * strain_rate_magnitude;
131  Real nu_turb_base_v = CsDeltaSqr_v * strain_rate_magnitude;
132 
133  Real stability_factor = 1.0;
134 
135  if (l_use_Ri_corr && l_has_xvel && l_has_yvel) {
136  Real N2 = ComputeN2(i, j, k, dzInv, l_abs_g, cell_data, moisture_indices);
137  Real S2_vert = ComputeVerticalShear2(i, j, k, dzInv, u_arr, v_arr);
138  Real Ri = ComputeRichardson(N2, S2_vert);
139  stability_factor = StabilityFunction(Ri, l_Ri_crit);
140  }
141 
142  if (isotropic) {
143  mu_turb(i, j, k, EddyDiff::Mom_h) = rho * nu_turb_base_h * stability_factor;
144  mu_turb(i, j, k, EddyDiff::Mom_v) = rho * nu_turb_base_v * stability_factor;
145  } else {
146  mu_turb(i, j, k, EddyDiff::Mom_h) = rho * nu_turb_base_h;
147  mu_turb(i, j, k, EddyDiff::Mom_v) = rho * nu_turb_base_v * stability_factor;
148  }
149 
150  Real dtheta_dz = 0.5 * ( cell_data(i,j,k+1,RhoTheta_comp)/cell_data(i,j,k+1,Rho_comp)
151  - cell_data(i,j,k-1,RhoTheta_comp)/cell_data(i,j,k-1,Rho_comp) )*dzInv;
152 
153  hfx_x(i,j,k) = 0.0;
154  hfx_y(i,j,k) = 0.0;
155  hfx_z(i,j,k) = -inv_Pr_t * mu_turb(i,j,k,EddyDiff::Mom_v) * dtheta_dz;
156  });
157  }
158  }
159  // DEARDORFF: Fill Kturb for momentum in horizontal and vertical
160  //***********************************************************************************
161  else if (turbChoice.les_type == LESType::Deardorff)
162  {
163  const Real l_C_k = turbChoice.Ck;
164  const Real l_C_e = turbChoice.Ce;
165  const Real l_C_e_wall = turbChoice.Ce_wall;
166  const Real Ce_lcoeff = amrex::max(0.0, l_C_e - 1.9*l_C_k);
167  const Real l_abs_g = const_grav;
168 
169  const bool use_ref_theta = (turbChoice.theta_ref > 0);
170  const Real l_inv_theta0 = (use_ref_theta) ? 1.0 / turbChoice.theta_ref : 1.0;
171 
172 #ifdef _OPENMP
173 #pragma omp parallel if (Gpu::notInLaunchRegion())
174 #endif
175  for ( MFIter mfi(eddyViscosity,TilingIfNotGPU()); mfi.isValid(); ++mfi)
176  {
177  Box bxcc = mfi.tilebox();
178 
179  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
180  const Array4<Real>& hfx_x = Hfx1.array(mfi);
181  const Array4<Real>& hfx_y = Hfx2.array(mfi);
182  const Array4<Real>& hfx_z = Hfx3.array(mfi);
183  const Array4<Real>& diss = Diss.array(mfi);
184 
185  const Array4<Real const > &cell_data = cons_in.array(mfi);
186 
187  Array4<Real const> mf_u = mapfac[MapFacType::u_x]->const_array(mfi);
188  Array4<Real const> mf_v = mapfac[MapFacType::v_y]->const_array(mfi);
189 
190  Array4<Real const> z_nd_arr = z_phys_nd->const_array(mfi);
191 
192  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
193  {
194  Real dxInv = cellSizeInv[0];
195  Real dyInv = cellSizeInv[1];
196  Real dzInv = cellSizeInv[2];
197  if (use_terrain_fitted_coords) {
198  // the terrain grid is only deformed in z for now
199  dzInv /= Compute_h_zeta_AtCellCenter(i,j,k, cellSizeInv, z_nd_arr);
200  }
201  Real Delta;
202  if (isotropic) {
203  Real cellVolMsf = 1.0 / (dxInv * mf_u(i,j,0) * dyInv * mf_v(i,j,0) * dzInv);
204  Delta = std::cbrt(cellVolMsf);
205  } else {
206  Delta = 1.0 / dzInv;
207  }
208 
209  Real dtheta_dz;
210  if (use_thetav_grad) {
211  dtheta_dz = 0.5 * ( GetThetav(i, j, k+1, cell_data, moisture_indices)
212  - GetThetav(i, j, k-1, cell_data, moisture_indices) )*dzInv;
213  } else if (use_thetal_grad) {
214  dtheta_dz = 0.5 * ( GetThetal(i, j, k+1, cell_data, moisture_indices)
215  - GetThetal(i, j, k-1, cell_data, moisture_indices) )*dzInv;
216  } else {
217  dtheta_dz = 0.5 * ( cell_data(i, j, k+1, RhoTheta_comp) / cell_data(i, j, k+1, Rho_comp)
218  - cell_data(i, j, k-1, RhoTheta_comp) / cell_data(i, j, k-1, Rho_comp) )*dzInv;
219  }
220 
221  // Calculate stratification-dependent mixing length (Deardorff 1980, Eqn. 10a)
222  Real E = amrex::max(cell_data(i,j,k,RhoKE_comp)/cell_data(i,j,k,Rho_comp),Real(0.0));
223  Real stratification = l_abs_g * dtheta_dz * l_inv_theta0;
224  if (!use_ref_theta) {
225  // l_inv_theta0 == 1, divide by actual theta
226  stratification *= cell_data(i,j,k,Rho_comp) /
227  cell_data(i,j,k,RhoTheta_comp);
228  }
229 
230  // Following WRF, the stratification effects are applied to the vertical length scales
231  // in the case of anisotropic mixing
232  Real length;
234  if (stratification <= eps) {
235  length = Delta; // cbrt(dx*dy*dz) -or- dz
236  } else {
237  length = 0.76 * std::sqrt(E / amrex::max(stratification,eps));
238  // mixing length should be _reduced_ for stable stratification
239  length = amrex::min(length, Delta);
240  // following WRF, make sure the mixing length isn't too small
241  length = amrex::max(length, 0.001 * Delta);
242  }
243 
244  Real DeltaH = (isotropic) ? length : std::sqrt(1.0 / (dxInv * mf_u(i,j,0) * dyInv * mf_v(i,j,0)));
245 
246  Real Pr_inv_v = (1. + 2.*length/Delta);
247  Real Pr_inv_h = (isotropic) ? Pr_inv_v : inv_Pr_t;
248 
249  // Calculate eddy diffusivities
250  // K = rho * C_k * l * KE^(1/2)
251  mu_turb(i,j,k,EddyDiff::Mom_h) = cell_data(i,j,k,Rho_comp) * l_C_k * DeltaH * std::sqrt(E);
252  mu_turb(i,j,k,EddyDiff::Mom_v) = cell_data(i,j,k,Rho_comp) * l_C_k * length * std::sqrt(E);
253  // KH = (1 + 2*l/delta) * mu_turb
254  mu_turb(i,j,k,EddyDiff::Theta_h) = Pr_inv_h * mu_turb(i,j,k,EddyDiff::Mom_h);
255  mu_turb(i,j,k,EddyDiff::Theta_v) = Pr_inv_v * mu_turb(i,j,k,EddyDiff::Mom_v);
256  // Store lengthscale for TKE source terms
257  mu_turb(i,j,k,EddyDiff::Turb_lengthscale) = length;
258 
259  // Calculate SFS quantities
260  // - dissipation
261  Real Ce;
262  if ((l_C_e_wall > 0) && (k==0)) {
263  Ce = l_C_e_wall;
264  } else {
265  Ce = 1.9*l_C_k + Ce_lcoeff*length / Delta;
266  }
267  diss(i,j,k) = cell_data(i,j,k,Rho_comp) * Ce * std::pow(E,1.5) / length;
268 
269  // - heat flux
270  // (Note: If using SurfaceLayer, the value at k=0 will
271  // be overwritten)
272  hfx_x(i,j,k) = 0.0;
273  hfx_y(i,j,k) = 0.0;
274  hfx_z(i,j,k) = -mu_turb(i,j,k,EddyDiff::Theta_v) * dtheta_dz; // (rho*w)' theta' [kg m^-2 s^-1 K]
275  });
276  }
277  }
278 
279  // Extrapolate Kturb in x/y, fill remaining elements (relevant to lev==0)
280  //***********************************************************************************
281  int ngc(1);
282  // EddyDiff mapping : Theta_h KE_h Scalar_h Q_h
283  Vector<Real> Factors = {inv_Pr_t, inv_sigma_k, inv_Sc_t, inv_Sc_t}; // alpha = mu/Pr
284  Gpu::AsyncVector<Real> d_Factors; d_Factors.resize(Factors.size());
285  Gpu::copy(Gpu::hostToDevice, Factors.begin(), Factors.end(), d_Factors.begin());
286  Real* fac_ptr = d_Factors.data();
287 
288  const bool use_KE = ( turbChoice.les_type == LESType::Deardorff );
289 
290 #ifdef _OPENMP
291 #pragma omp parallel if (Gpu::notInLaunchRegion())
292 #endif
293  for ( MFIter mfi(eddyViscosity,TilingIfNotGPU()); mfi.isValid(); ++mfi)
294  {
295  Box bxcc = mfi.tilebox();
296  Box planex = bxcc; planex.setSmall(0, 1); planex.setBig(0, ngc); planex.grow(1,1);
297  Box planey = bxcc; planey.setSmall(1, 1); planey.setBig(1, ngc); planey.grow(0,1);
298  bxcc.growLo(0,ngc); bxcc.growHi(0,ngc);
299  bxcc.growLo(1,ngc); bxcc.growHi(1,ngc);
300 
301  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
302 
303  for (auto n = 1; n < (EddyDiff::NumDiffs-1)/2; ++n) {
304  int offset = (EddyDiff::NumDiffs-1)/2;
305  switch (n)
306  {
307  case EddyDiff::KE_h:
308  if (use_KE) {
309  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
310  {
311  int indx = n;
312  int indx_v = indx + offset;
313  mu_turb(i,j,k,indx) = mu_turb(i,j,k,EddyDiff::Mom_h) * fac_ptr[indx-1];
314  mu_turb(i,j,k,indx_v) = mu_turb(i,j,k,indx);
315  });
316  }
317  break;
318  default:
319  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
320  {
321  int indx = n;
322  int indx_v = indx + offset;
323 
324  // NOTE: Theta_h, Theta_v have already been set for Deardorff
325  if (!(indx_v == EddyDiff::Theta_v && use_KE)) {
326  mu_turb(i,j,k,indx) = mu_turb(i,j,k,EddyDiff::Mom_h) * fac_ptr[indx-1];
327  mu_turb(i,j,k,indx_v) = mu_turb(i,j,k,indx);
328  }
329  });
330  break;
331  }
332  }
333  }
334 }
if(l_use_mynn &&start_comp<=RhoKE_comp &&end_comp >=RhoKE_comp)
Definition: ERF_AddQKESources.H:2
@ tau12
Definition: ERF_DataStruct.H:32
@ tau23
Definition: ERF_DataStruct.H:32
@ tau33
Definition: ERF_DataStruct.H:32
@ tau22
Definition: ERF_DataStruct.H:32
@ tau11
Definition: ERF_DataStruct.H:32
@ tau13
Definition: ERF_DataStruct.H:32
@ v_y
Definition: ERF_DataStruct.H:25
@ u_x
Definition: ERF_DataStruct.H:24
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeSmnSmn2D(int &i, int &j, int &k, const amrex::Array4< amrex::Real const > &tau11, const amrex::Array4< amrex::Real const > &tau22, const amrex::Array4< amrex::Real const > &tau12)
Definition: ERF_EddyViscosity.H:204
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeSmnSmn(int &i, int &j, int &k, const amrex::Array4< amrex::Real const > &tau11, const amrex::Array4< amrex::Real const > &tau22, const amrex::Array4< amrex::Real const > &tau33, const amrex::Array4< amrex::Real const > &tau12, const amrex::Array4< amrex::Real const > &tau13, const amrex::Array4< amrex::Real const > &tau23)
Definition: ERF_EddyViscosity.H:83
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:16
rho
Definition: ERF_InitCustomPert_Bubble.H:106
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real GetThetal(int i, int j, int k, const amrex::Array4< amrex::Real const > &cell_data, const MoistureComponentIndices &moisture_indices)
Definition: ERF_MoistUtils.H:101
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real GetThetav(int i, int j, int k, const amrex::Array4< amrex::Real const > &cell_data, const MoistureComponentIndices &moisture_indices)
Definition: ERF_MoistUtils.H:72
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real StabilityFunction(amrex::Real Ri, amrex::Real Ri_crit)
Definition: ERF_RichardsonNumber.H:386
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeVerticalShear2(int i, int j, int k, amrex::Real dzInv, const amrex::Array4< const amrex::Real > &u, const amrex::Array4< const amrex::Real > &v)
Definition: ERF_RichardsonNumber.H:281
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeRichardson(amrex::Real N2_moist, amrex::Real S2_vert)
Definition: ERF_RichardsonNumber.H:370
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeN2(int i, int j, int k, amrex::Real dzInv, amrex::Real const_grav, const amrex::Array4< const amrex::Real > &cell_data, const MoistureComponentIndices &moisture_indices)
Definition: ERF_RichardsonNumber.H:126
const bool use_ref_theta
Definition: ERF_SetupDiff.H:16
const Real l_inv_theta0
Definition: ERF_SetupDiff.H:17
AMREX_FORCE_INLINE AMREX_GPU_DEVICE amrex::Real Compute_h_zeta_AtCellCenter(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:53
@ Theta_v
Definition: ERF_IndexDefines.H:192
@ Turb_lengthscale
Definition: ERF_IndexDefines.H:196
@ Mom_h
Definition: ERF_IndexDefines.H:186
@ Mom_v
Definition: ERF_IndexDefines.H:191
@ NumDiffs
Definition: ERF_IndexDefines.H:197
@ Theta_h
Definition: ERF_IndexDefines.H:187
@ KE_h
Definition: ERF_IndexDefines.H:188
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12
bool smag2d
Definition: ERF_TurbStruct.H:383
amrex::Real sigma_k
Definition: ERF_TurbStruct.H:399
amrex::Real Sc_t_inv
Definition: ERF_TurbStruct.H:379
amrex::Real Ri_crit
Definition: ERF_TurbStruct.H:410
StratType strat_type
Definition: ERF_TurbStruct.H:404
amrex::Real Ck
Definition: ERF_TurbStruct.H:388
bool use_Ri_correction
Definition: ERF_TurbStruct.H:409
amrex::Real Cs
Definition: ERF_TurbStruct.H:382
amrex::Real Pr_t_inv
Definition: ERF_TurbStruct.H:375
amrex::Real Ce_wall
Definition: ERF_TurbStruct.H:387
amrex::Real Ce
Definition: ERF_TurbStruct.H:386
amrex::Real theta_ref
Definition: ERF_TurbStruct.H:402
bool mix_isotropic
Definition: ERF_TurbStruct.H:407

Referenced by ComputeTurbulentViscosity().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeTurbulentViscosityLES_EB()

void ComputeTurbulentViscosityLES_EB ( Vector< std::unique_ptr< MultiFab >> &  Tau_lev,
const MultiFab &  cons_in,
MultiFab &  eddyViscosity,
MultiFab &  Hfx1,
MultiFab &  Hfx2,
MultiFab &  Hfx3,
const Geometry &  geom,
Vector< std::unique_ptr< MultiFab >> &  mapfac,
const TurbChoice turbChoice,
const Real  const_grav,
[[maybe_unused] ] const SolverChoice solverChoice,
std::unique_ptr< SurfaceLayer > &  ,
const MoistureComponentIndices moisture_indices,
const eb_ ebfact,
const MultiFab *  xvel,
const MultiFab *  yvel 
)
348 {
349  const GpuArray<Real, AMREX_SPACEDIM> cellSizeInv = geom.InvCellSizeArray();
350  Real dxInv = cellSizeInv[0];
351  Real dyInv = cellSizeInv[1];
352  Real dzInv = cellSizeInv[2];
353 
354  const Box& domain = geom.Domain();
355 
356  Real inv_Pr_t = turbChoice.Pr_t_inv;
357  Real inv_Sc_t = turbChoice.Sc_t_inv;
358  Real inv_sigma_k = 1.0 / turbChoice.sigma_k;
359 
360  bool isotropic = turbChoice.mix_isotropic;
361 
362  AMREX_ASSERT(turbChoice.les_type == LESType::Smagorinsky);
363 
364  // SMAGORINSKY: Fill Kturb for momentum in horizontal and vertical
365  //***********************************************************************************
366  if (turbChoice.les_type == LESType::Smagorinsky)
367  {
368  Real Cs = turbChoice.Cs;
369 
370  // Define variables required inside device lambdas (scalars only)
371  Real l_abs_g = const_grav;
372  Real l_Ri_crit = turbChoice.Ri_crit;
373  bool l_use_Ri_corr = turbChoice.use_Ri_correction;
374  bool l_has_xvel = (xvel != nullptr);
375  bool l_has_yvel = (yvel != nullptr);
376 
377 #ifdef _OPENMP
378 #pragma omp parallel if (Gpu::notInLaunchRegion())
379 #endif
380  for (MFIter mfi(eddyViscosity,TilingIfNotGPU()); mfi.isValid(); ++mfi)
381  {
382  Box bxcc = mfi.growntilebox(1) & domain;
383  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
384  const Array4<Real>& hfx_x = Hfx1.array(mfi);
385  const Array4<Real>& hfx_y = Hfx2.array(mfi);
386  const Array4<Real>& hfx_z = Hfx3.array(mfi);
387  const Array4<Real const >& cell_data = cons_in.array(mfi);
388  Array4<Real const> tau11 = Tau_lev[TauType::tau11]->array(mfi);
389  Array4<Real const> tau22 = Tau_lev[TauType::tau22]->array(mfi);
390  Array4<Real const> tau33 = Tau_lev[TauType::tau33]->array(mfi);
391  Array4<Real const> tau12 = Tau_lev[TauType::tau12]->array(mfi);
392  Array4<Real const> tau13 = Tau_lev[TauType::tau13]->array(mfi);
393  Array4<Real const> tau23 = Tau_lev[TauType::tau23]->array(mfi);
394  Array4<Real const> mf_u = mapfac[MapFacType::u_x]->const_array(mfi);
395  Array4<Real const> mf_v = mapfac[MapFacType::v_y]->const_array(mfi);
396 
397  Array4<Real const> u_arr = (l_has_xvel) ? xvel->const_array(mfi) : Array4<Real const>{};
398  Array4<Real const> v_arr = (l_has_yvel) ? yvel->const_array(mfi) : Array4<Real const>{};
399 
400  Array4<const EBCellFlag> c_cflag = (ebfact.get_const_factory())->getMultiEBCellFlagFab()[mfi].const_array();
401  Array4<const EBCellFlag> u_cflag = (ebfact.get_u_const_factory())->getMultiEBCellFlagFab()[mfi].const_array();
402  Array4<const EBCellFlag> v_cflag = (ebfact.get_v_const_factory())->getMultiEBCellFlagFab()[mfi].const_array();
403  Array4<const EBCellFlag> w_cflag = (ebfact.get_w_const_factory())->getMultiEBCellFlagFab()[mfi].const_array();
404  Array4<const Real > u_vfrac = (ebfact.get_u_const_factory())->getVolFrac().const_array(mfi);
405  Array4<const Real > v_vfrac = (ebfact.get_v_const_factory())->getVolFrac().const_array(mfi);
406 
407  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
408  {
409  // =====================================================================
410  // 1. STRAIN RATE MAGNITUDE CALCULATION
411  // =====================================================================
412  Real SmnSmn=0.0;
413  if (c_cflag(i,j,k).isRegular()) {
414  SmnSmn = ComputeSmnSmn(i,j,k,tau11,tau22,tau33,tau12,tau13,tau23);
415  } else if (c_cflag(i,j,k).isSingleValued()) {
416  SmnSmn = ComputeSmnSmn_EB(i,j,k,tau11,tau22,tau33,tau12,tau13,tau23,c_cflag,u_cflag,v_cflag,w_cflag);
417  }
418  Real strain_rate_magnitude = std::sqrt(2.0 * SmnSmn);
419 
420  // =====================================================================
421  // 2. GRID SCALE CALCULATION (filter width Δ)
422  // =====================================================================
423  Real Delta;
424  Real DeltaH;
425  if (isotropic) {
426  Real cellVolMsf = 1.0 / (dxInv * mf_u(i,j,0) * dyInv * mf_v(i,j,0) * dzInv);
427  Delta = std::cbrt(cellVolMsf);
428  DeltaH = Delta;
429  } else {
430  Delta = 1.0 / dzInv;
431  DeltaH = std::sqrt(1.0 / (dxInv * mf_u(i,j,0) * dyInv * mf_v(i,j,0)));
432  }
433 
434  Real rho = cell_data(i, j, k, Rho_comp);
435  Real CsDeltaSqr_h = Cs * Cs * DeltaH * DeltaH;
436  Real CsDeltaSqr_v = Cs * Cs * Delta * Delta;
437 
438  Real nu_turb_base_h = CsDeltaSqr_h * strain_rate_magnitude;
439  Real nu_turb_base_v = CsDeltaSqr_v * strain_rate_magnitude;
440 
441  Real stability_factor = 1.0;
442 
443  if (l_use_Ri_corr && l_has_xvel && l_has_yvel) {
444  Real N2 = 0.0;
445  Real S2_vert = 0.0;
446  if (c_cflag(i,j,k).isRegular()) {
447  N2 = ComputeN2(i, j, k, dzInv, l_abs_g, cell_data, moisture_indices);
448  S2_vert = ComputeVerticalShear2(i, j, k, dzInv, u_arr, v_arr);
449  } else if (c_cflag(i,j,k).isSingleValued()) {
450  N2 = ComputeN2_EB(i, j, k, c_cflag, dzInv, l_abs_g, cell_data, moisture_indices);
451  S2_vert = ComputeVerticalShear2_EB(i, j, k, c_cflag, u_vfrac, v_vfrac, dzInv, u_arr, v_arr);
452  }
453  Real Ri = ComputeRichardson(N2, S2_vert);
454  stability_factor = StabilityFunction(Ri, l_Ri_crit);
455  }
456 
457  if (isotropic) {
458  mu_turb(i, j, k, EddyDiff::Mom_h) = rho * nu_turb_base_h * stability_factor;
459  mu_turb(i, j, k, EddyDiff::Mom_v) = rho * nu_turb_base_v * stability_factor;
460  } else {
461  mu_turb(i, j, k, EddyDiff::Mom_h) = rho * nu_turb_base_h;
462  mu_turb(i, j, k, EddyDiff::Mom_v) = rho * nu_turb_base_v * stability_factor;
463  }
464 
465  amrex::Real theta_km1 = cell_data(i,j,k-1,RhoTheta_comp)/cell_data(i,j,k-1,Rho_comp);
466  amrex::Real theta_kp1 = cell_data(i,j,k+1,RhoTheta_comp)/cell_data(i,j,k+1,Rho_comp);
467  Real dtheta_dz = 0.0;
468 
469  if (c_cflag(i,j,k).isRegular()) {
470  dtheta_dz = 0.5 * ( theta_kp1 - theta_km1 )*dzInv;
471  } else if (c_cflag(i,j,k).isSingleValued()) {
472 
473  amrex::Real theta = cell_data(i, j, k, RhoTheta_comp) / rho;
474  if (c_cflag(i,j,k+1).isCovered()) {
475  amrex::Real theta_km2 = cell_data(i, j, k-2, RhoTheta_comp) / cell_data(i, j, k-2, Rho_comp);
476  dtheta_dz = (3.*theta - 4.*theta_km1 + theta_km2) * 0.5 * dzInv;
477  } else if (c_cflag(i,j,k-1).isCovered()) {
478  amrex::Real theta_kp2 = cell_data(i, j, k+2, RhoTheta_comp) / cell_data(i, j, k+2, Rho_comp);
479  dtheta_dz = (-theta_kp2 + 4.*theta_kp1 - 3.*theta) * 0.5 * dzInv;
480  } else {
481  dtheta_dz = 0.5 * (theta_kp1 - theta_km1) * dzInv;
482  }
483  }
484 
485  hfx_x(i,j,k) = 0.0;
486  hfx_y(i,j,k) = 0.0;
487  hfx_z(i,j,k) = -inv_Pr_t * mu_turb(i,j,k,EddyDiff::Mom_v) * dtheta_dz;
488  });
489  }
490  }
491 
492  // Extrapolate Kturb in x/y, fill remaining elements (relevant to lev==0)
493  //***********************************************************************************
494  int ngc(1);
495  // EddyDiff mapping : Theta_h KE_h Scalar_h Q_h
496  Vector<Real> Factors = {inv_Pr_t, inv_sigma_k, inv_Sc_t, inv_Sc_t}; // alpha = mu/Pr
497  Gpu::AsyncVector<Real> d_Factors; d_Factors.resize(Factors.size());
498  Gpu::copy(Gpu::hostToDevice, Factors.begin(), Factors.end(), d_Factors.begin());
499  Real* fac_ptr = d_Factors.data();
500 
501  const bool use_KE = ( turbChoice.les_type == LESType::Deardorff );
502 
503 #ifdef _OPENMP
504 #pragma omp parallel if (Gpu::notInLaunchRegion())
505 #endif
506  for ( MFIter mfi(eddyViscosity,TilingIfNotGPU()); mfi.isValid(); ++mfi)
507  {
508  Box bxcc = mfi.tilebox();
509  Box planex = bxcc; planex.setSmall(0, 1); planex.setBig(0, ngc); planex.grow(1,1);
510  Box planey = bxcc; planey.setSmall(1, 1); planey.setBig(1, ngc); planey.grow(0,1);
511  bxcc.growLo(0,ngc); bxcc.growHi(0,ngc);
512  bxcc.growLo(1,ngc); bxcc.growHi(1,ngc);
513 
514  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
515 
516  for (auto n = 1; n < (EddyDiff::NumDiffs-1)/2; ++n) {
517  int offset = (EddyDiff::NumDiffs-1)/2;
518  switch (n)
519  {
520  case EddyDiff::KE_h:
521  if (use_KE) {
522  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
523  {
524  int indx = n;
525  int indx_v = indx + offset;
526  mu_turb(i,j,k,indx) = mu_turb(i,j,k,EddyDiff::Mom_h) * fac_ptr[indx-1];
527  mu_turb(i,j,k,indx_v) = mu_turb(i,j,k,indx);
528  });
529  }
530  break;
531  default:
532  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
533  {
534  int indx = n;
535  int indx_v = indx + offset;
536 
537  // NOTE: Theta_h, Theta_v have already been set for Deardorff
538  if (!(indx_v == EddyDiff::Theta_v && use_KE)) {
539  mu_turb(i,j,k,indx) = mu_turb(i,j,k,EddyDiff::Mom_h) * fac_ptr[indx-1];
540  mu_turb(i,j,k,indx_v) = mu_turb(i,j,k,indx);
541  }
542  });
543  break;
544  }
545  }
546  }
547 }
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeSmnSmn_EB(int &i, int &j, int &k, const amrex::Array4< amrex::Real const > &tau11, const amrex::Array4< amrex::Real const > &tau22, const amrex::Array4< amrex::Real const > &tau33, const amrex::Array4< amrex::Real const > &tau12, const amrex::Array4< amrex::Real const > &tau13, const amrex::Array4< amrex::Real const > &tau23, const amrex::Array4< const amrex::EBCellFlag > &c_cflag, const amrex::Array4< const amrex::EBCellFlag > &u_cflag, const amrex::Array4< const amrex::EBCellFlag > &v_cflag, const amrex::Array4< const amrex::EBCellFlag > &w_cflag)
Definition: ERF_EddyViscosity.H:110
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeVerticalShear2_EB(int i, int j, int k, const amrex::Array4< const amrex::EBCellFlag > &c_cflag, const amrex::Array4< const amrex::Real > &u_vfrac, const amrex::Array4< const amrex::Real > &v_vfrac, amrex::Real dzInv, const amrex::Array4< const amrex::Real > &u, const amrex::Array4< const amrex::Real > &v)
Definition: ERF_RichardsonNumber.H:303
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeN2_EB(int i, int j, int k, const amrex::Array4< const amrex::EBCellFlag > &c_cflag, amrex::Real dzInv, amrex::Real const_grav, const amrex::Array4< const amrex::Real > &cell_data, const MoistureComponentIndices &moisture_indices)
Definition: ERF_RichardsonNumber.H:175
eb_aux_ const * get_w_const_factory() const noexcept
Definition: ERF_EB.H:52
const std::unique_ptr< amrex::EBFArrayBoxFactory > & get_const_factory() const noexcept
Definition: ERF_EB.H:46
eb_aux_ const * get_v_const_factory() const noexcept
Definition: ERF_EB.H:51
eb_aux_ const * get_u_const_factory() const noexcept
Definition: ERF_EB.H:50
@ theta
Definition: ERF_MM5.H:20

Referenced by ComputeTurbulentViscosity().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeTurbulentViscosityRANS()

void ComputeTurbulentViscosityRANS ( Vector< std::unique_ptr< MultiFab >> &  ,
const MultiFab &  cons_in,
const MultiFab &  wdist,
MultiFab &  eddyViscosity,
MultiFab &  Hfx1,
MultiFab &  Hfx2,
MultiFab &  Hfx3,
MultiFab &  Diss,
const Geometry &  geom,
bool  use_terrain_fitted_coords,
Vector< std::unique_ptr< MultiFab >> &  ,
const std::unique_ptr< MultiFab > &  z_phys_nd,
const TurbChoice turbChoice,
const Real  const_grav,
std::unique_ptr< SurfaceLayer > &  SurfLayer,
const MultiFab *  z_0 
)

Function for computing the eddy viscosity with RANS.

Parameters
[in]Tau_levstrain at this level
[in]cons_incell center conserved quantities
[out]eddyViscosityturbulent viscosity
[in]Hfx1heat flux in x-dir
[in]Hfx2heat flux in y-dir
[in]Hfx3heat flux in z-dir
[in]Dissdissipation of turbulent kinetic energy
[in]geomproblem geometry
[in]mapfacmap factor
[in]turbChoicecontainer with turbulence parameters
579 {
580  const GpuArray<Real, AMREX_SPACEDIM> cellSizeInv = geom.InvCellSizeArray();
581  const bool use_SurfLayer = (SurfLayer != nullptr);
582 
583  Real inv_Pr_t = turbChoice.Pr_t_inv;
584  Real inv_Sc_t = turbChoice.Sc_t_inv;
585  Real inv_sigma_k = 1.0 / turbChoice.sigma_k;
586 
587  // One-Equation k model (Axell & Liungman 2001, Environ Fluid Mech)
588  //***********************************************************************************
589  if (turbChoice.rans_type == RANSType::kEqn)
590  {
591  const Real Cmu0 = turbChoice.Cmu0;
592  const Real Cmu0_pow3 = Cmu0 * Cmu0 * Cmu0;
593  const Real inv_Cb_sq = 1.0 / (turbChoice.Cb * turbChoice.Cb);
594  const Real Rt_crit = turbChoice.Rt_crit;
595  const Real Rt_min = turbChoice.Rt_min;
596  const Real l_g_max = turbChoice.l_g_max;
597  const Real abs_g = const_grav;
598 
599  const bool use_ref_theta = (turbChoice.theta_ref > 0);
600  const Real inv_theta0 = (use_ref_theta) ? 1.0 / turbChoice.theta_ref : 1.0;
601 
602 #ifdef _OPENMP
603 #pragma omp parallel if (Gpu::notInLaunchRegion())
604 #endif
605  for ( MFIter mfi(eddyViscosity,TilingIfNotGPU()); mfi.isValid(); ++mfi)
606  {
607  Box bxcc = mfi.tilebox();
608 
609  const Array4<Real const>& d_arr = wdist.const_array(mfi);
610  const Array4<Real const>& z0_arr = (use_SurfLayer) ? z_0->const_array(mfi) : Array4<Real const>{};
611 
612  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
613  const Array4<Real>& hfx_x = Hfx1.array(mfi);
614  const Array4<Real>& hfx_y = Hfx2.array(mfi);
615  const Array4<Real>& hfx_z = Hfx3.array(mfi);
616  const Array4<Real>& diss = Diss.array(mfi);
617 
618  const Array4<Real const>& cell_data = cons_in.array(mfi);
619 
620  const Array4<Real const>& z_nd_arr = z_phys_nd->const_array(mfi);
621 
622  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
623  {
625  Real tke = amrex::max(cell_data(i,j,k,RhoKE_comp)/cell_data(i,j,k,Rho_comp), eps);
626 
627  // Estimate stratification
628  Real dzInv = cellSizeInv[2];
629  if (use_terrain_fitted_coords) {
630  // the terrain grid is only deformed in z for now
631  dzInv /= Compute_h_zeta_AtCellCenter(i,j,k, cellSizeInv, z_nd_arr);
632  }
633  Real dtheta_dz = 0.5 * ( cell_data(i,j,k+1,RhoTheta_comp)/cell_data(i,j,k+1,Rho_comp)
634  - cell_data(i,j,k-1,RhoTheta_comp)/cell_data(i,j,k-1,Rho_comp) )*dzInv;
635  Real N2 = abs_g * inv_theta0 * dtheta_dz; // Brunt–Väisälä frequency squared
636  if (!use_ref_theta) {
637  // inv_theta0 == 1, divide by actual theta
638  N2 *= cell_data(i,j,k,Rho_comp) /
639  cell_data(i,j,k,RhoTheta_comp);
640  }
641 
642  // Geometric length scale (AL01, Eqn. 22)
643  Real l_g = (z0_arr) ? KAPPA * (d_arr(i, j, k) + z0_arr(i, j, 0))
644  : KAPPA * d_arr(i, j, k);
645 
646  // Enforce a maximum value
647  l_g = l_g_max * l_g / (l_g_max + l_g);
648 
649  // Turbulent length scale
650  Real length, Rt;
651  if (std::abs(N2) <= eps) {
652  length = l_g;
653  } else if (N2 > eps) {
654  // Stable (AL01, Eqn. 26)
655  length = std::sqrt(1.0 /
656  (1.0 / (l_g * l_g) + inv_Cb_sq * N2 / tke));
657  } else {
658  Real diss0 = Cmu0_pow3 * std::pow(tke,1.5) / l_g; // approx
659  Rt = tke*tke * N2 / (diss0*diss0);
660 
661  // Unstable (AL01, Eqn. 28)
662  // - predict
663  length = l_g * std::sqrt(1.0 - Cmu0_pow3*Cmu0_pow3 * inv_Cb_sq * Rt);
664  // - correct
665  diss0 = Cmu0_pow3 * std::pow(tke,1.5) / length;
666  Rt = tke*tke * N2 / (diss0*diss0);
667  length = l_g * std::sqrt(1.0 - Cmu0_pow3*Cmu0_pow3 * inv_Cb_sq * Rt);
668  }
669  mu_turb(i, j, k, EddyDiff::Turb_lengthscale) = length;
670 
671  // Dissipation rate (AL01, Eqn. 19)
672  diss(i, j, k) = cell_data(i, j, k, Rho_comp) * Cmu0_pow3 * std::pow(tke,1.5) / length;
673 
674  // Turbulent Richardson number (AL01, Eqn. 29)
675  //Real Rt = tke*tke * N2 / (diss(i,j,k)*diss(i,j,k));
676  Rt = length*length * N2 / (tke * Cmu0_pow3 * Cmu0_pow3); // combined with Eqn. 19
677 
678  // Burchard & Petersen smoothing function
679  Rt = (Rt >= Rt_crit) ? Rt
680  : std::max(Rt,
681  Rt - std::pow(Rt - Rt_crit, 2) /
682  (Rt + Rt_min - 2*Rt_crit));
683 
684  // Stability functions
685  // Note: These use the smoothed turbulent Richardson number
686  Real cmu = (Cmu0 + 0.108*Rt)
687  / (1.0 + 0.308*Rt + 0.00837*Rt*Rt); // (AL01, Eqn. 31)
688  Real cmu_prime = Cmu0 / (1 + 0.277*Rt); // (AL01, Eqn. 32)
689 
690  // Calculate eddy diffusivities
691  // K = rho * nu_t = rho * c_mu * tke^(1/2) * length
692  Real nut = cmu * std::sqrt(tke) * length; // eddy viscosity
693  Real nut_prime = cmu_prime / cmu * nut; // eddy diffusivity
694  mu_turb(i, j, k, EddyDiff::Mom_h) = cell_data(i, j, k, Rho_comp) * nut;
695  mu_turb(i, j, k, EddyDiff::Mom_v) = mu_turb(i, j, k, EddyDiff::Mom_h);
696  mu_turb(i, j, k, EddyDiff::Theta_v) = cell_data(i, j, k, Rho_comp) * nut_prime;
697 
698  // Calculate heat flux
699  // - If using SurfaceLayer, the value at k=0 will be overwritten
700  hfx_x(i, j, k) = 0.0;
701  hfx_y(i, j, k) = 0.0;
702  // Note: buoyant production = g/theta0 * hfx == -nut_prime * N^2 (c.f. AL01 Eqn. 15)
703  // = nut_prime * g/theta0 * dtheta/dz
704  // ==> hfx = nut_prime * dtheta/dz
705  // Our convention is such that dtheta/dz < 0 gives a positive
706  // (upward) heat flux.
707  hfx_z(i, j, k) = -mu_turb(i, j, k, EddyDiff::Theta_v) * dtheta_dz; // (rho*w)' theta' [kg m^-2 s^-1 K]
708  });
709  }
710  }
711 
712  // Extrapolate Kturb in x/y, fill remaining elements (relevant to lev==0)
713  //***********************************************************************************
714  int ngc(1);
715  // EddyDiff mapping : Theta_h KE_h Scalar_h Q_h
716  Vector<Real> Factors = {inv_Pr_t, inv_sigma_k, inv_Sc_t, inv_Sc_t}; // alpha = mu/Pr
717  Gpu::AsyncVector<Real> d_Factors; d_Factors.resize(Factors.size());
718  Gpu::copy(Gpu::hostToDevice, Factors.begin(), Factors.end(), d_Factors.begin());
719  Real* fac_ptr = d_Factors.data();
720 
721  const bool use_KE = ( turbChoice.rans_type == RANSType::kEqn );
722 
723 #ifdef _OPENMP
724 #pragma omp parallel if (Gpu::notInLaunchRegion())
725 #endif
726  for ( MFIter mfi(eddyViscosity,TilingIfNotGPU()); mfi.isValid(); ++mfi)
727  {
728  Box bxcc = mfi.tilebox();
729  Box planex = bxcc; planex.setSmall(0, 1); planex.setBig(0, ngc); planex.grow(1,1);
730  Box planey = bxcc; planey.setSmall(1, 1); planey.setBig(1, ngc); planey.grow(0,1);
731  bxcc.growLo(0,ngc); bxcc.growHi(0,ngc);
732  bxcc.growLo(1,ngc); bxcc.growHi(1,ngc);
733 
734  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
735 
736  for (auto n = 1; n < (EddyDiff::NumDiffs-1)/2; ++n) {
737  int offset = (EddyDiff::NumDiffs-1)/2;
738  switch (n)
739  {
740  case EddyDiff::KE_h:
741  if (use_KE) {
742  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
743  {
744  int indx = n;
745  int indx_v = indx + offset;
746  mu_turb(i,j,k,indx) = mu_turb(i,j,k,EddyDiff::Mom_h) * fac_ptr[indx-1];
747  mu_turb(i,j,k,indx_v) = mu_turb(i,j,k,indx);
748  });
749  }
750  break;
751  default:
752  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
753  {
754  int indx = n;
755  int indx_v = indx + offset;
756 
757  mu_turb(i,j,k,indx) = mu_turb(i,j,k,EddyDiff::Mom_h) * fac_ptr[indx-1];
758 
759  // NOTE: Theta_v has already been set for Deardorff
760  if (!(indx_v == EddyDiff::Theta_v && use_KE)) {
761  mu_turb(i,j,k,indx_v) = mu_turb(i,j,k,indx);
762  }
763  });
764  break;
765  }
766  }
767  }
768 }
constexpr amrex::Real KAPPA
Definition: ERF_Constants.H:20
amrex::Real Rt_min
Definition: ERF_TurbStruct.H:394
amrex::Real Cmu0
Definition: ERF_TurbStruct.H:391
amrex::Real Cb
Definition: ERF_TurbStruct.H:392
amrex::Real Rt_crit
Definition: ERF_TurbStruct.H:393
amrex::Real l_g_max
Definition: ERF_TurbStruct.H:395

Referenced by ComputeTurbulentViscosity().

Here is the call graph for this function:
Here is the caller graph for this function: