ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_PBLModels.H
Go to the documentation of this file.
1 #ifndef ERF_PBLMODELS_H_
2 #define ERF_PBLMODELS_H_
3 
4 #include <ERF_Thetav.H>
5 
6 /**
7  * Compute eddy diffusivities of momentum (eddy viscosity) and heat using the
8  * Mellor-Yamada-Nakanishi-Niino Level 2.5 PBL scheme
9  *
10  * @param[in] xvel Velocity in x-dir
11  * @param[in] yvel Velocity in y-dir
12  * @param[in] cons_in Cell center conserved quantities
13  * @param[out] eddyViscosity Holds turbulent viscosity
14  * @param[in] geom Problem geometry
15  * @param[in] turbChoice Container with turbulence parameters
16  * @param[in] most Pointer to Monin-Obukhov class if instantiated
17  * @param[in] use_moisture If we have microphysics enabled
18  * @param[in] level Current level
19  * @param[in] bc_ptr Pointer to array with boundary condition info
20  * @param[in] vert_only Only compute vertical eddy diffusivities
21  * @param[in] z_phys_nd Physical location of grid nodes, if terrain (or grid stretching) is enabled
22  */
23 void
24 ComputeDiffusivityMYNN25 (const amrex::MultiFab& xvel,
25  const amrex::MultiFab& yvel,
26  const amrex::MultiFab& cons_in,
27  amrex::MultiFab& eddyViscosity,
28  const amrex::Geometry& geom,
29  const TurbChoice& turbChoice,
30  std::unique_ptr<ABLMost>& most,
31  bool use_moisture,
32  int level,
33  const amrex::BCRec* bc_ptr,
34  bool /*vert_only*/,
35  const std::unique_ptr<amrex::MultiFab>& z_phys_nd,
36  const int RhoQv_comp,
37  const int RhoQc_comp,
38  const int RhoQr_comp);
39 
40 /**
41  * Compute eddy diffusivities of momentum (eddy viscosity) and heat using the
42  * Yonsei University PBL scheme
43  *
44  * @param[in] xvel Velocity in x-dir
45  * @param[in] yvel Velocity in y-dir
46  * @param[in] cons_in Cell center conserved quantities
47  * @param[out] eddyViscosity Holds turbulent viscosity
48  * @param[in] geom Problem geometry
49  * @param[in] turbChoice Container with turbulence parameters
50  * @param[in] most Pointer to Monin-Obukhov class if instantiated
51  * @param[in] use_moisture If we have microphysics enabled
52  * @param[in] level Current level
53  * @param[in] bc_ptr Pointer to array with boundary condition info
54  * @param[in] vert_only Only compute vertical eddy diffusivities
55  * @param[in] z_phys_nd Physical location of grid nodes, if terrain (or grid stretching) is enabled
56  */
57 void
58 ComputeDiffusivityYSU (const amrex::MultiFab& xvel,
59  const amrex::MultiFab& yvel,
60  const amrex::MultiFab& cons_in,
61  amrex::MultiFab& eddyViscosity,
62  const amrex::Geometry& geom,
63  const TurbChoice& turbChoice,
64  std::unique_ptr<ABLMost>& most,
65  bool use_moisture,
66  int level,
67  const amrex::BCRec* bc_ptr,
68  bool /*vert_only*/,
69  const std::unique_ptr<amrex::MultiFab>& z_phys_nd);
70 
71 
72 /**
73  * Function for computing vertical derivatives for use in PBL model
74  *
75  * @param[in] u velocity in x-dir
76  * @param[in] v velocity in y-dir
77  * @param[in] cell_data conserved cell center vars
78  */
79 AMREX_GPU_DEVICE
80 AMREX_FORCE_INLINE
81 void
82 ComputeVerticalDerivativesPBL (int i, int j, int k,
83  const amrex::Array4<const amrex::Real>& uvel,
84  const amrex::Array4<const amrex::Real>& vvel,
85  const amrex::Array4<const amrex::Real>& cell_data,
86  const int izmin,
87  const int izmax,
88  const amrex::Real& dz_inv,
89  const bool c_ext_dir_on_zlo,
90  const bool c_ext_dir_on_zhi,
91  const bool u_ext_dir_on_zlo,
92  const bool u_ext_dir_on_zhi,
93  const bool v_ext_dir_on_zlo,
94  const bool v_ext_dir_on_zhi,
95  amrex::Real& dthetadz,
96  amrex::Real& dudz,
97  amrex::Real& dvdz,
98  const int RhoQv_comp,
99  const int RhoQc_comp,
100  const int RhoQr_comp,
101  const bool use_most=true)
102 {
103  if ( k==izmax && c_ext_dir_on_zhi ) {
104  dthetadz = (1.0/3.0)*( -Thetav(i,j,k-1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
105  - 3.0 * Thetav(i,j,k ,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
106  + 4.0 * Thetav(i,j,k+1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp) )*dz_inv;
107  } else if ( k==izmin && c_ext_dir_on_zlo ) {
108  dthetadz = (1.0/3.0)*( Thetav(i,j,k+1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
109  + 3.0 * Thetav(i,j,k ,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
110  - 4.0 * Thetav(i,j,k-1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp) )*dz_inv;
111  } else if ( k==izmin && use_most ) {
112  dthetadz = ( Thetav(i,j,k+1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
113  - Thetav(i,j,k ,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp) )*dz_inv;
114  } else {
115  dthetadz = 0.5*( Thetav(i,j,k+1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
116  - Thetav(i,j,k-1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp) )*dz_inv;
117  }
118 
119  if ( k==izmax && u_ext_dir_on_zhi ) {
120  dudz = (1.0/6.0)*( (-uvel(i ,j,k-1) - 3.0 * uvel(i ,j,k ) + 4.0 * uvel(i ,j,k+1))
121  + (-uvel(i+1,j,k-1) - 3.0 * uvel(i+1,j,k ) + 4.0 * uvel(i+1,j,k+1)) )*dz_inv;
122  } else if ( k==izmin && u_ext_dir_on_zlo ) {
123  dudz = (1.0/6.0)*( (uvel(i ,j,k+1) + 3.0 * uvel(i ,j,k ) - 4.0 * uvel(i ,j,k-1))
124  + (uvel(i+1,j,k+1) + 3.0 * uvel(i+1,j,k ) - 4.0 * uvel(i+1,j,k-1)) )*dz_inv;
125  } else if ( k==izmin && use_most ) {
126  dudz = 0.50*( uvel(i,j,k+1) - uvel(i,j,k ) + uvel(i+1,j,k+1) - uvel(i+1,j,k ) )*dz_inv;
127  } else {
128  dudz = 0.25*( uvel(i,j,k+1) - uvel(i,j,k-1) + uvel(i+1,j,k+1) - uvel(i+1,j,k-1) )*dz_inv;
129  }
130 
131  if ( k==izmax && v_ext_dir_on_zhi ) {
132  dvdz = (1.0/6.0)*( (-vvel(i,j ,k-1) - 3.0 * vvel(i,j ,k ) + 4.0 * vvel(i,j ,k+1))
133  + (-vvel(i,j+1,k-1) - 3.0 * vvel(i,j+1,k ) + 4.0 * vvel(i,j+1,k+1)) )*dz_inv;
134  } else if ( k==izmin && v_ext_dir_on_zlo ) {
135  dvdz = (1.0/6.0)*( (vvel(i,j ,k+1) + 3.0 * vvel(i,j ,k ) - 4.0 * vvel(i,j ,k-1))
136  + (vvel(i,j+1,k+1) + 3.0 * vvel(i,j+1,k ) - 4.0 * vvel(i,j+1,k-1)) )*dz_inv;
137  } else if ( k==izmin && use_most ) {
138  dvdz = 0.50*( vvel(i,j,k+1) - vvel(i,j,k ) + vvel(i,j+1,k+1) - vvel(i,j+1,k ) )*dz_inv;
139  } else {
140  dvdz = 0.25*( vvel(i,j,k+1) - vvel(i,j,k-1) + vvel(i,j+1,k+1) - vvel(i,j+1,k-1) )*dz_inv;
141  }
142 }
143 
144 /**
145  * Function for computing the QKE source terms (NN09, Eqn. 5).
146  *
147  * @param[in] u velocity in x-dir
148  * @param[in] v velocity in y-dir
149  * @param[in] cell_data conserved cell center vars
150  * @param[in] cell_prim primitive cell center vars
151  * @param[in] K_turb turbulent viscosity
152  * @param[in] cellSizeInv inverse cell size array
153  * @param[in] domain box of the whole domain
154  * @param[in] pbl_mynn_B1_l a parameter
155  * @param[in] theta_mean average theta
156  */
157 AMREX_GPU_DEVICE
158 AMREX_FORCE_INLINE
159 amrex::Real
160 ComputeQKESourceTerms (int i, int j, int k,
161  const amrex::Array4<const amrex::Real>& uvel,
162  const amrex::Array4<const amrex::Real>& vvel,
163  const amrex::Array4<const amrex::Real>& cell_data,
164  const amrex::Array4<const amrex::Real>& cell_prim,
165  const amrex::Array4<const amrex::Real>& K_turb,
166  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
167  const amrex::Box& domain,
168  amrex::Real pbl_mynn_B1_l,
169  const amrex::Real theta_mean,
170  const int RhoQv_comp,
171  const int RhoQc_comp,
172  const int RhoQr_comp,
173  bool c_ext_dir_on_zlo,
174  bool c_ext_dir_on_zhi,
175  bool u_ext_dir_on_zlo,
176  bool u_ext_dir_on_zhi,
177  bool v_ext_dir_on_zlo,
178  bool v_ext_dir_on_zhi,
179  const bool use_most=false,
180  const amrex::Real met_h_zeta=1.0)
181 {
182  // Compute some relevant derivatives
183  amrex::Real dthetadz, dudz, dvdz;
184  amrex::Real source_term = 0.0;
185 
186  amrex::Real dz_inv = cellSizeInv[2];
187  int izmin = domain.smallEnd(2);
188  int izmax = domain.bigEnd(2);
189 
190  // NOTE: With MOST, the ghost cells are filled AFTER k_turb is computed
191  // so that the non-explicit pathway works. Therefore, at this
192  // point we DO have valid ghost cells from MOST. We are passing
193  // the MOST flag to use one-sided diffs here to be consistent with
194  // the explicit pathway.
195 
197  uvel, vvel, cell_data, izmin, izmax, dz_inv/met_h_zeta,
198  c_ext_dir_on_zlo, c_ext_dir_on_zhi,
199  u_ext_dir_on_zlo, u_ext_dir_on_zhi,
200  v_ext_dir_on_zlo, v_ext_dir_on_zhi,
201  dthetadz, dudz, dvdz,
202  RhoQv_comp, RhoQc_comp, RhoQr_comp, use_most);
203 
204  // Note: Transport terms due to turbulence and pressure are included when
205  // DiffusionSrcForState_* is called from ERF_slow_rhs_post.
206 
207  // Shear Production
208  source_term += K_turb(i,j,k,EddyDiff::Mom_v) * (dudz*dudz + dvdz*dvdz);
209 
210  // Buoyancy Production
211  source_term -= (CONST_GRAV/theta_mean)*K_turb(i,j,k,EddyDiff::Theta_v)*dthetadz;
212 
213  // Dissipation
214  amrex::Real qke = 2.0 * cell_prim(i,j,k,PrimKE_comp);
215  if (std::abs(qke) > 0.0) {
216  source_term -= cell_data(i,j,k,Rho_comp) * std::pow(qke,1.5) /
217  (pbl_mynn_B1_l * K_turb(i,j,k,EddyDiff::Turb_lengthscale));
218  }
219 
220  return source_term;
221 }
222 #endif
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:21
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define PrimKE_comp
Definition: ERF_IndexDefines.H:51
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeQKESourceTerms(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 amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< const amrex::Real > &K_turb, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Box &domain, amrex::Real pbl_mynn_B1_l, const amrex::Real theta_mean, const int RhoQv_comp, const int RhoQc_comp, const int RhoQr_comp, bool c_ext_dir_on_zlo, bool c_ext_dir_on_zhi, bool u_ext_dir_on_zlo, bool u_ext_dir_on_zhi, bool v_ext_dir_on_zlo, bool v_ext_dir_on_zhi, const bool use_most=false, const amrex::Real met_h_zeta=1.0)
Definition: ERF_PBLModels.H:160
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:82
void ComputeDiffusivityYSU(const amrex::MultiFab &xvel, const amrex::MultiFab &yvel, const amrex::MultiFab &cons_in, amrex::MultiFab &eddyViscosity, const amrex::Geometry &geom, const TurbChoice &turbChoice, std::unique_ptr< ABLMost > &most, bool use_moisture, int level, const amrex::BCRec *bc_ptr, bool, const std::unique_ptr< amrex::MultiFab > &z_phys_nd)
void ComputeDiffusivityMYNN25(const amrex::MultiFab &xvel, const amrex::MultiFab &yvel, const amrex::MultiFab &cons_in, amrex::MultiFab &eddyViscosity, const amrex::Geometry &geom, const TurbChoice &turbChoice, std::unique_ptr< ABLMost > &most, bool use_moisture, int level, const amrex::BCRec *bc_ptr, bool, const std::unique_ptr< amrex::MultiFab > &z_phys_nd, const int RhoQv_comp, const int RhoQc_comp, const int RhoQr_comp)
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Thetav(int i, int j, int k, const amrex::Array4< const amrex::Real > &cell_data, const int RhoQv_comp, const int RhoQc_comp, const int RhoQr_comp)
Definition: ERF_Thetav.H:10
@ Theta_v
Definition: ERF_IndexDefines.H:157
@ Turb_lengthscale
Definition: ERF_IndexDefines.H:161
@ Mom_v
Definition: ERF_IndexDefines.H:156
@ xvel
Definition: ERF_IndexDefines.H:130
@ yvel
Definition: ERF_IndexDefines.H:131
Definition: ERF_TurbStruct.H:31