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
810 {
811  BL_PROFILE_VAR("ComputeTurbulentViscosity()",ComputeTurbulentViscosity);
812  //
813  // In LES mode, the turbulent viscosity is isotropic (unless mix_isotropic is set to false), so
814  // the LES model sets both horizontal and vertical viscosities
815  //
816  // In PBL mode, the primary purpose of the PBL model is to control vertical transport, so the PBL model sets the vertical viscosity.
817  // Optionally, the PBL model can be run in conjunction with an LES model that sets the horizontal viscosity
818  // (this isn’t truly LES, but the model form is the same as Smagorinsky).
819  //
820  // ComputeTurbulentViscosityLES populates the LES viscosity for both horizontal and vertical components.
821  // ComputeTurbulentViscosityPBL computes the PBL viscosity just for the vertical component.
822  //
823 
824  TurbChoice turbChoice = solverChoice.turbChoice[level];
825  const Real const_grav = solverChoice.gravity;
826 
827  if (!SurfLayer) {
828  AMREX_ALWAYS_ASSERT(!vert_only);
829  }
830 
831  bool impose_phys_bcs = false;
832 
833  if (turbChoice.les_type != LESType::None) {
834  impose_phys_bcs = true;
835  if (solverChoice.terrain_type == TerrainType::EB) {
837  cons_in, eddyViscosity,
838  Hfx1, Hfx2, Hfx3,
839  geom,
840  mapfac, turbChoice, const_grav,
841  solverChoice,
842  SurfLayer, solverChoice.moisture_indices,
843  ebfact, &xvel, &yvel);
844  } else {
846  cons_in, eddyViscosity,
847  Hfx1, Hfx2, Hfx3, Diss,
848  geom, use_terrain_fitted_coords,
849  mapfac, z_phys_nd, turbChoice, const_grav,
850  SurfLayer, solverChoice.moisture_indices,
851  &xvel, &yvel);
852  }
853  }
854 
855  if (turbChoice.rans_type != RANSType::None) {
856  impose_phys_bcs = true;
858  cons_in, wdist,
859  eddyViscosity,
860  Hfx1, Hfx2, Hfx3, Diss,
861  geom, use_terrain_fitted_coords,
862  mapfac, z_phys_nd, turbChoice, const_grav,
863  SurfLayer, z_0);
864  }
865 
866  if (turbChoice.pbl_type == PBLType::MYJ) {
867  ComputeDiffusivityMYJ(dt, xvel, yvel, cons_in, eddyViscosity,
868  geom, turbChoice, SurfLayer,
869  use_terrain_fitted_coords, use_moisture,
870  level, bc_ptr, vert_only, z_phys_nd,
871  solverChoice.moisture_indices);
872  } else if (turbChoice.pbl_type == PBLType::MYNN25) {
873  ComputeDiffusivityMYNN25(xvel, yvel, cons_in, eddyViscosity,
874  geom, turbChoice, SurfLayer,
875  use_terrain_fitted_coords, use_moisture,
876  level, bc_ptr, vert_only, z_phys_nd,
877  solverChoice.moisture_indices);
878  } else if (turbChoice.pbl_type == PBLType::MYNNEDMF) {
879  ComputeDiffusivityMYNNEDMF(xvel, yvel, cons_in, eddyViscosity,
880  geom, turbChoice, SurfLayer,
881  use_terrain_fitted_coords, use_moisture,
882  level, bc_ptr, vert_only, z_phys_nd,
883  solverChoice.moisture_indices);
884  } else if (turbChoice.pbl_type == PBLType::YSU) {
885  ComputeDiffusivityYSU(xvel, yvel, cons_in, eddyViscosity,
886  geom, turbChoice, SurfLayer,
887  use_terrain_fitted_coords, use_moisture,
888  level, bc_ptr, vert_only, z_phys_nd,
889  solverChoice.moisture_indices);
890  } else if (turbChoice.pbl_type == PBLType::MRF) {
891  ComputeDiffusivityMRF(xvel, yvel, cons_in, eddyViscosity,
892  geom, turbChoice, SurfLayer,
893  use_terrain_fitted_coords, use_moisture,
894  level, bc_ptr, vert_only, z_phys_nd,
895  solverChoice.moisture_indices);
896  } else if (turbChoice.pbl_type == PBLType::SHOC) {
897  // NOTE: Nothing to do here. The SHOC class handles setting the vertical
898  // components of eddyDiffs in slow RHS pre.
899  }
900 
901  //
902  // At all levels we need to fill values outside the physical boundary for the LES coeffs.
903  // In addition, for all cases, if at level > 0, we want to fill fine ghost cell values that
904  // overlie coarse grid cells (and that are not in another fine valid region) with
905  // extrapolated values from the interior, rather than interpolating from the coarser level,
906  // since we may be using a different turbulence model there.
907  //
908  // Note: here "covered" refers to "covered by valid region of another grid at this level"
909  // Note: here "physbnd" refers to "cells outside the domain if not periodic"
910  // Note: here "interior" refers to "valid cells, i.e. inside 'my' grid"
911  //
912  {
913  int is_covered = 0;
914  int is_notcovered = 1;
915  int is_physbnd = 2;
916  int is_interior = 3;
917  iMultiFab cc_mask(eddyViscosity.boxArray(),eddyViscosity.DistributionMap(),1,1);
918  cc_mask.BuildMask(geom.Domain(), geom.periodicity(), is_covered, is_notcovered, is_physbnd, is_interior);
919 
920  Box domain = geom.Domain();
921  for (int i = 0; i < AMREX_SPACEDIM; ++i) {
922  if (geom.isPeriodic(i)) {
923  domain.grow(i,1);
924  }
925  }
926 
927  eddyViscosity.FillBoundary(geom.periodicity());
928 
929  int ncomp = eddyViscosity.nComp();
930 
931 #ifdef _OPENMP
932 #pragma omp parallel if (Gpu::notInLaunchRegion())
933 #endif
934  for (MFIter mfi(eddyViscosity); mfi.isValid(); ++mfi)
935  {
936  Box vbx = mfi.validbox();
937 
938  Box planex_lo = mfi.growntilebox(1); planex_lo.setBig(0, vbx.smallEnd(0)-1);
939  Box planey_lo = mfi.growntilebox(1); planey_lo.setBig(1, vbx.smallEnd(1)-1);
940  Box planez_lo = mfi.growntilebox(1); planez_lo.setBig(2, vbx.smallEnd(2)-1);
941 
942  Box planex_hi = mfi.growntilebox(1); planex_hi.setSmall(0, vbx.bigEnd(0)+1);
943  Box planey_hi = mfi.growntilebox(1); planey_hi.setSmall(1, vbx.bigEnd(1)+1);
944  Box planez_hi = mfi.growntilebox(1); planez_hi.setSmall(2, vbx.bigEnd(2)+1);
945 
946  int i_lo = vbx.smallEnd(0); int i_hi = vbx.bigEnd(0);
947  int j_lo = vbx.smallEnd(1); int j_hi = vbx.bigEnd(1);
948  int k_lo = vbx.smallEnd(2); int k_hi = vbx.bigEnd(2);
949 
950  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
951  const Array4<int>& mask_arr = cc_mask.array(mfi);
952 
953  auto domlo = lbound(domain);
954  auto domhi = ubound(domain);
955 
956  ParallelFor(planex_lo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
957  {
958  int lj = amrex::min(amrex::max(j, domlo.y), domhi.y);
959  int lk = amrex::min(amrex::max(k, domlo.z), domhi.z);
960  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i_lo,lj,lk) != is_notcovered) ||
961  (mask_arr(i,j,k) == is_physbnd && i < domlo.x && impose_phys_bcs)) {
962  for (int n = 0; n < ncomp; n++) {
963  mu_turb(i,j,k,n) = mu_turb(i_lo,lj,lk,n);
964  }
965  }
966  });
967  ParallelFor(planex_hi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
968  {
969  int lj = amrex::min(amrex::max(j, domlo.y), domhi.y);
970  int lk = amrex::min(amrex::max(k, domlo.z), domhi.z);
971  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i_hi,lj,lk) != is_notcovered) ||
972  (mask_arr(i,j,k) == is_physbnd && i > domhi.x && impose_phys_bcs)) {
973  for (int n = 0; n < ncomp; n++) {
974  mu_turb(i,j,k,n) = mu_turb(i_hi,lj,lk,n);
975  }
976  }
977  });
978  ParallelFor(planey_lo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
979  {
980  int lk = amrex::min(amrex::max(k, domlo.z), domhi.z);
981  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i,j_lo,lk) != is_notcovered) ||
982  (mask_arr(i,j,k) == is_physbnd && j < domlo.y && impose_phys_bcs)) {
983  for (int n = 0; n < ncomp; n++) {
984  mu_turb(i,j,k,n) = mu_turb(i,j_lo,lk,n);
985  }
986  }
987  });
988  ParallelFor(planey_hi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
989  {
990  int lk = amrex::min(amrex::max(k, domlo.z), domhi.z);
991  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i,j_hi,lk) != is_notcovered)||
992  (mask_arr(i,j,k) == is_physbnd && j > domhi.y && impose_phys_bcs)) {
993  for (int n = 0; n < ncomp; n++) {
994  mu_turb(i,j,k,n) = mu_turb(i,j_hi,lk,n);
995  }
996  }
997  });
998  ParallelFor(planez_lo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
999  {
1000  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i,j,k_lo) != is_notcovered) ||
1001  (mask_arr(i,j,k) == is_physbnd && k < domlo.z && impose_phys_bcs)) {
1002  if (mask_arr(i,j,k) == is_physbnd) {
1003  }
1004  for (int n = 0; n < ncomp; n++) {
1005  mu_turb(i,j,k,n) = mu_turb(i,j,k_lo,n);
1006  }
1007  }
1008  });
1009  ParallelFor(planez_hi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
1010  {
1011  if ((mask_arr(i,j,k) == is_notcovered && mask_arr(i,j,k_hi) != is_notcovered) ||
1012  (mask_arr(i,j,k) == is_physbnd && k > domhi.z && impose_phys_bcs)) {
1013  for (int n = 0; n < ncomp; n++) {
1014  mu_turb(i,j,k,n) = mu_turb(i,j,k_hi,n);
1015  }
1016  }
1017  });
1018  } // mfi
1019  eddyViscosity.FillBoundary(geom.periodicity());
1020  }
1021 }
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:569
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:791
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+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);} } })
for(int i=0;i< m_num_species;i++)
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:48
TurbChoice turbChoice
Definition: ERF_SetupVertDiff.H:6
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Real z_0
Definition: ERF_UpdateWSubsidence_Bomex.H:10
@ xvel
Definition: ERF_IndexDefines.H:159
@ yvel
Definition: ERF_IndexDefines.H:160
amrex::Real gravity
Definition: ERF_DataStruct.H:1154
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:1091
static TerrainType terrain_type
Definition: ERF_DataStruct.H:1067
MoistureComponentIndices moisture_indices
Definition: ERF_DataStruct.H:1223
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 = one / 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  // one 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(two * SmnSmn);
104 
105  // =====================================================================
106  // two 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 = one / (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 = one / dzInv;
123  DeltaH = std::sqrt(one / (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 = one;
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 = myhalf * ( 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) = zero;
154  hfx_y(i,j,k) = zero;
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(zero, l_C_e - Real(1.9)*l_C_k);
167  const Real l_abs_g = const_grav;
168 
169  const bool use_ref_theta = (turbChoice.theta_ref > 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 = one / (dxInv * mf_u(i,j,0) * dyInv * mf_v(i,j,0) * dzInv);
204  Delta = std::cbrt(cellVolMsf);
205  } else {
206  Delta = one / dzInv;
207  }
208 
209  Real dtheta_dz;
210  if (use_thetav_grad) {
211  dtheta_dz = myhalf * ( 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 = myhalf * ( 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 = myhalf * ( 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 = Real(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, Real(0.001) * Delta);
242  }
243 
244  Real DeltaH = (isotropic) ? length : std::sqrt(one / (dxInv * mf_u(i,j,0) * dyInv * mf_v(i,j,0)));
245 
246  Real Pr_inv_v = (one + two*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 = Real(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,Real(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) = zero;
273  hfx_y(i,j,k) = zero;
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
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
@ 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:17
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:55
@ Theta_v
Definition: ERF_IndexDefines.H:194
@ Turb_lengthscale
Definition: ERF_IndexDefines.H:198
@ Mom_h
Definition: ERF_IndexDefines.H:188
@ Mom_v
Definition: ERF_IndexDefines.H:193
@ NumDiffs
Definition: ERF_IndexDefines.H:199
@ Theta_h
Definition: ERF_IndexDefines.H:189
@ KE_h
Definition: ERF_IndexDefines.H:190
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 = one / 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  if (!c_cflag(i,j,k).isCovered()) {
410  // =====================================================================
411  // 1. STRAIN RATE MAGNITUDE CALCULATION
412  // =====================================================================
413  Real SmnSmn = zero;
414  if (c_cflag(i,j,k).isRegular()) {
415  SmnSmn = ComputeSmnSmn(i,j,k,tau11,tau22,tau33,tau12,tau13,tau23);
416  } else if (c_cflag(i,j,k).isSingleValued()) {
417  SmnSmn = ComputeSmnSmn_EB(i,j,k,tau11,tau22,tau33,tau12,tau13,tau23,c_cflag,u_cflag,v_cflag,w_cflag);
418  }
419  Real strain_rate_magnitude = std::sqrt(two * SmnSmn);
420 
421  // =====================================================================
422  // 2. GRID SCALE CALCULATION (filter width Δ)
423  // =====================================================================
424  Real Delta;
425  Real DeltaH;
426  if (isotropic) {
427  Real cellVolMsf = one / (dxInv * mf_u(i,j,0) * dyInv * mf_v(i,j,0) * dzInv);
428  Delta = std::cbrt(cellVolMsf);
429  DeltaH = Delta;
430  } else {
431  Delta = one / dzInv;
432  DeltaH = std::sqrt(one / (dxInv * mf_u(i,j,0) * dyInv * mf_v(i,j,0)));
433  }
434 
435  Real rho = cell_data(i, j, k, Rho_comp);
436  Real CsDeltaSqr_h = Cs * Cs * DeltaH * DeltaH;
437  Real CsDeltaSqr_v = Cs * Cs * Delta * Delta;
438 
439  Real nu_turb_base_h = CsDeltaSqr_h * strain_rate_magnitude;
440  Real nu_turb_base_v = CsDeltaSqr_v * strain_rate_magnitude;
441 
442  Real stability_factor = one;
443 
444  if (l_use_Ri_corr && l_has_xvel && l_has_yvel) {
445  Real N2 = zero;
446  Real S2_vert = zero;
447  if (c_cflag(i,j,k).isRegular()) {
448  N2 = ComputeN2(i, j, k, dzInv, l_abs_g, cell_data, moisture_indices);
449  S2_vert = ComputeVerticalShear2(i, j, k, dzInv, u_arr, v_arr);
450  } else if (c_cflag(i,j,k).isSingleValued()) {
451  N2 = ComputeN2_EB(i, j, k, c_cflag, dzInv, l_abs_g, cell_data, moisture_indices);
452  S2_vert = ComputeVerticalShear2_EB(i, j, k, c_cflag, u_vfrac, v_vfrac, dzInv, u_arr, v_arr);
453  }
454  Real Ri = ComputeRichardson(N2, S2_vert);
455  stability_factor = StabilityFunction(Ri, l_Ri_crit);
456  }
457 
458  if (isotropic) {
459  mu_turb(i, j, k, EddyDiff::Mom_h) = rho * nu_turb_base_h * stability_factor;
460  mu_turb(i, j, k, EddyDiff::Mom_v) = rho * nu_turb_base_v * stability_factor;
461  } else {
462  mu_turb(i, j, k, EddyDiff::Mom_h) = rho * nu_turb_base_h;
463  mu_turb(i, j, k, EddyDiff::Mom_v) = rho * nu_turb_base_v * stability_factor;
464  }
465 
466  Real dtheta_dz = zero;
467 
468  if (c_cflag(i,j,k).isSingleValued() && c_cflag(i,j,k+1).isCovered()) {
469  amrex::Real theta = cell_data(i, j, k, RhoTheta_comp) / rho;
470  amrex::Real theta_km1 = cell_data(i,j,k-1,RhoTheta_comp)/cell_data(i,j,k-1,Rho_comp);
471  amrex::Real theta_km2 = cell_data(i,j,k-2,RhoTheta_comp)/cell_data(i,j,k-2,Rho_comp);
472  dtheta_dz = (three*theta - Real(4.)*theta_km1 + theta_km2) * myhalf * dzInv;
473  } else if (c_cflag(i,j,k).isSingleValued() && c_cflag(i,j,k-1).isCovered()) {
474  amrex::Real theta = cell_data(i, j, k, RhoTheta_comp) / rho;
475  amrex::Real theta_kp1 = cell_data(i,j,k+1,RhoTheta_comp)/cell_data(i,j,k+1,Rho_comp);
476  amrex::Real theta_kp2 = cell_data(i,j,k+2,RhoTheta_comp)/cell_data(i,j,k+2,Rho_comp);
477  dtheta_dz = (-theta_kp2 + Real(4.)*theta_kp1 - three*theta) * myhalf * dzInv;
478  } else if (!c_cflag(i,j,k).isCovered()) {
479  amrex::Real theta_km1 = cell_data(i,j,k-1,RhoTheta_comp)/cell_data(i,j,k-1,Rho_comp);
480  amrex::Real theta_kp1 = cell_data(i,j,k+1,RhoTheta_comp)/cell_data(i,j,k+1,Rho_comp);
481  dtheta_dz = myhalf * (theta_kp1 - theta_km1) * dzInv;
482  }
483 
484  hfx_x(i,j,k) = zero;
485  hfx_y(i,j,k) = zero;
486  hfx_z(i,j,k) = -inv_Pr_t * mu_turb(i,j,k,EddyDiff::Mom_v) * dtheta_dz;
487 
488  } else {
489 
490  hfx_x(i,j,k) = zero;
491  hfx_y(i,j,k) = zero;
492  hfx_z(i,j,k) = zero;
493  }// End of if not covered
494  });
495  }
496  }
497 
498  // Extrapolate Kturb in x/y, fill remaining elements (relevant to lev==0)
499  //***********************************************************************************
500  int ngc(1);
501  // EddyDiff mapping : Theta_h KE_h Scalar_h Q_h
502  Vector<Real> Factors = {inv_Pr_t, inv_sigma_k, inv_Sc_t, inv_Sc_t}; // alpha = mu/Pr
503  Gpu::AsyncVector<Real> d_Factors; d_Factors.resize(Factors.size());
504  Gpu::copy(Gpu::hostToDevice, Factors.begin(), Factors.end(), d_Factors.begin());
505  Real* fac_ptr = d_Factors.data();
506 
507  const bool use_KE = ( turbChoice.les_type == LESType::Deardorff );
508 
509 #ifdef _OPENMP
510 #pragma omp parallel if (Gpu::notInLaunchRegion())
511 #endif
512  for ( MFIter mfi(eddyViscosity,TilingIfNotGPU()); mfi.isValid(); ++mfi)
513  {
514  Box bxcc = mfi.tilebox();
515  Box planex = bxcc; planex.setSmall(0, 1); planex.setBig(0, ngc); planex.grow(1,1);
516  Box planey = bxcc; planey.setSmall(1, 1); planey.setBig(1, ngc); planey.grow(0,1);
517  bxcc.growLo(0,ngc); bxcc.growHi(0,ngc);
518  bxcc.growLo(1,ngc); bxcc.growHi(1,ngc);
519 
520  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
521 
522  for (auto n = 1; n < (EddyDiff::NumDiffs-1)/2; ++n) {
523  int offset = (EddyDiff::NumDiffs-1)/2;
524  switch (n)
525  {
526  case EddyDiff::KE_h:
527  if (use_KE) {
528  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
529  {
530  int indx = n;
531  int indx_v = indx + offset;
532  mu_turb(i,j,k,indx) = mu_turb(i,j,k,EddyDiff::Mom_h) * fac_ptr[indx-1];
533  mu_turb(i,j,k,indx_v) = mu_turb(i,j,k,indx);
534  });
535  }
536  break;
537  default:
538  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
539  {
540  int indx = n;
541  int indx_v = indx + offset;
542 
543  // NOTE: Theta_h, Theta_v have already been set for Deardorff
544  if (!(indx_v == EddyDiff::Theta_v && use_KE)) {
545  mu_turb(i,j,k,indx) = mu_turb(i,j,k,EddyDiff::Mom_h) * fac_ptr[indx-1];
546  mu_turb(i,j,k,indx_v) = mu_turb(i,j,k,indx);
547  }
548  });
549  break;
550  }
551  }
552  }
553 }
constexpr amrex::Real three
Definition: ERF_Constants.H:9
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
585 {
586  const GpuArray<Real, AMREX_SPACEDIM> cellSizeInv = geom.InvCellSizeArray();
587  const bool use_SurfLayer = (SurfLayer != nullptr);
588 
589  Real inv_Pr_t = turbChoice.Pr_t_inv;
590  Real inv_Sc_t = turbChoice.Sc_t_inv;
591  Real inv_sigma_k = one / turbChoice.sigma_k;
592 
593  // One-Equation k model (Axell & Liungman 2001, Environ Fluid Mech)
594  //***********************************************************************************
595  if (turbChoice.rans_type == RANSType::kEqn)
596  {
597  const Real Cmu0 = turbChoice.Cmu0;
598  const Real Cmu0_pow3 = Cmu0 * Cmu0 * Cmu0;
599  const Real inv_Cb_sq = one / (turbChoice.Cb * turbChoice.Cb);
600  const Real Rt_crit = turbChoice.Rt_crit;
601  const Real Rt_min = turbChoice.Rt_min;
602  const Real l_g_max = turbChoice.l_g_max;
603  const Real abs_g = const_grav;
604 
605  const bool use_ref_theta = (turbChoice.theta_ref > 0);
606  const Real inv_theta0 = (use_ref_theta) ? one / turbChoice.theta_ref : one;
607 
608 #ifdef _OPENMP
609 #pragma omp parallel if (Gpu::notInLaunchRegion())
610 #endif
611  for ( MFIter mfi(eddyViscosity,TilingIfNotGPU()); mfi.isValid(); ++mfi)
612  {
613  Box bxcc = mfi.tilebox();
614 
615  const Array4<Real const>& d_arr = wdist.const_array(mfi);
616  const Array4<Real const>& z0_arr = (use_SurfLayer) ? z_0->const_array(mfi) : Array4<Real const>{};
617 
618  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
619  const Array4<Real>& hfx_x = Hfx1.array(mfi);
620  const Array4<Real>& hfx_y = Hfx2.array(mfi);
621  const Array4<Real>& hfx_z = Hfx3.array(mfi);
622  const Array4<Real>& diss = Diss.array(mfi);
623 
624  const Array4<Real const>& cell_data = cons_in.array(mfi);
625 
626  const Array4<Real const>& z_nd_arr = z_phys_nd->const_array(mfi);
627 
628  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
629  {
631  Real tke = amrex::max(cell_data(i,j,k,RhoKE_comp)/cell_data(i,j,k,Rho_comp), eps);
632 
633  // Estimate stratification
634  Real dzInv = cellSizeInv[2];
635  if (use_terrain_fitted_coords) {
636  // the terrain grid is only deformed in z for now
637  dzInv /= Compute_h_zeta_AtCellCenter(i,j,k, cellSizeInv, z_nd_arr);
638  }
639  Real dtheta_dz = myhalf * ( cell_data(i,j,k+1,RhoTheta_comp)/cell_data(i,j,k+1,Rho_comp)
640  - cell_data(i,j,k-1,RhoTheta_comp)/cell_data(i,j,k-1,Rho_comp) )*dzInv;
641  Real N2 = abs_g * inv_theta0 * dtheta_dz; // Brunt–Väisälä frequency squared
642  if (!use_ref_theta) {
643  // inv_theta0 == 1, divide by actual theta
644  N2 *= cell_data(i,j,k,Rho_comp) /
645  cell_data(i,j,k,RhoTheta_comp);
646  }
647 
648  // Geometric length scale (AL01, Eqn. 22)
649  Real l_g = (z0_arr) ? KAPPA * (d_arr(i, j, k) + z0_arr(i, j, 0))
650  : KAPPA * d_arr(i, j, k);
651 
652  // Enforce a maximum value
653  l_g = l_g_max * l_g / (l_g_max + l_g);
654 
655  // Turbulent length scale
656  Real length, Rt;
657  if (std::abs(N2) <= eps) {
658  length = l_g;
659  } else if (N2 > eps) {
660  // Stable (AL01, Eqn. 26)
661  length = std::sqrt(one /
662  (one / (l_g * l_g) + inv_Cb_sq * N2 / tke));
663  } else {
664  Real diss0 = Cmu0_pow3 * std::pow(tke,Real(1.5)) / l_g; // approx
665  Rt = tke*tke * N2 / (diss0*diss0);
666 
667  // Unstable (AL01, Eqn. 28)
668  // - predict
669  length = l_g * std::sqrt(one - Cmu0_pow3*Cmu0_pow3 * inv_Cb_sq * Rt);
670  // - correct
671  diss0 = Cmu0_pow3 * std::pow(tke,Real(1.5)) / length;
672  Rt = tke*tke * N2 / (diss0*diss0);
673  length = l_g * std::sqrt(one - Cmu0_pow3*Cmu0_pow3 * inv_Cb_sq * Rt);
674  }
675  mu_turb(i, j, k, EddyDiff::Turb_lengthscale) = length;
676 
677  // Dissipation rate (AL01, Eqn. 19)
678  diss(i, j, k) = cell_data(i, j, k, Rho_comp) * Cmu0_pow3 * std::pow(tke,Real(1.5)) / length;
679 
680  // Turbulent Richardson number (AL01, Eqn. 29)
681  //Real Rt = tke*tke * N2 / (diss(i,j,k)*diss(i,j,k));
682  Rt = length*length * N2 / (tke * Cmu0_pow3 * Cmu0_pow3); // combined with Eqn. 19
683 
684  // Burchard & Petersen smoothing function
685  Rt = (Rt >= Rt_crit) ? Rt : std::max(Rt, Rt - (Rt - Rt_crit)*(Rt - Rt_crit) / (Rt + Rt_min - 2*Rt_crit));
686 
687  // Stability functions
688  // Note: These use the smoothed turbulent Richardson number
689  Real cmu = (Cmu0 + Real(0.108)*Rt)
690  / (one + Real(0.308)*Rt + Real(0.00837)*Rt*Rt); // (AL01, Eqn. 31)
691  Real cmu_prime = Cmu0 / (1 + Real(0.277)*Rt); // (AL01, Eqn. 32)
692 
693  // Calculate eddy diffusivities
694  // K = rho * nu_t = rho * c_mu * tke^(1/2) * length
695  Real nut = cmu * std::sqrt(tke) * length; // eddy viscosity
696  Real nut_prime = cmu_prime / cmu * nut; // eddy diffusivity
697  mu_turb(i, j, k, EddyDiff::Mom_h) = cell_data(i, j, k, Rho_comp) * nut;
698  mu_turb(i, j, k, EddyDiff::Mom_v) = mu_turb(i, j, k, EddyDiff::Mom_h);
699  mu_turb(i, j, k, EddyDiff::Theta_v) = cell_data(i, j, k, Rho_comp) * nut_prime;
700 
701  // Calculate heat flux
702  // - If using SurfaceLayer, the value at k=0 will be overwritten
703  hfx_x(i, j, k) = zero;
704  hfx_y(i, j, k) = zero;
705  // Note: buoyant production = g/theta0 * hfx == -nut_prime * N^2 (c.f. AL01 Eqn. 15)
706  // = nut_prime * g/theta0 * dtheta/dz
707  // ==> hfx = nut_prime * dtheta/dz
708  // Our convention is such that dtheta/dz < 0 gives a positive
709  // (upward) heat flux.
710  hfx_z(i, j, k) = -mu_turb(i, j, k, EddyDiff::Theta_v) * dtheta_dz; // (rho*w)' theta' [kg m^-2 s^-1 K]
711  });
712  }
713  }
714 
715  // Extrapolate Kturb in x/y, fill remaining elements (relevant to lev==0)
716  //***********************************************************************************
717  int ngc(1);
718  // EddyDiff mapping : Theta_h KE_h Scalar_h Q_h
719  Vector<Real> Factors = {inv_Pr_t, inv_sigma_k, inv_Sc_t, inv_Sc_t}; // alpha = mu/Pr
720  Gpu::AsyncVector<Real> d_Factors; d_Factors.resize(Factors.size());
721  Gpu::copy(Gpu::hostToDevice, Factors.begin(), Factors.end(), d_Factors.begin());
722  Real* fac_ptr = d_Factors.data();
723 
724  const bool use_KE = ( turbChoice.rans_type == RANSType::kEqn );
725 
726 #ifdef _OPENMP
727 #pragma omp parallel if (Gpu::notInLaunchRegion())
728 #endif
729  for ( MFIter mfi(eddyViscosity,TilingIfNotGPU()); mfi.isValid(); ++mfi)
730  {
731  Box bxcc = mfi.tilebox();
732  Box planex = bxcc; planex.setSmall(0, 1); planex.setBig(0, ngc); planex.grow(1,1);
733  Box planey = bxcc; planey.setSmall(1, 1); planey.setBig(1, ngc); planey.grow(0,1);
734  bxcc.growLo(0,ngc); bxcc.growHi(0,ngc);
735  bxcc.growLo(1,ngc); bxcc.growHi(1,ngc);
736 
737  const Array4<Real>& mu_turb = eddyViscosity.array(mfi);
738 
739  for (auto n = 1; n < (EddyDiff::NumDiffs-1)/2; ++n) {
740  int offset = (EddyDiff::NumDiffs-1)/2;
741  switch (n)
742  {
743  case EddyDiff::KE_h:
744  if (use_KE) {
745  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
746  {
747  int indx = n;
748  int indx_v = indx + offset;
749  mu_turb(i,j,k,indx) = mu_turb(i,j,k,EddyDiff::Mom_h) * fac_ptr[indx-1];
750  mu_turb(i,j,k,indx_v) = mu_turb(i,j,k,indx);
751  });
752  }
753  break;
754  default:
755  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
756  {
757  int indx = n;
758  int indx_v = indx + offset;
759 
760  mu_turb(i,j,k,indx) = mu_turb(i,j,k,EddyDiff::Mom_h) * fac_ptr[indx-1];
761 
762  // NOTE: Theta_v has already been set for Deardorff
763  if (!(indx_v == EddyDiff::Theta_v && use_KE)) {
764  mu_turb(i,j,k,indx_v) = mu_turb(i,j,k,indx);
765  }
766  });
767  break;
768  }
769  }
770  }
771 }
constexpr amrex::Real KAPPA
Definition: ERF_Constants.H:30
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: