1 #ifndef ERF_PBLMODELS_H_
2 #define ERF_PBLMODELS_H_
25 const amrex::MultiFab&
yvel,
26 const amrex::MultiFab& cons_in,
27 amrex::MultiFab& eddyViscosity,
28 const amrex::Geometry& geom,
30 std::unique_ptr<ABLMost>& most,
31 bool use_terrain_fitted_coords,
34 const amrex::BCRec* bc_ptr,
36 const std::unique_ptr<amrex::MultiFab>& z_phys_nd,
39 const int RhoQr_comp);
60 const amrex::MultiFab&
yvel,
61 const amrex::MultiFab& cons_in,
62 amrex::MultiFab& eddyViscosity,
63 const amrex::Geometry& geom,
65 std::unique_ptr<ABLMost>& most,
66 bool use_terrain_fitted_coords,
69 const amrex::BCRec* bc_ptr,
71 const std::unique_ptr<amrex::MultiFab>& z_phys_nd,
74 const int RhoQr_comp);
95 const amrex::MultiFab&
yvel,
96 const amrex::MultiFab& cons_in,
97 amrex::MultiFab& eddyViscosity,
98 const amrex::Geometry& geom,
100 std::unique_ptr<ABLMost>& most,
101 bool use_terrain_fitted_coords,
104 const amrex::BCRec* bc_ptr,
106 const std::unique_ptr<amrex::MultiFab>& z_phys_nd);
120 const amrex::Array4<const amrex::Real>& uvel,
121 const amrex::Array4<const amrex::Real>& vvel,
122 const amrex::Array4<const amrex::Real>& cell_data,
125 const amrex::Real& dz_inv,
126 const bool c_ext_dir_on_zlo,
127 const bool c_ext_dir_on_zhi,
128 const bool u_ext_dir_on_zlo,
129 const bool u_ext_dir_on_zhi,
130 const bool v_ext_dir_on_zlo,
131 const bool v_ext_dir_on_zhi,
132 amrex::Real& dthetadz,
135 const int RhoQv_comp,
136 const int RhoQc_comp,
137 const int RhoQr_comp,
138 const bool use_most=
true)
140 if ( k==izmax && c_ext_dir_on_zhi ) {
141 dthetadz = (1.0/3.0)*( -
Thetav(i,j,k-1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
142 - 3.0 *
Thetav(i,j,k ,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
143 + 4.0 *
Thetav(i,j,k+1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp) )*dz_inv;
144 }
else if ( k==izmin && c_ext_dir_on_zlo ) {
145 dthetadz = (1.0/3.0)*(
Thetav(i,j,k+1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
146 + 3.0 *
Thetav(i,j,k ,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
147 - 4.0 *
Thetav(i,j,k-1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp) )*dz_inv;
148 }
else if ( k==izmin && use_most ) {
149 dthetadz = (
Thetav(i,j,k+1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
150 -
Thetav(i,j,k ,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp) )*dz_inv;
152 dthetadz = 0.5*(
Thetav(i,j,k+1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp)
153 -
Thetav(i,j,k-1,cell_data,RhoQv_comp,RhoQc_comp,RhoQr_comp) )*dz_inv;
156 if ( k==izmax && u_ext_dir_on_zhi ) {
157 dudz = (1.0/6.0)*( (-uvel(i ,j,k-1) - 3.0 * uvel(i ,j,k ) + 4.0 * uvel(i ,j,k+1))
158 + (-uvel(i+1,j,k-1) - 3.0 * uvel(i+1,j,k ) + 4.0 * uvel(i+1,j,k+1)) )*dz_inv;
159 }
else if ( k==izmin && u_ext_dir_on_zlo ) {
160 dudz = (1.0/6.0)*( (uvel(i ,j,k+1) + 3.0 * uvel(i ,j,k ) - 4.0 * uvel(i ,j,k-1))
161 + (uvel(i+1,j,k+1) + 3.0 * uvel(i+1,j,k ) - 4.0 * uvel(i+1,j,k-1)) )*dz_inv;
162 }
else if ( k==izmin && use_most ) {
163 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;
165 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;
168 if ( k==izmax && v_ext_dir_on_zhi ) {
169 dvdz = (1.0/6.0)*( (-vvel(i,j ,k-1) - 3.0 * vvel(i,j ,k ) + 4.0 * vvel(i,j ,k+1))
170 + (-vvel(i,j+1,k-1) - 3.0 * vvel(i,j+1,k ) + 4.0 * vvel(i,j+1,k+1)) )*dz_inv;
171 }
else if ( k==izmin && v_ext_dir_on_zlo ) {
172 dvdz = (1.0/6.0)*( (vvel(i,j ,k+1) + 3.0 * vvel(i,j ,k ) - 4.0 * vvel(i,j ,k-1))
173 + (vvel(i,j+1,k+1) + 3.0 * vvel(i,j+1,k ) - 4.0 * vvel(i,j+1,k-1)) )*dz_inv;
174 }
else if ( k==izmin && use_most ) {
175 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;
177 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;
198 const amrex::Array4<const amrex::Real>& uvel,
199 const amrex::Array4<const amrex::Real>& vvel,
200 const amrex::Array4<const amrex::Real>& cell_data,
201 const amrex::Array4<const amrex::Real>& cell_prim,
202 const amrex::Array4<const amrex::Real>& K_turb,
203 const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
204 const amrex::Box& domain,
205 amrex::Real pbl_mynn_B1_l,
206 const amrex::Real theta_mean,
207 const int RhoQv_comp,
208 const int RhoQc_comp,
209 const int RhoQr_comp,
210 bool c_ext_dir_on_zlo,
211 bool c_ext_dir_on_zhi,
212 bool u_ext_dir_on_zlo,
213 bool u_ext_dir_on_zhi,
214 bool v_ext_dir_on_zlo,
215 bool v_ext_dir_on_zhi,
216 const bool use_most=
false,
217 const amrex::Real met_h_zeta=1.0)
220 amrex::Real dthetadz, dudz, dvdz;
221 amrex::Real source_term = 0.0;
223 amrex::Real dz_inv = cellSizeInv[2];
224 int izmin = domain.smallEnd(2);
225 int izmax = domain.bigEnd(2);
234 uvel, vvel, cell_data, izmin, izmax, dz_inv/met_h_zeta,
235 c_ext_dir_on_zlo, c_ext_dir_on_zhi,
236 u_ext_dir_on_zlo, u_ext_dir_on_zhi,
237 v_ext_dir_on_zlo, v_ext_dir_on_zhi,
238 dthetadz, dudz, dvdz,
239 RhoQv_comp, RhoQc_comp, RhoQr_comp, use_most);
245 source_term += K_turb(i,j,k,
EddyDiff::Mom_v) * (dudz*dudz + dvdz*dvdz);
251 amrex::Real qke = 2.0 * cell_prim(i,j,k,
PrimKE_comp);
252 if (std::abs(qke) > 0.0) {
253 source_term -= cell_data(i,j,k,
Rho_comp) * std::pow(qke,1.5) /
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
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_terrain_fitted_coords, bool use_moisture, int level, const amrex::BCRec *bc_ptr, bool, const std::unique_ptr< amrex::MultiFab > &z_phys_nd)
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:197
void ComputeDiffusivityMYNNEDMF(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_terrain_fitted_coords, 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)
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_terrain_fitted_coords, 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 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_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:168
@ Turb_lengthscale
Definition: ERF_IndexDefines.H:172
@ Mom_v
Definition: ERF_IndexDefines.H:167
@ xvel
Definition: ERF_IndexDefines.H:141
@ yvel
Definition: ERF_IndexDefines.H:142
Definition: ERF_TurbStruct.H:31