ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
TurbChoice Struct Reference

#include <ERF_TurbStruct.H>

Collaboration diagram for TurbChoice:

Public Member Functions

void init_params (int lev, int max_level, std::string pp_prefix)
 
void check_params (amrex::GpuArray< ERF_BC, AMREX_SPACEDIM *2 > &phys_bc_type)
 
void display (int lev)
 

Public Attributes

LESType les_type
 
amrex::Real Pr_t = one / three
 
amrex::Real Pr_t_inv = three
 
amrex::Real Sc_t = one
 
amrex::Real Sc_t_inv = one
 
amrex::Real Cs = zero
 
bool smag2d = false
 
amrex::Real Ce = amrex::Real(0.93)
 
amrex::Real Ce_wall = zero
 
amrex::Real Ck = amrex::Real(0.1)
 
amrex::Real Cmu0 = amrex::Real(0.5562)
 
amrex::Real Cb = amrex::Real(0.35)
 
amrex::Real Rt_crit = -one
 
amrex::Real Rt_min = -three
 
amrex::Real l_g_max = amrex::Real(30.0)
 
amrex::Real sigma_k = amrex::Real(0.5)
 
amrex::Real theta_ref = zero
 
StratType strat_type = StratType::theta
 
bool mix_isotropic = true
 
bool use_Ri_correction = true
 
amrex::Real Ri_crit = fourth
 
RANSType rans_type
 
bool dirichlet_k = false
 
PBLType pbl_type
 
MYNNLevel25 pbl_mynn
 
MYNNLevel2 pbl_mynn_level2
 
bool use_kturb = false
 
bool use_keqn
 
bool use_pbl_tke
 
bool use_tke = false
 
bool init_tke_from_ustar = false
 
amrex::Real pbl_ysu_coriolis_freq
 
bool pbl_ysu_use_consistent_coriolis
 
bool pbl_ysu_force_over_water
 
amrex::Real pbl_ysu_land_Ribcr
 
amrex::Real pbl_ysu_unst_Ribcr
 
amrex::Real pbl_mrf_coriolis_freq = amrex::Real(1.0e-4)
 
amrex::Real pbl_mrf_Ribcr = amrex::Real(0.5)
 
amrex::Real pbl_mrf_const_b = amrex::Real(7.8)
 
amrex::Real pbl_mrf_sf = amrex::Real(0.1)
 
bool mrf_moistvars = false
 
bool advect_tke = true
 
bool diffuse_tke_3D = true
 

Detailed Description

Container holding quantities related to turbulence parametrizations

Member Function Documentation

◆ check_params()

void TurbChoice::check_params ( amrex::GpuArray< ERF_BC, AMREX_SPACEDIM *2 > &  phys_bc_type)
inline
255  {
256  // BC compatibility
257  if ( ( (pbl_type == PBLType::MYNN25) ||
258  (pbl_type == PBLType::MYNNEDMF) ||
259  (pbl_type == PBLType::YSU) ||
260  (pbl_type == PBLType::MRF)
261  ) &&
262  phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::surface_layer ) {
263  amrex::Abort("MYNN2.5/MYNNEDMF/YSU/MRF PBL Model requires MOST at lower boundary");
264  }
265  if ( (les_type == LESType::Deardorff) && (Ce_wall > 0) &&
266  (phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::surface_layer) &&
267  (phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::slip_wall) &&
268  (phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::no_slip_wall) )
269  {
270  amrex::Warning("Deardorff LES assumes wall at zlo when applying Ce_wall");
271  }
272  }
@ no_slip_wall
@ surface_layer
PBLType pbl_type
Definition: ERF_TurbStruct.H:439
amrex::Real Ce_wall
Definition: ERF_TurbStruct.H:408
LESType les_type
Definition: ERF_TurbStruct.H:392

◆ display()

void TurbChoice::display ( int  lev)
inline
275  {
276  amrex::Print() << "Turbulence Settings at level " << lev << std::endl;
277 
278  if (
279  les_type == LESType::None && rans_type == RANSType::None &&
280  pbl_type == PBLType::None) {
281  amrex::Print() << " Using DNS model at level " << lev << std::endl;
282  } else if (les_type == LESType::Smagorinsky) {
283  if (smag2d) {
284  amrex::Print() << " Using 2D Smagorinsky LES model at level " << lev << std::endl;
285  } else {
286  amrex::Print() << " Using Smagorinsky LES model at level " << lev << std::endl;
287  }
288  if (use_Ri_correction) {
289  amrex::Print() << " Smagorinsky uses Richardson number correction with Ri_crit = "
290  << Ri_crit << std::endl;
291  }
292  } else if (les_type == LESType::Deardorff) {
293  amrex::Print() << " Using Deardorff LES model at level " << lev << std::endl;
294  } else if (rans_type == RANSType::kEqn) {
295  amrex::Print()
296  << " Using Axell & Liungman one-equation RANS k model at level " << lev
297  << std::endl;
298  } else if (pbl_type == PBLType::MYJ) {
299  amrex::Print() << " Using MYJ PBL model at level " << lev << std::endl;
300  } else if (pbl_type == PBLType::MYNN25) {
301  amrex::Print() << " Using MYNN2.5 PBL model at level " << lev << std::endl;
302  } else if (pbl_type == PBLType::MYNNEDMF) {
303  amrex::Print() << " Using MYNNEDMF PBL model at level " << lev << std::endl;
304  } else if (pbl_type == PBLType::YSU) {
305  amrex::Print() << " Using YSU PBL model at level " << lev << std::endl;
306  } else if (pbl_type == PBLType::MRF) {
307  amrex::Print() << " Using MRF PBL model at level " << lev << std::endl;
308  } else {
309  amrex::Error("Unknown turbulence model");
310  }
311 
312  if (les_type != LESType::None) {
313  if (les_type == LESType::Smagorinsky) {
314  amrex::Print() << " Cs : " << Cs << std::endl;
315  }
316  if (les_type == LESType::Deardorff) {
317  amrex::Print() << " Ce : " << Ce << std::endl;
318  amrex::Print() << " Ce at wall : " << Ce_wall << std::endl;
319  amrex::Print() << " Ck : " << Ck << std::endl;
320  amrex::Print() << " sigma_k : " << sigma_k << std::endl;
321 
322  // Sullivan et al 1994, Eqn 14
323  amrex::Real Cs_equiv = std::sqrt(Ck * std::sqrt(Ck / Ce));
324  amrex::Print() << " equivalent Cs : " << Cs_equiv
325  << std::endl;
326  }
327  amrex::Print() << " isotropic mixing : " << mix_isotropic
328  << std::endl;
329  }
330 
331  if (rans_type != RANSType::None) {
332  if (rans_type == RANSType::kEqn) {
333  amrex::Print() << "Cmu0 : " << Cmu0 << std::endl;
334  amrex::Print() << "sigma_k : " << sigma_k << std::endl;
335  amrex::Print() << "Cb : " << Cb << std::endl;
336  amrex::Print() << "Rt_crit : " << Rt_crit << std::endl;
337  amrex::Print() << "Rt_min : " << Rt_min << std::endl;
338  amrex::Print() << "max_geom_lscale : " << l_g_max << std::endl;
339  }
340  }
341 
342  if ((les_type == LESType::Deardorff) ||
343  (rans_type == RANSType::kEqn)) {
344  if (theta_ref > 0) {
345  amrex::Print() << " reference theta : " << theta_ref << std::endl;
346  } else {
347  amrex::Print() << " reference theta : n/a" << std::endl;
348  }
349  }
350 
351  if ((les_type != LESType::None) || (rans_type != RANSType::None)) {
352  amrex::Print() << " Pr_t : " << Pr_t << std::endl;
353  amrex::Print() << " Sc_t : " << Sc_t << std::endl;
354  }
355 
356  if (pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
357  amrex::Print() << " pbl_mynn_A1 : " << pbl_mynn.A1 << std::endl;
358  amrex::Print() << " pbl_mynn_A2 : " << pbl_mynn.A2 << std::endl;
359  amrex::Print() << " pbl_mynn_B1 : " << pbl_mynn.B1 << std::endl;
360  amrex::Print() << " pbl_mynn_B2 : " << pbl_mynn.B2 << std::endl;
361  amrex::Print() << " pbl_mynn_C1 : " << pbl_mynn.C1 << std::endl;
362  amrex::Print() << " pbl_mynn_C2 : " << pbl_mynn.C2 << std::endl;
363  amrex::Print() << " pbl_mynn_C3 : " << pbl_mynn.C3 << std::endl;
364  amrex::Print() << " pbl_mynn_C4 : " << pbl_mynn.C4 << std::endl;
365  amrex::Print() << " pbl_mynn_C5 : " << pbl_mynn.C5 << std::endl;
366  } else if (pbl_type == PBLType::YSU) {
367  amrex::Print() << " pbl_ysu_coriolis_freq : "
368  << pbl_ysu_coriolis_freq << std::endl;
369  amrex::Print() << " pbl_ysu_use_consistent_coriolis : "
370  << pbl_ysu_use_consistent_coriolis << std::endl;
371  amrex::Print() << " pbl_ysu_force_over_water : "
372  << pbl_ysu_force_over_water << std::endl;
373  amrex::Print() << " pbl_ysu_land_Ribcr : "
374  << pbl_ysu_land_Ribcr << std::endl;
375  amrex::Print() << " pbl_ysu_unst_Ribcr : "
376  << pbl_ysu_unst_Ribcr << std::endl;
377  } else if (pbl_type == PBLType::MRF) {
378  amrex::Print() << " pbl_mrf_coriolis_freq : " << pbl_mrf_coriolis_freq
379  << std::endl;
380  amrex::Print() << " pbl_mrf_Ribcr : " << pbl_mrf_Ribcr
381  << std::endl;
382  amrex::Print() << " pbl_mrf_const_b : " << pbl_mrf_const_b
383  << std::endl;
384  amrex::Print() << " pbl_mrf_sf : " << pbl_mrf_sf
385  << std::endl;
386  amrex::Print() << " mrf_moistvars : " << mrf_moistvars
387  << std::endl;
388  }
389  }
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Real C4
Definition: ERF_MYNNStruct.H:50
amrex::Real C1
Definition: ERF_MYNNStruct.H:47
amrex::Real C3
Definition: ERF_MYNNStruct.H:49
amrex::Real C2
Definition: ERF_MYNNStruct.H:48
amrex::Real A2
Definition: ERF_MYNNStruct.H:44
amrex::Real B1
Definition: ERF_MYNNStruct.H:45
amrex::Real B2
Definition: ERF_MYNNStruct.H:46
amrex::Real C5
Definition: ERF_MYNNStruct.H:51
amrex::Real A1
Definition: ERF_MYNNStruct.H:43
bool smag2d
Definition: ERF_TurbStruct.H:404
amrex::Real sigma_k
Definition: ERF_TurbStruct.H:420
MYNNLevel25 pbl_mynn
Definition: ERF_TurbStruct.H:441
amrex::Real pbl_mrf_const_b
Definition: ERF_TurbStruct.H:476
amrex::Real Rt_min
Definition: ERF_TurbStruct.H:415
amrex::Real pbl_mrf_coriolis_freq
Definition: ERF_TurbStruct.H:474
amrex::Real Ri_crit
Definition: ERF_TurbStruct.H:431
RANSType rans_type
Definition: ERF_TurbStruct.H:434
amrex::Real pbl_mrf_Ribcr
Definition: ERF_TurbStruct.H:475
amrex::Real Ck
Definition: ERF_TurbStruct.H:409
amrex::Real Cmu0
Definition: ERF_TurbStruct.H:412
bool pbl_ysu_use_consistent_coriolis
Definition: ERF_TurbStruct.H:463
amrex::Real Cb
Definition: ERF_TurbStruct.H:413
amrex::Real pbl_mrf_sf
Definition: ERF_TurbStruct.H:477
bool use_Ri_correction
Definition: ERF_TurbStruct.H:430
amrex::Real pbl_ysu_land_Ribcr
Definition: ERF_TurbStruct.H:469
bool mrf_moistvars
Definition: ERF_TurbStruct.H:478
amrex::Real Cs
Definition: ERF_TurbStruct.H:403
amrex::Real Ce
Definition: ERF_TurbStruct.H:407
amrex::Real pbl_ysu_coriolis_freq
Definition: ERF_TurbStruct.H:460
bool pbl_ysu_force_over_water
Definition: ERF_TurbStruct.H:466
amrex::Real Pr_t
Definition: ERF_TurbStruct.H:395
amrex::Real theta_ref
Definition: ERF_TurbStruct.H:423
amrex::Real Sc_t
Definition: ERF_TurbStruct.H:399
amrex::Real Rt_crit
Definition: ERF_TurbStruct.H:414
amrex::Real pbl_ysu_unst_Ribcr
Definition: ERF_TurbStruct.H:471
bool mix_isotropic
Definition: ERF_TurbStruct.H:428
amrex::Real l_g_max
Definition: ERF_TurbStruct.H:416

◆ init_params()

void TurbChoice::init_params ( int  lev,
int  max_level,
std::string  pp_prefix 
)
inline
69  {
70  amrex::ParmParse pp(pp_prefix);
71 
72  // Which LES closure?
73  query_one_or_per_level(pp, "les_type", les_type, lev, max_level);
74 
75  // Handle 2-D Smag
76  if (les_type == LESType::Smagorinsky2D) {
77  les_type = LESType::Smagorinsky;
78  smag2d = true;
79  }
80 
81  // Which RANS closure?
82  query_one_or_per_level(pp, "rans_type", rans_type, lev, max_level);
83 
84  if ((rans_type != RANSType::None) && (les_type != LESType::None)) {
85  amrex::Error("Hybrid RANS-LES not implemented");
86  }
87 
88  // Which PBL Closure
89  query_one_or_per_level(pp, "pbl_type", pbl_type, lev, max_level);
90 
91  // Do some more stuff for PBL Modeling
92  if (pbl_type != PBLType::None) {
93  // Check for compatibility between PBL, LES, Molec Transport
94  if (les_type != LESType::None) {
95  amrex::Print() << "Selected a PBL model and an LES model: "
96  << "Using PBL for vertical transport, LES for horizontal"
97  << std::endl;
98  }
99  if (les_type == LESType::Smagorinsky) {
100  if (!smag2d)
101  amrex::Error("If using Smagorinsky with a PBL model, the 2-D "
102  "formulation should be used");
103  } else if (les_type == LESType::Deardorff) {
104  amrex::Error(
105  "It is not recommended to use Deardorff LES and a PBL model");
106  }
107 
108  if (pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
109  query_one_or_per_level(pp, "pbl_mynn_A1", pbl_mynn.A1, lev, max_level);
110  query_one_or_per_level(pp, "pbl_mynn_A2", pbl_mynn.A2, lev, max_level);
111  query_one_or_per_level(pp, "pbl_mynn_B1", pbl_mynn.B1, lev, max_level);
112  query_one_or_per_level(pp, "pbl_mynn_B2", pbl_mynn.B2, lev, max_level);
113  query_one_or_per_level(pp, "pbl_mynn_C1", pbl_mynn.C1, lev, max_level);
114  query_one_or_per_level(pp, "pbl_mynn_C2", pbl_mynn.C2, lev, max_level);
115  query_one_or_per_level(pp, "pbl_mynn_C3", pbl_mynn.C3, lev, max_level);
116  query_one_or_per_level(pp, "pbl_mynn_C4", pbl_mynn.C4, lev, max_level);
117  query_one_or_per_level(pp, "pbl_mynn_C5", pbl_mynn.C5, lev, max_level);
122  pp, "pbl_mynn_diffuse_moistvars", pbl_mynn.diffuse_moistvars, lev,
123  max_level);
125  pp, "pbl_mynn_SMmin", pbl_mynn.SMmin, lev, max_level);
127  pp, "pbl_mynn_SMmax", pbl_mynn.SMmax, lev, max_level);
129  pp, "pbl_mynn_SHmin", pbl_mynn.SHmin, lev, max_level);
131  pp, "pbl_mynn_SHmax", pbl_mynn.SHmax, lev, max_level);
133  pp, "pbl_mynn_SQfactor", pbl_mynn.SQfac, lev, max_level);
134  } else if (pbl_type == PBLType::YSU) {
136  pp, "pbl_ysu_coriolis_freq", pbl_ysu_coriolis_freq, lev, max_level);
138  pp, "pbl_ysu_use_consistent_coriolis",
139  pbl_ysu_use_consistent_coriolis, lev, max_level);
141  pp, "pbl_ysu_force_over_water", pbl_ysu_force_over_water, lev,
142  max_level);
144  pp, "pbl_ysu_land_Ribcr", pbl_ysu_land_Ribcr, lev, max_level);
146  pp, "pbl_ysu_unst_Ribcr", pbl_ysu_unst_Ribcr, lev, max_level);
147  } else if (pbl_type == PBLType::MRF) {
149  pp, "pbl_mrf_coriolis_freq", pbl_mrf_coriolis_freq, lev, max_level);
151  pp, "pbl_mrf_Ribcr", pbl_mrf_Ribcr, lev, max_level);
153  pp, "pbl_mrf_const_b", pbl_mrf_const_b, lev, max_level);
154  query_one_or_per_level(pp, "pbl_mrf_sf", pbl_mrf_sf, lev, max_level);
156  pp, "mrf_moistvars", mrf_moistvars, lev, max_level);
157  } else if (pbl_type == PBLType::SHOC) {
158 #ifndef ERF_USE_SHOC
159  amrex::Abort("You set use_shoc to true but didn't build with SHOC; you must rebuild the executable");
160 #endif
161  std::string zlo_bc = "none";
162  amrex::ParmParse pp_bc("zlo");
163  pp_bc.get("type",zlo_bc);
164  if (amrex::toLower(zlo_bc) != "surface_layer") {
165  amrex::Abort("You must use the surface_layer BC at zlo with SHOC.");
166  }
167  }
168  }
169 
170  // Flags for QKE/TKE equation
171  if (pbl_type == PBLType::MYJ || pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
172  // Add sources/sinks to QKE/TKE? (MYJ does this inline)
173  if (pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
174  use_pbl_tke = true;
175  }
176  // Advect QKE/TKE?
177  query_one_or_per_level(pp, "advect_tke" , advect_tke , lev, max_level);
178  // Apply numerical diffusion to QKE/TKE?
179  query_one_or_per_level(pp, "diffuse_tke_3D", diffuse_tke_3D, lev, max_level);
180  }
181 
182  // LES constants...
183  query_one_or_per_level(pp, "Cs", Cs, lev, max_level);
184 
185  query_one_or_per_level(pp, "Pr_t", Pr_t, lev, max_level);
186  query_one_or_per_level(pp, "Sc_t", Sc_t, lev, max_level);
187 
188  // Compute relevant forms of diffusion parameters
189  Pr_t_inv = one / Pr_t;
190  Sc_t_inv = one / Sc_t;
191 
192  if (les_type == LESType::Deardorff) {
193  query_one_or_per_level(pp, "Ck", Ck, lev, max_level);
194  query_one_or_per_level(pp, "Ce", Ce, lev, max_level);
195  query_one_or_per_level(pp, "Ce_wall", Ce_wall, lev, max_level);
196  }
197 
198  // To quantify atmospheric stability for subgrid modeling
199  query_one_or_per_level(pp, "thermal_stratification", strat_type, lev, max_level);
200  if (strat_type == StratType::theta) {
201  amrex::Print() << "Thermal stratification based on gradient of potential temperature" << std::endl;
202  } else if (strat_type == StratType::thetav) {
203  amrex::Print() << "Thermal stratification based on gradient of virtual potential temperature" << std::endl;
204  } else if (strat_type == StratType::thetal) {
205  amrex::Print() << "Thermal stratification based on gradient of linearized liquid-water potential temperature" << std::endl;
206  }
207 
208  // k-eqn constants
209  query_one_or_per_level(pp, "Cmu0", Cmu0, lev, max_level);
210  query_one_or_per_level(pp, "Cb", Cb, lev, max_level);
211  query_one_or_per_level(pp, "Rt_crit", Rt_crit, lev, max_level);
212  query_one_or_per_level(pp, "Rt_min", Rt_min, lev, max_level);
213  query_one_or_per_level(pp, "max_geom_lscale", l_g_max, lev, max_level);
214  query_one_or_per_level(pp, "dirichlet_k", dirichlet_k, lev, max_level);
215 
216  // Common inputs (LES or RANS)
217  if (!query_one_or_per_level(pp, "sigma_k", sigma_k, lev, max_level) && rans_type == RANSType::kEqn) {
218  amrex::Print() << "Overriding default sigma_k for k-eqn RANS" << std::endl;
219  sigma_k = one;
220  };
221  query_one_or_per_level(pp, "theta_ref", theta_ref, lev, max_level);
222 
223  query_one_or_per_level(pp, "mix_isotropic", mix_isotropic, lev, max_level);
224  query_one_or_per_level(pp, "use_Ri_correction", use_Ri_correction, lev, max_level);
225  query_one_or_per_level(pp, "Ri_crit", Ri_crit, lev, max_level);
226 
227  // Set common flags
228  use_kturb =
229  ((les_type != LESType::None) || (rans_type != RANSType::None) ||
230  (pbl_type != PBLType::None));
231  use_keqn =
232  ((les_type == LESType::Deardorff) || (rans_type == RANSType::kEqn));
233  use_tke =
234  ((les_type == LESType::Deardorff) || (rans_type == RANSType::kEqn) ||
235  (pbl_type == PBLType::MYJ) || (pbl_type == PBLType::MYNN25) ||
236  (pbl_type == PBLType::MYNNEDMF) || (pbl_type == PBLType::SHOC));
237 
238  if (use_tke) {
239  query_one_or_per_level(pp, "init_tke_from_ustar", init_tke_from_ustar, lev, max_level);
240  }
241 
242  // Validate inputs
243  if (les_type == LESType::Smagorinsky) {
244  if (Cs == 0) {
245  amrex::Error("Need to specify Cs for Smagorsinky LES");
246  }
247  if (smag2d && mix_isotropic) {
248  amrex::Print() << "Turning off mix_isotropic for 2-D Smagorinsky" << std::endl;
249  mix_isotropic = false;
250  }
251  }
252  }
constexpr amrex::Real one
Definition: ERF_Constants.H:7
ParmParse pp("prob")
int query_one_or_per_level(const amrex::ParmParse &pp, const char *query_string, T &query_var, const int lev, const int maxlev)
Definition: ERF_TurbStruct.H:16
@ theta
Definition: ERF_MM5.H:20
amrex::Real SMmax
Definition: ERF_MYNNStruct.H:58
amrex::Real SHmax
Definition: ERF_MYNNStruct.H:60
amrex::Real SQfac
Definition: ERF_MYNNStruct.H:54
amrex::Real SHmin
Definition: ERF_MYNNStruct.H:59
amrex::Real SMmin
Definition: ERF_MYNNStruct.H:57
bool diffuse_moistvars
Definition: ERF_MYNNStruct.H:65
void init_coeffs(amrex::Real A1_lvl25, amrex::Real A2_lvl25, amrex::Real B1, amrex::Real B2, amrex::Real C1, amrex::Real C2, amrex::Real C3, amrex::Real, amrex::Real C5)
Definition: ERF_MYNNStruct.H:69
bool use_pbl_tke
Definition: ERF_TurbStruct.H:449
amrex::Real Sc_t_inv
Definition: ERF_TurbStruct.H:400
bool use_keqn
Definition: ERF_TurbStruct.H:446
MYNNLevel2 pbl_mynn_level2
Definition: ERF_TurbStruct.H:442
StratType strat_type
Definition: ERF_TurbStruct.H:425
bool dirichlet_k
Definition: ERF_TurbStruct.H:436
bool advect_tke
Definition: ERF_TurbStruct.H:481
bool init_tke_from_ustar
Definition: ERF_TurbStruct.H:456
bool use_tke
Definition: ERF_TurbStruct.H:451
bool diffuse_tke_3D
Definition: ERF_TurbStruct.H:483
amrex::Real Pr_t_inv
Definition: ERF_TurbStruct.H:396
bool use_kturb
Definition: ERF_TurbStruct.H:445
Here is the call graph for this function:

Member Data Documentation

◆ advect_tke

bool TurbChoice::advect_tke = true

Referenced by erf_slow_rhs_post(), and init_params().

◆ Cb

amrex::Real TurbChoice::Cb = amrex::Real(0.35)

◆ Ce

amrex::Real TurbChoice::Ce = amrex::Real(0.93)

◆ Ce_wall

◆ Ck

◆ Cmu0

◆ Cs

◆ diffuse_tke_3D

bool TurbChoice::diffuse_tke_3D = true

Referenced by init_params(), and make_sources().

◆ dirichlet_k

bool TurbChoice::dirichlet_k = false

◆ init_tke_from_ustar

bool TurbChoice::init_tke_from_ustar = false

Referenced by init_params().

◆ l_g_max

amrex::Real TurbChoice::l_g_max = amrex::Real(30.0)

◆ les_type

◆ mix_isotropic

bool TurbChoice::mix_isotropic = true

◆ mrf_moistvars

bool TurbChoice::mrf_moistvars = false

Referenced by display(), and init_params().

◆ pbl_mrf_const_b

amrex::Real TurbChoice::pbl_mrf_const_b = amrex::Real(7.8)

◆ pbl_mrf_coriolis_freq

amrex::Real TurbChoice::pbl_mrf_coriolis_freq = amrex::Real(1.0e-4)

Referenced by display(), and init_params().

◆ pbl_mrf_Ribcr

amrex::Real TurbChoice::pbl_mrf_Ribcr = amrex::Real(0.5)

◆ pbl_mrf_sf

amrex::Real TurbChoice::pbl_mrf_sf = amrex::Real(0.1)

◆ pbl_mynn

◆ pbl_mynn_level2

MYNNLevel2 TurbChoice::pbl_mynn_level2

◆ pbl_type

◆ pbl_ysu_coriolis_freq

amrex::Real TurbChoice::pbl_ysu_coriolis_freq
Initial value:
=
amrex::Real(1.0e-4)

Referenced by ComputeDiffusivityYSU(), display(), and init_params().

◆ pbl_ysu_force_over_water

bool TurbChoice::pbl_ysu_force_over_water
Initial value:
=
false

Referenced by ComputeDiffusivityYSU(), display(), and init_params().

◆ pbl_ysu_land_Ribcr

amrex::Real TurbChoice::pbl_ysu_land_Ribcr
Initial value:
=
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12

Referenced by ComputeDiffusivityYSU(), display(), and init_params().

◆ pbl_ysu_unst_Ribcr

amrex::Real TurbChoice::pbl_ysu_unst_Ribcr
Initial value:
=
constexpr amrex::Real zero
Definition: ERF_Constants.H:6

Referenced by ComputeDiffusivityYSU(), display(), and init_params().

◆ pbl_ysu_use_consistent_coriolis

bool TurbChoice::pbl_ysu_use_consistent_coriolis
Initial value:
=
false

Referenced by display(), and init_params().

◆ Pr_t

amrex::Real TurbChoice::Pr_t = one / three

Referenced by display(), and init_params().

◆ Pr_t_inv

◆ rans_type

◆ Ri_crit

◆ Rt_crit

amrex::Real TurbChoice::Rt_crit = -one

◆ Rt_min

amrex::Real TurbChoice::Rt_min = -three

◆ Sc_t

amrex::Real TurbChoice::Sc_t = one

Referenced by display(), and init_params().

◆ Sc_t_inv

◆ sigma_k

◆ smag2d

bool TurbChoice::smag2d = false

◆ strat_type

StratType TurbChoice::strat_type = StratType::theta

◆ theta_ref

◆ use_keqn

bool TurbChoice::use_keqn
Initial value:
=
false

Referenced by erf_slow_rhs_pre(), and init_params().

◆ use_kturb

◆ use_pbl_tke

bool TurbChoice::use_pbl_tke
Initial value:
=
false

Referenced by init_params().

◆ use_Ri_correction

bool TurbChoice::use_Ri_correction = true

◆ use_tke

bool TurbChoice::use_tke = false

The documentation for this struct was generated from the following file: