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 = LESType::None
 
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 = RANSType::None
 
bool dirichlet_k = false
 
PBLType pbl_type = PBLType::None
 
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 tke_min = amrex::Real(1.e-6)
 
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
275  {
276  // BC compatibility
277  if ( ( (pbl_type == PBLType::MYNN25) ||
278  (pbl_type == PBLType::MYNNEDMF) ||
279  (pbl_type == PBLType::YSU) ||
280  (pbl_type == PBLType::MRF)
281  ) &&
282  phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::surface_layer ) {
283  amrex::Abort("MYNN2.5/MYNNEDMF/YSU/MRF PBL Model requires MOST at lower boundary");
284  }
285  if ( (les_type == LESType::Deardorff) && (Ce_wall > 0) &&
286  (phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::surface_layer) &&
287  (phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::slip_wall) &&
288  (phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::no_slip_wall) )
289  {
290  amrex::Warning("Deardorff LES assumes wall at zlo when applying Ce_wall");
291  }
292  }
@ no_slip_wall
@ surface_layer
PBLType pbl_type
Definition: ERF_TurbStruct.H:459
amrex::Real Ce_wall
Definition: ERF_TurbStruct.H:428
LESType les_type
Definition: ERF_TurbStruct.H:412

◆ display()

void TurbChoice::display ( int  lev)
inline
295  {
296  amrex::Print() << "Turbulence Settings at level " << lev << std::endl;
297 
298  if (
299  les_type == LESType::None && rans_type == RANSType::None &&
300  pbl_type == PBLType::None) {
301  amrex::Print() << " Using DNS model at level " << lev << std::endl;
302  } else if (les_type == LESType::Smagorinsky) {
303  if (smag2d) {
304  amrex::Print() << " Using 2D Smagorinsky LES model at level " << lev << std::endl;
305  } else {
306  amrex::Print() << " Using Smagorinsky LES model at level " << lev << std::endl;
307  }
308  if (use_Ri_correction) {
309  amrex::Print() << " Smagorinsky uses Richardson number correction with Ri_crit = "
310  << Ri_crit << std::endl;
311  }
312  } else if (les_type == LESType::Deardorff) {
313  amrex::Print() << " Using Deardorff LES model at level " << lev << std::endl;
314  } else if (rans_type == RANSType::kEqn) {
315  amrex::Print()
316  << " Using Axell & Liungman one-equation RANS k model at level " << lev
317  << std::endl;
318  } else if (pbl_type == PBLType::MYJ) {
319  amrex::Print() << " Using MYJ PBL model at level " << lev << std::endl;
320  } else if (pbl_type == PBLType::MYNN25) {
321  amrex::Print() << " Using MYNN2.5 PBL model at level " << lev << std::endl;
322  } else if (pbl_type == PBLType::MYNNEDMF) {
323  amrex::Print() << " Using MYNNEDMF PBL model at level " << lev << std::endl;
324  } else if (pbl_type == PBLType::YSU) {
325  amrex::Print() << " Using YSU PBL model at level " << lev << std::endl;
326  } else if (pbl_type == PBLType::MRF) {
327  amrex::Print() << " Using MRF PBL model at level " << lev << std::endl;
328  } else {
329  amrex::Error("Unknown turbulence model");
330  }
331 
332  if (les_type != LESType::None) {
333  if (les_type == LESType::Smagorinsky) {
334  amrex::Print() << " Cs : " << Cs << std::endl;
335  }
336  if (les_type == LESType::Deardorff) {
337  amrex::Print() << " Ce : " << Ce << std::endl;
338  amrex::Print() << " Ce at wall : " << Ce_wall << std::endl;
339  amrex::Print() << " Ck : " << Ck << std::endl;
340  amrex::Print() << " sigma_k : " << sigma_k << std::endl;
341 
342  // Sullivan et al 1994, Eqn 14
343  amrex::Real Cs_equiv = std::sqrt(Ck * std::sqrt(Ck / Ce));
344  amrex::Print() << " equivalent Cs : " << Cs_equiv
345  << std::endl;
346  }
347  amrex::Print() << " isotropic mixing : " << mix_isotropic
348  << std::endl;
349  }
350 
351  if (rans_type != RANSType::None) {
352  if (rans_type == RANSType::kEqn) {
353  amrex::Print() << "Cmu0 : " << Cmu0 << std::endl;
354  amrex::Print() << "sigma_k : " << sigma_k << std::endl;
355  amrex::Print() << "Cb : " << Cb << std::endl;
356  amrex::Print() << "Rt_crit : " << Rt_crit << std::endl;
357  amrex::Print() << "Rt_min : " << Rt_min << std::endl;
358  amrex::Print() << "max_geom_lscale : " << l_g_max << std::endl;
359  }
360  }
361 
362  if ((les_type == LESType::Deardorff) ||
363  (rans_type == RANSType::kEqn)) {
364  if (theta_ref > 0) {
365  amrex::Print() << " reference theta : " << theta_ref << std::endl;
366  } else {
367  amrex::Print() << " reference theta : n/a" << std::endl;
368  }
369  }
370 
371  if ((les_type != LESType::None) || (rans_type != RANSType::None)) {
372  amrex::Print() << " Pr_t : " << Pr_t << std::endl;
373  amrex::Print() << " Sc_t : " << Sc_t << std::endl;
374  }
375 
376  if (pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
377  amrex::Print() << " pbl_mynn_A1 : " << pbl_mynn.A1 << std::endl;
378  amrex::Print() << " pbl_mynn_A2 : " << pbl_mynn.A2 << std::endl;
379  amrex::Print() << " pbl_mynn_B1 : " << pbl_mynn.B1 << std::endl;
380  amrex::Print() << " pbl_mynn_B2 : " << pbl_mynn.B2 << std::endl;
381  amrex::Print() << " pbl_mynn_C1 : " << pbl_mynn.C1 << std::endl;
382  amrex::Print() << " pbl_mynn_C2 : " << pbl_mynn.C2 << std::endl;
383  amrex::Print() << " pbl_mynn_C3 : " << pbl_mynn.C3 << std::endl;
384  amrex::Print() << " pbl_mynn_C4 : " << pbl_mynn.C4 << std::endl;
385  amrex::Print() << " pbl_mynn_C5 : " << pbl_mynn.C5 << std::endl;
386  } else if (pbl_type == PBLType::YSU) {
387  amrex::Print() << " pbl_ysu_coriolis_freq : "
388  << pbl_ysu_coriolis_freq << std::endl;
389  amrex::Print() << " pbl_ysu_use_consistent_coriolis : "
390  << pbl_ysu_use_consistent_coriolis << std::endl;
391  amrex::Print() << " pbl_ysu_force_over_water : "
392  << pbl_ysu_force_over_water << std::endl;
393  amrex::Print() << " pbl_ysu_land_Ribcr : "
394  << pbl_ysu_land_Ribcr << std::endl;
395  amrex::Print() << " pbl_ysu_unst_Ribcr : "
396  << pbl_ysu_unst_Ribcr << std::endl;
397  } else if (pbl_type == PBLType::MRF) {
398  amrex::Print() << " pbl_mrf_coriolis_freq : " << pbl_mrf_coriolis_freq
399  << std::endl;
400  amrex::Print() << " pbl_mrf_Ribcr : " << pbl_mrf_Ribcr
401  << std::endl;
402  amrex::Print() << " pbl_mrf_const_b : " << pbl_mrf_const_b
403  << std::endl;
404  amrex::Print() << " pbl_mrf_sf : " << pbl_mrf_sf
405  << std::endl;
406  amrex::Print() << " mrf_moistvars : " << mrf_moistvars
407  << std::endl;
408  }
409  }
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:424
amrex::Real sigma_k
Definition: ERF_TurbStruct.H:440
MYNNLevel25 pbl_mynn
Definition: ERF_TurbStruct.H:461
amrex::Real pbl_mrf_const_b
Definition: ERF_TurbStruct.H:499
amrex::Real Rt_min
Definition: ERF_TurbStruct.H:435
amrex::Real pbl_mrf_coriolis_freq
Definition: ERF_TurbStruct.H:497
amrex::Real Ri_crit
Definition: ERF_TurbStruct.H:451
RANSType rans_type
Definition: ERF_TurbStruct.H:454
amrex::Real pbl_mrf_Ribcr
Definition: ERF_TurbStruct.H:498
amrex::Real Ck
Definition: ERF_TurbStruct.H:429
amrex::Real Cmu0
Definition: ERF_TurbStruct.H:432
bool pbl_ysu_use_consistent_coriolis
Definition: ERF_TurbStruct.H:486
amrex::Real Cb
Definition: ERF_TurbStruct.H:433
amrex::Real pbl_mrf_sf
Definition: ERF_TurbStruct.H:500
bool use_Ri_correction
Definition: ERF_TurbStruct.H:450
amrex::Real pbl_ysu_land_Ribcr
Definition: ERF_TurbStruct.H:492
bool mrf_moistvars
Definition: ERF_TurbStruct.H:501
amrex::Real Cs
Definition: ERF_TurbStruct.H:423
amrex::Real Ce
Definition: ERF_TurbStruct.H:427
amrex::Real pbl_ysu_coriolis_freq
Definition: ERF_TurbStruct.H:483
bool pbl_ysu_force_over_water
Definition: ERF_TurbStruct.H:489
amrex::Real Pr_t
Definition: ERF_TurbStruct.H:415
amrex::Real theta_ref
Definition: ERF_TurbStruct.H:443
amrex::Real Sc_t
Definition: ERF_TurbStruct.H:419
amrex::Real Rt_crit
Definition: ERF_TurbStruct.H:434
amrex::Real pbl_ysu_unst_Ribcr
Definition: ERF_TurbStruct.H:494
bool mix_isotropic
Definition: ERF_TurbStruct.H:448
amrex::Real l_g_max
Definition: ERF_TurbStruct.H:436

◆ init_params()

void TurbChoice::init_params ( int  lev,
int  max_level,
std::string  pp_prefix 
)
inline
85  {
86  amrex::ParmParse pp(pp_prefix);
87 
88  // Which LES closure?
89  query_one_or_per_level(pp, "les_type", les_type, lev, max_level);
90 
91  // Handle 2-D Smag
92  if (les_type == LESType::Smagorinsky2D) {
93  les_type = LESType::Smagorinsky;
94  smag2d = true;
95  }
96 
97  // Which RANS closure?
98  query_one_or_per_level(pp, "rans_type", rans_type, lev, max_level);
99 
100  if ((rans_type != RANSType::None) && (les_type != LESType::None)) {
101  amrex::Error("Hybrid RANS-LES not implemented");
102  }
103 
104  // Which PBL Closure
105  query_one_or_per_level(pp, "pbl_type", pbl_type, lev, max_level);
106 
107  // Do some more stuff for PBL Modeling
108  if (pbl_type != PBLType::None) {
109  // Check for compatibility between PBL, LES, Molec Transport
110  if (les_type != LESType::None) {
111  amrex::Print() << "Selected a PBL model and an LES model: "
112  << "Using PBL for vertical transport, LES for horizontal"
113  << std::endl;
114  }
115  if (les_type == LESType::Smagorinsky) {
116  if (!smag2d)
117  amrex::Error("If using Smagorinsky with a PBL model, the 2-D "
118  "formulation should be used");
119  } else if (les_type == LESType::Deardorff) {
120  amrex::Error(
121  "It is not recommended to use Deardorff LES and a PBL model");
122  }
123 
124  if (pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
125  query_one_or_per_level(pp, "pbl_mynn_A1", pbl_mynn.A1, lev, max_level);
126  query_one_or_per_level(pp, "pbl_mynn_A2", pbl_mynn.A2, lev, max_level);
127  query_one_or_per_level(pp, "pbl_mynn_B1", pbl_mynn.B1, lev, max_level);
128  query_one_or_per_level(pp, "pbl_mynn_B2", pbl_mynn.B2, lev, max_level);
129  query_one_or_per_level(pp, "pbl_mynn_C1", pbl_mynn.C1, lev, max_level);
130  query_one_or_per_level(pp, "pbl_mynn_C2", pbl_mynn.C2, lev, max_level);
131  query_one_or_per_level(pp, "pbl_mynn_C3", pbl_mynn.C3, lev, max_level);
132  query_one_or_per_level(pp, "pbl_mynn_C4", pbl_mynn.C4, lev, max_level);
133  query_one_or_per_level(pp, "pbl_mynn_C5", pbl_mynn.C5, lev, max_level);
138  pp, "pbl_mynn_diffuse_moistvars", pbl_mynn.diffuse_moistvars, lev,
139  max_level);
141  pp, "pbl_mynn_SMmin", pbl_mynn.SMmin, lev, max_level);
143  pp, "pbl_mynn_SMmax", pbl_mynn.SMmax, lev, max_level);
145  pp, "pbl_mynn_SHmin", pbl_mynn.SHmin, lev, max_level);
147  pp, "pbl_mynn_SHmax", pbl_mynn.SHmax, lev, max_level);
149  pp, "pbl_mynn_SQfactor", pbl_mynn.SQfac, lev, max_level);
150  } else if (pbl_type == PBLType::YSU) {
152  pp, "pbl_ysu_coriolis_freq", pbl_ysu_coriolis_freq, lev, max_level);
154  pp, "pbl_ysu_use_consistent_coriolis",
155  pbl_ysu_use_consistent_coriolis, lev, max_level);
157  pp, "pbl_ysu_force_over_water", pbl_ysu_force_over_water, lev,
158  max_level);
160  pp, "pbl_ysu_land_Ribcr", pbl_ysu_land_Ribcr, lev, max_level);
162  pp, "pbl_ysu_unst_Ribcr", pbl_ysu_unst_Ribcr, lev, max_level);
163  } else if (pbl_type == PBLType::MRF) {
165  pp, "pbl_mrf_coriolis_freq", pbl_mrf_coriolis_freq, lev, max_level);
167  pp, "pbl_mrf_Ribcr", pbl_mrf_Ribcr, lev, max_level);
169  pp, "pbl_mrf_const_b", pbl_mrf_const_b, lev, max_level);
170  query_one_or_per_level(pp, "pbl_mrf_sf", pbl_mrf_sf, lev, max_level);
172  pp, "mrf_moistvars", mrf_moistvars, lev, max_level);
173  } else if (pbl_type == PBLType::SHOC) {
174 #ifndef ERF_USE_SHOC
175  amrex::Abort("You set use_shoc to true but didn't build with SHOC; you must rebuild the executable");
176 #endif
177  std::string zlo_bc = "none";
178  amrex::ParmParse pp_bc("zlo");
179  pp_bc.get("type",zlo_bc);
180  if (amrex::toLower(zlo_bc) != "surface_layer") {
181  amrex::Abort("You must use the surface_layer BC at zlo with SHOC.");
182  }
183  }
184  }
185 
186  // Flags for QKE/TKE equation
187  if (pbl_type == PBLType::MYJ || pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
188  // Add sources/sinks to QKE/TKE? (MYJ does this inline)
189  if (pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
190  use_pbl_tke = true;
191  }
192  // Advect QKE/TKE?
193  query_one_or_per_level(pp, "advect_tke" , advect_tke , lev, max_level);
194  // Apply numerical diffusion to QKE/TKE?
195  query_one_or_per_level(pp, "diffuse_tke_3D", diffuse_tke_3D, lev, max_level);
196  }
197 
198  // There is a default value of 1e-6 but the user can override the value here,
199  // and in addition can add perturbational values. and in addition can add perturbational values.
200  pp.query("tke_min",tke_min);
201 
202  // LES constants...
203  query_one_or_per_level(pp, "Cs", Cs, lev, max_level);
204 
205  query_one_or_per_level(pp, "Pr_t", Pr_t, lev, max_level);
206  query_one_or_per_level(pp, "Sc_t", Sc_t, lev, max_level);
207 
208  // Compute relevant forms of diffusion parameters
209  Pr_t_inv = one / Pr_t;
210  Sc_t_inv = one / Sc_t;
211 
212  if (les_type == LESType::Deardorff) {
213  query_one_or_per_level(pp, "Ck", Ck, lev, max_level);
214  query_one_or_per_level(pp, "Ce", Ce, lev, max_level);
215  query_one_or_per_level(pp, "Ce_wall", Ce_wall, lev, max_level);
216  }
217 
218  // To quantify atmospheric stability for subgrid modeling
219  query_one_or_per_level(pp, "thermal_stratification", strat_type, lev, max_level);
220  if (strat_type == StratType::theta) {
221  amrex::Print() << "Thermal stratification based on gradient of potential temperature" << std::endl;
222  } else if (strat_type == StratType::thetav) {
223  amrex::Print() << "Thermal stratification based on gradient of virtual potential temperature" << std::endl;
224  } else if (strat_type == StratType::thetal) {
225  amrex::Print() << "Thermal stratification based on gradient of linearized liquid-water potential temperature" << std::endl;
226  }
227 
228  // k-eqn constants
229  query_one_or_per_level(pp, "Cmu0", Cmu0, lev, max_level);
230  query_one_or_per_level(pp, "Cb", Cb, lev, max_level);
231  query_one_or_per_level(pp, "Rt_crit", Rt_crit, lev, max_level);
232  query_one_or_per_level(pp, "Rt_min", Rt_min, lev, max_level);
233  query_one_or_per_level(pp, "max_geom_lscale", l_g_max, lev, max_level);
234  query_one_or_per_level(pp, "dirichlet_k", dirichlet_k, lev, max_level);
235 
236  // Common inputs (LES or RANS)
237  if (!query_one_or_per_level(pp, "sigma_k", sigma_k, lev, max_level) && rans_type == RANSType::kEqn) {
238  amrex::Print() << "Overriding default sigma_k for k-eqn RANS" << std::endl;
239  sigma_k = one;
240  };
241  query_one_or_per_level(pp, "theta_ref", theta_ref, lev, max_level);
242 
243  query_one_or_per_level(pp, "mix_isotropic", mix_isotropic, lev, max_level);
244  query_one_or_per_level(pp, "use_Ri_correction", use_Ri_correction, lev, max_level);
245  query_one_or_per_level(pp, "Ri_crit", Ri_crit, lev, max_level);
246 
247  // Set common flags
248  use_kturb =
249  ((les_type != LESType::None) || (rans_type != RANSType::None) ||
250  (pbl_type != PBLType::None));
251  use_keqn =
252  ((les_type == LESType::Deardorff) || (rans_type == RANSType::kEqn));
253  use_tke =
254  ((les_type == LESType::Deardorff) || (rans_type == RANSType::kEqn) ||
255  (pbl_type == PBLType::MYJ) || (pbl_type == PBLType::MYNN25) ||
256  (pbl_type == PBLType::MYNNEDMF) || (pbl_type == PBLType::SHOC));
257 
258  if (use_tke) {
259  query_one_or_per_level(pp, "init_tke_from_ustar", init_tke_from_ustar, lev, max_level);
260  }
261 
262  // Validate inputs
263  if (les_type == LESType::Smagorinsky) {
264  if (Cs == 0) {
265  amrex::Error("Need to specify Cs for Smagorsinky LES");
266  }
267  if (smag2d && mix_isotropic) {
268  amrex::Print() << "Turning off mix_isotropic for 2-D Smagorinsky" << std::endl;
269  mix_isotropic = false;
270  }
271  }
272  }
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
amrex::Real tke_min
Definition: ERF_TurbStruct.H:479
bool use_pbl_tke
Definition: ERF_TurbStruct.H:469
amrex::Real Sc_t_inv
Definition: ERF_TurbStruct.H:420
bool use_keqn
Definition: ERF_TurbStruct.H:466
MYNNLevel2 pbl_mynn_level2
Definition: ERF_TurbStruct.H:462
StratType strat_type
Definition: ERF_TurbStruct.H:445
bool dirichlet_k
Definition: ERF_TurbStruct.H:456
bool advect_tke
Definition: ERF_TurbStruct.H:504
bool init_tke_from_ustar
Definition: ERF_TurbStruct.H:476
bool use_tke
Definition: ERF_TurbStruct.H:471
bool diffuse_tke_3D
Definition: ERF_TurbStruct.H:506
amrex::Real Pr_t_inv
Definition: ERF_TurbStruct.H:416
bool use_kturb
Definition: ERF_TurbStruct.H:465
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

amrex::Real TurbChoice::Cmu0 = amrex::Real(0.5562)

◆ Cs

◆ diffuse_tke_3D

bool TurbChoice::diffuse_tke_3D = true

Referenced by init_params(), and make_sources().

◆ dirichlet_k

bool TurbChoice::dirichlet_k = false

Referenced by init_params().

◆ 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

◆ tke_min

amrex::Real TurbChoice::tke_min = amrex::Real(1.e-6)

Referenced by init_params().

◆ 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: