ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ERF_ComputeDiffusivityYSU.cpp File Reference
#include "ERF_ABLMost.H"
#include "ERF_DirectionSelector.H"
#include "ERF_Diffusion.H"
#include "ERF_Constants.H"
#include "ERF_TurbStruct.H"
#include "ERF_PBLModels.H"
Include dependency graph for ERF_ComputeDiffusivityYSU.cpp:

Functions

void ComputeDiffusivityYSU (const MultiFab &xvel, const MultiFab &yvel, const MultiFab &cons_in, MultiFab &eddyViscosity, const Geometry &geom, const TurbChoice &turbChoice, std::unique_ptr< ABLMost > &most, bool use_terrain_fitted_coords, bool, int level, const BCRec *bc_ptr, bool, const std::unique_ptr< MultiFab > &z_phys_nd)
 

Function Documentation

◆ ComputeDiffusivityYSU()

void ComputeDiffusivityYSU ( const MultiFab &  xvel,
const MultiFab &  yvel,
const MultiFab &  cons_in,
MultiFab &  eddyViscosity,
const Geometry &  geom,
const TurbChoice turbChoice,
std::unique_ptr< ABLMost > &  most,
bool  use_terrain_fitted_coords,
bool  ,
int  level,
const BCRec *  bc_ptr,
bool  ,
const std::unique_ptr< MultiFab > &  z_phys_nd 
)
24 {
25  {
26  /*
27  YSU PBL initially introduced by S.-Y. Hong, Y. Noh, and J. Dudhia, MWR, 2006 [HND06]
28 
29  Further Modifications from S.-Y. Hong, Q. J. R. Meteorol. Soc., 2010 [H10]
30 
31  Implementation follows WRF as of early 2024 with some simplifications
32  */
33 
34 #ifdef _OPENMP
35 #pragma omp parallel if (Gpu::notInLaunchRegion())
36 #endif
37  for ( MFIter mfi(eddyViscosity,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
38 
39  // Pull out the box we're working on, make sure it covers full domain in z-direction
40  const Box &bx = mfi.growntilebox(1);
41  const Box &dbx = geom.Domain();
42  Box sbx(bx.smallEnd(), bx.bigEnd());
43  sbx.grow(2,-1);
44  AMREX_ALWAYS_ASSERT(sbx.smallEnd(2) == dbx.smallEnd(2) && sbx.bigEnd(2) == dbx.bigEnd(2));
45 
46  // Get some data in arrays
47  const auto& cell_data = cons_in.const_array(mfi);
48  const auto& uvel = xvel.const_array(mfi);
49  const auto& vvel = yvel.const_array(mfi);
50 
51  const auto& z0_arr = most->get_z0(level)->const_array();
52  const auto& ws10av_arr = most->get_mac_avg(level,5)->const_array(mfi);
53  const auto& t10av_arr = most->get_mac_avg(level,2)->const_array(mfi);
54  const auto& t_surf_arr = most->get_t_surf(level)->const_array(mfi);
55  const auto& over_land_arr = (most->get_lmask(level)) ? most->get_lmask(level)->const_array(mfi) : Array4<int> {};
56  const Array4<Real const> z_nd_arr = z_phys_nd->array(mfi);
57  const Real most_zref = most->get_zref();
58 
59  // Require that MOST zref is 10 m so we get the wind speed at 10 m from most
60  bool invalid_zref = false;
61  if (use_terrain_fitted_coords) {
62  invalid_zref = most_zref != 10.0;
63  } else {
64  // zref gets reset to nearest cell center, so assert that zref is in the same cell as the 10m point
65  Real dz = geom.CellSize(2);
66  invalid_zref = int((most_zref - 0.5*dz)/dz) != int((10.0 - 0.5*dz)/dz);
67  }
68  if (invalid_zref) {
69  Print() << "most_zref = " << most_zref << std::endl;
70  Abort("MOST Zref must be 10m for YSU PBL scheme");
71  }
72 
73  // create flattened boxes to store PBL height
74  const GeometryData gdata = geom.data();
75  const Box xybx = PerpendicularBox<ZDir>(bx, IntVect{0,0,0});
76  FArrayBox pbl_height(xybx,1);
77  IArrayBox pbl_index(xybx,1);
78  const auto& pblh_arr = pbl_height.array();
79  const auto& pbli_arr = pbl_index.array();
80 
81  // -- Diagnose PBL height - starting out assuming non-moist --
82  // loop is only over i,j in order to find height at each x,y
83  const Real f0 = turbChoice.pbl_ysu_coriolis_freq;
84  const bool force_over_water = turbChoice.pbl_ysu_force_over_water;
85  const Real land_Ribcr = turbChoice.pbl_ysu_land_Ribcr;
86  const Real unst_Ribcr = turbChoice.pbl_ysu_unst_Ribcr;
87  ParallelFor(xybx, [=] AMREX_GPU_DEVICE (int i, int j, int) noexcept
88  {
89  // Reconstruct a surface bulk Richardson number from the surface layer model
90  // In WRF, this value is supplied to YSU by the MM5 surface layer model
91  const Real t_surf = t_surf_arr(i,j,0);
92  const Real t_layer = t10av_arr(i,j,0);
93  const Real ws_layer = ws10av_arr(i,j,0);
94  const Real Rib_layer = CONST_GRAV * most_zref / (ws_layer*ws_layer) * (t_layer - t_surf)/(t_layer);
95 
96  // For now, we only support stable boundary layers
97  if (Rib_layer < unst_Ribcr) {
98  Abort("For now, YSU PBL only supports stable conditions");
99  }
100 
101  // TODO: unstable BLs
102 
103  // PBL Height: Stable Conditions
104  Real Rib_cr;
105  bool over_land = (over_land_arr) ? over_land_arr(i,j,0) : 1;
106  if (over_land && !force_over_water) {
107  Rib_cr = land_Ribcr;
108  } else { // over water
109  // Velocity at z=10 m comes from MOST -> currently the average using whatever averaging MOST uses.
110  // TODO: Revisit this calculation with local ws10?
111  const Real z0 = z0_arr(i,j,0);
112  const Real Rossby = ws_layer/(f0*z0);
113  Rib_cr = min(0.16*std::pow(1.0e-7*Rossby,-0.18),0.3); // Note: upper bound in WRF code, but not H10 paper
114  }
115 
116  bool above_critical = false;
117  int kpbl = 0;
118  Real Rib_up = Rib_layer, Rib_dn;
119  const Real base_theta = cell_data(i,j,0,RhoTheta_comp) / cell_data(i,j,0,Rho_comp);
120  while (!above_critical and bx.contains(i,j,kpbl+1)) {
121  kpbl += 1;
122  const Real zval = use_terrain_fitted_coords ?
123  Compute_Zrel_AtCellCenter(i,j,kpbl,z_nd_arr) : gdata.ProbLo(2) + (kpbl + 0.5)*gdata.CellSize(2);
124  const Real ws2_level = 0.25*( (uvel(i,j,kpbl)+uvel(i+1,j ,kpbl))*(uvel(i,j,kpbl)+uvel(i+1,j ,kpbl))
125  + (vvel(i,j,kpbl)+vvel(i ,j+1,kpbl))*(vvel(i,j,kpbl)+vvel(i ,j+1,kpbl)) );
126  const Real theta = cell_data(i,j,kpbl,RhoTheta_comp) / cell_data(i,j,kpbl,Rho_comp);
127  Rib_dn = Rib_up;
128  Rib_up = (theta-base_theta)/base_theta * CONST_GRAV * zval / ws2_level;
129  above_critical = Rib_up >= Rib_cr;
130  }
131 
132  Real interp_fact;
133  if (Rib_dn >= Rib_cr) {
134  interp_fact = 0.0;
135  } else if (Rib_up <= Rib_cr)
136  interp_fact = 1.0;
137  else {
138  interp_fact = (Rib_cr - Rib_dn) / (Rib_up - Rib_dn);
139  }
140 
141  const Real zval_up = use_terrain_fitted_coords ?
142  Compute_Zrel_AtCellCenter(i,j,kpbl,z_nd_arr) : gdata.ProbLo(2) + (kpbl + 0.5)*gdata.CellSize(2);
143  const Real zval_dn = use_terrain_fitted_coords ?
144  Compute_Zrel_AtCellCenter(i,j,kpbl-1,z_nd_arr) : gdata.ProbLo(2) + (kpbl-1 + 0.5)*gdata.CellSize(2);
145  pblh_arr(i,j,0) = zval_dn + interp_fact*(zval_up-zval_dn);
146 
147  const Real zval_0 = use_terrain_fitted_coords ?
148  Compute_Zrel_AtCellCenter(i,j,0,z_nd_arr) : gdata.ProbLo(2) + (0.5)*gdata.CellSize(2);
149  const Real zval_1 = use_terrain_fitted_coords ?
150  Compute_Zrel_AtCellCenter(i,j,1,z_nd_arr) : gdata.ProbLo(2) + (1.5)*gdata.CellSize(2);
151  if (pblh_arr(i,j,0) < 0.5*(zval_0+zval_1) ) {
152  kpbl = 0;
153  }
154  pbli_arr(i,j,0) = kpbl;
155  });
156 
157  // -- Compute nonlocal/countergradient mixing parameters --
158  // Not included for stable so nothing to do until unstable treatment is added
159 
160  // -- Compute entrainment parameters --
161  // 0 for stable so nothing to do?
162 
163  // -- Compute diffusion coefficients --
164 
165  const auto& u_star_arr = most->get_u_star(level)->const_array(mfi);
166  const auto& l_obuk_arr = most->get_olen(level)->const_array(mfi);
167  const Array4<Real > &K_turb = eddyViscosity.array(mfi);
168 
169  // Dirichlet flags to switch derivative stencil
170  bool c_ext_dir_on_zlo = ( (bc_ptr[BCVars::cons_bc].lo(2) == ERFBCType::ext_dir) );
171  bool c_ext_dir_on_zhi = ( (bc_ptr[BCVars::cons_bc].lo(5) == ERFBCType::ext_dir) );
172  bool u_ext_dir_on_zlo = ( (bc_ptr[BCVars::xvel_bc].lo(2) == ERFBCType::ext_dir) );
173  bool u_ext_dir_on_zhi = ( (bc_ptr[BCVars::xvel_bc].lo(5) == ERFBCType::ext_dir) );
174  bool v_ext_dir_on_zlo = ( (bc_ptr[BCVars::yvel_bc].lo(2) == ERFBCType::ext_dir) );
175  bool v_ext_dir_on_zhi = ( (bc_ptr[BCVars::yvel_bc].lo(5) == ERFBCType::ext_dir) );
176 
177  const auto& dxInv = geom.InvCellSizeArray();
178  const Real dz_inv = geom.InvCellSize(2);
179  const int izmin = geom.Domain().smallEnd(2);
180  const int izmax = geom.Domain().bigEnd(2);
181 
182  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
183  {
184  const Real zval = use_terrain_fitted_coords ?
185  Compute_Zrel_AtCellCenter(i,j,k,z_nd_arr) : gdata.ProbLo(2) + (k + 0.5)*gdata.CellSize(2);
186  const Real rho = cell_data(i,j,k,Rho_comp);
187  const Real met_h_zeta = use_terrain_fitted_coords ? Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_nd_arr) : 1.0;
188  const Real dz_terrain = met_h_zeta/dz_inv;
189  if (k < pbli_arr(i,j,0)) {
190  // -- Compute diffusion coefficients within PBL
191  constexpr Real zfacmin = 1e-8; // value from WRF
192  constexpr Real phifac = 8.0; // value from H10 and WRF
193  constexpr Real wstar3 = 0.0; // only nonzero for unstable
194  constexpr Real pfac = 2.0; // profile exponent
195  const Real zfac = std::min(std::max(1 - zval / pblh_arr(i,j,0), zfacmin ), 1.0);
196  // Not including YSU top down PBL term (not in H10, added to WRF later)
197  const Real ust3 = u_star_arr(i,j,0) * u_star_arr(i,j,0) * u_star_arr(i,j,0);
198  Real wscalek = ust3 + phifac * KAPPA * wstar3 * (1.0 - zfac);
199  wscalek = std::pow(wscalek, 1.0/3.0);
200  // stable only
201  const Real phi_term = 1 + 5 * zval / l_obuk_arr(i,j,0); // phi_term appears in WRF but not papers
202  wscalek = std::max(u_star_arr(i,j,0) / phi_term, 0.001); // 0.001 limit appears in WRF but not papers
203  K_turb(i,j,k,EddyDiff::Mom_v) = rho * wscalek * KAPPA * zval * std::pow(zfac, pfac);
204  K_turb(i,j,k,EddyDiff::Theta_v) = K_turb(i,j,k,EddyDiff::Mom_v);
205  } else {
206  // -- Compute coefficients in free stream above PBL
207  constexpr Real lam0 = 30.0;
208  constexpr Real min_richardson = -100.0;
209  constexpr Real prandtl_max = 4.0;
210  Real dthetadz, dudz, dvdz;
211  const int RhoQv_comp = -1;
212  const int RhoQc_comp = -1;
213  const int RhoQr_comp = -1;
215  uvel, vvel, cell_data, izmin, izmax, 1.0/dz_terrain,
216  c_ext_dir_on_zlo, c_ext_dir_on_zhi,
217  u_ext_dir_on_zlo, u_ext_dir_on_zhi,
218  v_ext_dir_on_zlo, v_ext_dir_on_zhi,
219  dthetadz, dudz, dvdz, RhoQv_comp, RhoQc_comp, RhoQr_comp);
220  const Real shear_squared = dudz*dudz + dvdz*dvdz + 1.0e-9; // 1.0e-9 from WRF to avoid divide by zero
221  const Real theta = cell_data(i,j,k,RhoTheta_comp) / cell_data(i,j,k,Rho_comp);
222  Real richardson = CONST_GRAV / theta * dthetadz / shear_squared;
223  const Real lambdadz = std::min(std::max(0.1*dz_terrain , lam0), 300.0); // in WRF, H10 paper just says use lam0
224  const Real lengthscale = lambdadz * KAPPA * zval / (lambdadz + KAPPA * zval);
225  const Real turbfact = lengthscale * lengthscale * std::sqrt(shear_squared);
226 
227  if (richardson < 0) {
228  richardson = max(richardson, min_richardson);
229  Real sqrt_richardson = std::sqrt(-richardson);
230  K_turb(i,j,k,EddyDiff::Mom_v) = rho * turbfact * (1.0 - 8.0 * richardson / (1.0 + 1.746 * sqrt_richardson));
231  K_turb(i,j,k,EddyDiff::Theta_v) = rho * turbfact * (1.0 - 8.0 * richardson / (1.0 + 1.286 * sqrt_richardson));
232  } else {
233  const Real oneplus5ri = 1.0 + 5.0 * richardson;
234  K_turb(i,j,k,EddyDiff::Theta_v) = rho * turbfact / (oneplus5ri * oneplus5ri);
235  const Real prandtl = std::min(1.0+2.1*richardson, prandtl_max); // limit from WRF
236  K_turb(i,j,k,EddyDiff::Mom_v) = K_turb(i,j,k,EddyDiff::Theta_v) * prandtl;
237  }
238  }
239 
240  // limit both diffusion coefficients - from WRF, not documented in papers
241  constexpr Real ckz = 0.001;
242  constexpr Real Kmax = 1000.0;
243  const Real rhoKmin = ckz * dz_terrain * rho;
244  const Real rhoKmax = rho * Kmax;
245  K_turb(i,j,k,EddyDiff::Mom_v) = std::max(std::min(K_turb(i,j,k,EddyDiff::Mom_v) ,rhoKmax), rhoKmin);
246  K_turb(i,j,k,EddyDiff::Theta_v) = std::max(std::min(K_turb(i,j,k,EddyDiff::Theta_v) ,rhoKmax), rhoKmin);
247  K_turb(i,j,k,EddyDiff::Turb_lengthscale) = pblh_arr(i,j,0);
248  });
249 
250  // HACK set bottom ghost cell to 1st cell
251  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
252  {
253  if (k==-1) {
254  K_turb(i,j,k,EddyDiff::Mom_v) = K_turb(i,j,0,EddyDiff::Mom_v);
255  K_turb(i,j,k,EddyDiff::Theta_v) = K_turb(i,j,0,EddyDiff::Theta_v);
256  }
257  });
258  }
259  }
260 }
constexpr amrex::Real KAPPA
Definition: ERF_Constants.H:20
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:21
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void ComputeVerticalDerivativesPBL(int i, int j, int k, const amrex::Array4< const amrex::Real > &uvel, const amrex::Array4< const amrex::Real > &vvel, const amrex::Array4< const amrex::Real > &cell_data, const int izmin, const int izmax, const amrex::Real &dz_inv, const bool c_ext_dir_on_zlo, const bool c_ext_dir_on_zhi, const bool u_ext_dir_on_zlo, const bool u_ext_dir_on_zhi, const bool v_ext_dir_on_zlo, const bool v_ext_dir_on_zhi, amrex::Real &dthetadz, amrex::Real &dudz, amrex::Real &dvdz, const int RhoQv_comp, const int RhoQc_comp, const int RhoQr_comp, const bool use_most=true)
Definition: ERF_PBLModels.H:119
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:46
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Zrel_AtCellCenter(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:359
@ yvel_bc
Definition: ERF_IndexDefines.H:88
@ cons_bc
Definition: ERF_IndexDefines.H:76
@ xvel_bc
Definition: ERF_IndexDefines.H:87
@ ext_dir
Definition: ERF_IndexDefines.H:191
@ Theta_v
Definition: ERF_IndexDefines.H:168
@ Turb_lengthscale
Definition: ERF_IndexDefines.H:172
@ Mom_v
Definition: ERF_IndexDefines.H:167
@ theta
Definition: ERF_MM5.H:20
@ rho
Definition: ERF_Kessler.H:22
@ xvel
Definition: ERF_IndexDefines.H:141
@ yvel
Definition: ERF_IndexDefines.H:142
amrex::Real pbl_ysu_land_Ribcr
Definition: ERF_TurbStruct.H:250
amrex::Real pbl_ysu_coriolis_freq
Definition: ERF_TurbStruct.H:247
bool pbl_ysu_force_over_water
Definition: ERF_TurbStruct.H:249
amrex::Real pbl_ysu_unst_Ribcr
Definition: ERF_TurbStruct.H:251

Referenced by ComputeTurbulentViscosity().

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