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 = amrex::Real(1.0) / amrex::Real(3.0)
 
amrex::Real Pr_t_inv = amrex::Real(3.0)
 
amrex::Real Sc_t = 1.0
 
amrex::Real Sc_t_inv = 1.0
 
amrex::Real Cs = 0.0
 
bool smag2d = false
 
amrex::Real Ce = 0.93
 
amrex::Real Ce_wall = 0.0
 
amrex::Real Ck = 0.1
 
amrex::Real Cmu0 = 0.5562
 
amrex::Real Cb = 0.35
 
amrex::Real Rt_crit = -1.0
 
amrex::Real Rt_min = -3.0
 
amrex::Real l_g_max = 30.0
 
amrex::Real sigma_k = 0.5
 
amrex::Real theta_ref = 0.0
 
StratType strat_type = StratType::theta
 
bool mix_isotropic = true
 
bool use_Ri_correction = true
 
amrex::Real Ri_crit = 0.25
 
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 = 1.0e-4
 
amrex::Real pbl_mrf_Ribcr
 
amrex::Real pbl_mrf_const_b
 
amrex::Real pbl_mrf_sf
 
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
234  {
235  // BC compatibility
236  if ( ( (pbl_type == PBLType::MYNN25) ||
237  (pbl_type == PBLType::MYNNEDMF) ||
238  (pbl_type == PBLType::YSU) ||
239  (pbl_type == PBLType::MRF)
240  ) &&
241  phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::surface_layer ) {
242  amrex::Abort("MYNN2.5/MYNNEDMF/YSU/MRF PBL Model requires MOST at lower boundary");
243  }
244  if ( (les_type == LESType::Deardorff) && (Ce_wall > 0) &&
245  (phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::surface_layer) &&
246  (phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::slip_wall) &&
247  (phys_bc_type[amrex::Orientation(amrex::Direction::z,amrex::Orientation::low)] != ERF_BC::no_slip_wall) )
248  {
249  amrex::Warning("Deardorff LES assumes wall at zlo when applying Ce_wall");
250  }
251  }
@ no_slip_wall
@ surface_layer
PBLType pbl_type
Definition: ERF_TurbStruct.H:418
amrex::Real Ce_wall
Definition: ERF_TurbStruct.H:387
LESType les_type
Definition: ERF_TurbStruct.H:371

◆ display()

void TurbChoice::display ( int  lev)
inline
254  {
255  amrex::Print() << "Turbulence Settings at level " << lev << std::endl;
256 
257  if (
258  les_type == LESType::None && rans_type == RANSType::None &&
259  pbl_type == PBLType::None) {
260  amrex::Print() << " Using DNS model at level " << lev << std::endl;
261  } else if (les_type == LESType::Smagorinsky) {
262  if (smag2d) {
263  amrex::Print() << " Using 2D Smagorinsky LES model at level " << lev << std::endl;
264  } else {
265  amrex::Print() << " Using Smagorinsky LES model at level " << lev << std::endl;
266  }
267  if (use_Ri_correction) {
268  amrex::Print() << " Smagorinsky uses Richardson number correction with Ri_crit = "
269  << Ri_crit << std::endl;
270  }
271  } else if (les_type == LESType::Deardorff) {
272  amrex::Print() << " Using Deardorff LES model at level " << lev << std::endl;
273  } else if (rans_type == RANSType::kEqn) {
274  amrex::Print()
275  << " Using Axell & Liungman one-equation RANS k model at level " << lev
276  << std::endl;
277  } else if (pbl_type == PBLType::MYJ) {
278  amrex::Print() << " Using MYJ PBL model at level " << lev << std::endl;
279  } else if (pbl_type == PBLType::MYNN25) {
280  amrex::Print() << " Using MYNN2.5 PBL model at level " << lev << std::endl;
281  } else if (pbl_type == PBLType::MYNNEDMF) {
282  amrex::Print() << " Using MYNNEDMF PBL model at level " << lev << std::endl;
283  } else if (pbl_type == PBLType::YSU) {
284  amrex::Print() << " Using YSU PBL model at level " << lev << std::endl;
285  } else if (pbl_type == PBLType::MRF) {
286  amrex::Print() << " Using MRF PBL model at level " << lev << std::endl;
287  } else {
288  amrex::Error("Unknown turbulence model");
289  }
290 
291  if (les_type != LESType::None) {
292  if (les_type == LESType::Smagorinsky) {
293  amrex::Print() << " Cs : " << Cs << std::endl;
294  }
295  if (les_type == LESType::Deardorff) {
296  amrex::Print() << " Ce : " << Ce << std::endl;
297  amrex::Print() << " Ce at wall : " << Ce_wall << std::endl;
298  amrex::Print() << " Ck : " << Ck << std::endl;
299  amrex::Print() << " sigma_k : " << sigma_k << std::endl;
300 
301  // Sullivan et al 1994, Eqn 14
302  amrex::Real Cs_equiv = std::sqrt(Ck * std::sqrt(Ck / Ce));
303  amrex::Print() << " equivalent Cs : " << Cs_equiv
304  << std::endl;
305  }
306  amrex::Print() << " isotropic mixing : " << mix_isotropic
307  << std::endl;
308  }
309 
310  if (rans_type != RANSType::None) {
311  if (rans_type == RANSType::kEqn) {
312  amrex::Print() << "Cmu0 : " << Cmu0 << std::endl;
313  amrex::Print() << "sigma_k : " << sigma_k << std::endl;
314  amrex::Print() << "Cb : " << Cb << std::endl;
315  amrex::Print() << "Rt_crit : " << Rt_crit << std::endl;
316  amrex::Print() << "Rt_min : " << Rt_min << std::endl;
317  amrex::Print() << "max_geom_lscale : " << l_g_max << std::endl;
318  }
319  }
320 
321  if ((les_type == LESType::Deardorff) ||
322  (rans_type == RANSType::kEqn)) {
323  if (theta_ref > 0) {
324  amrex::Print() << " reference theta : " << theta_ref << std::endl;
325  } else {
326  amrex::Print() << " reference theta : n/a" << std::endl;
327  }
328  }
329 
330  if ((les_type != LESType::None) || (rans_type != RANSType::None)) {
331  amrex::Print() << " Pr_t : " << Pr_t << std::endl;
332  amrex::Print() << " Sc_t : " << Sc_t << std::endl;
333  }
334 
335  if (pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
336  amrex::Print() << " pbl_mynn_A1 : " << pbl_mynn.A1 << std::endl;
337  amrex::Print() << " pbl_mynn_A2 : " << pbl_mynn.A2 << std::endl;
338  amrex::Print() << " pbl_mynn_B1 : " << pbl_mynn.B1 << std::endl;
339  amrex::Print() << " pbl_mynn_B2 : " << pbl_mynn.B2 << std::endl;
340  amrex::Print() << " pbl_mynn_C1 : " << pbl_mynn.C1 << std::endl;
341  amrex::Print() << " pbl_mynn_C2 : " << pbl_mynn.C2 << std::endl;
342  amrex::Print() << " pbl_mynn_C3 : " << pbl_mynn.C3 << std::endl;
343  amrex::Print() << " pbl_mynn_C4 : " << pbl_mynn.C4 << std::endl;
344  amrex::Print() << " pbl_mynn_C5 : " << pbl_mynn.C5 << std::endl;
345  } else if (pbl_type == PBLType::YSU) {
346  amrex::Print() << " pbl_ysu_coriolis_freq : "
347  << pbl_ysu_coriolis_freq << std::endl;
348  amrex::Print() << " pbl_ysu_use_consistent_coriolis : "
349  << pbl_ysu_use_consistent_coriolis << std::endl;
350  amrex::Print() << " pbl_ysu_force_over_water : "
351  << pbl_ysu_force_over_water << std::endl;
352  amrex::Print() << " pbl_ysu_land_Ribcr : "
353  << pbl_ysu_land_Ribcr << std::endl;
354  amrex::Print() << " pbl_ysu_unst_Ribcr : "
355  << pbl_ysu_unst_Ribcr << std::endl;
356  } else if (pbl_type == PBLType::MRF) {
357  amrex::Print() << " pbl_mrf_coriolis_freq : " << pbl_mrf_coriolis_freq
358  << std::endl;
359  amrex::Print() << " pbl_mrf_Ribcr : " << pbl_mrf_Ribcr
360  << std::endl;
361  amrex::Print() << " pbl_mrf_const_b : " << pbl_mrf_const_b
362  << std::endl;
363  amrex::Print() << " pbl_mrf_sf : " << pbl_mrf_sf
364  << std::endl;
365  amrex::Print() << " mrf_moistvars : " << mrf_moistvars
366  << std::endl;
367  }
368  }
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Real C4
Definition: ERF_MYNNStruct.H:48
amrex::Real C1
Definition: ERF_MYNNStruct.H:45
amrex::Real C3
Definition: ERF_MYNNStruct.H:47
amrex::Real C2
Definition: ERF_MYNNStruct.H:46
amrex::Real A2
Definition: ERF_MYNNStruct.H:42
amrex::Real B1
Definition: ERF_MYNNStruct.H:43
amrex::Real B2
Definition: ERF_MYNNStruct.H:44
amrex::Real C5
Definition: ERF_MYNNStruct.H:49
amrex::Real A1
Definition: ERF_MYNNStruct.H:41
bool smag2d
Definition: ERF_TurbStruct.H:383
amrex::Real sigma_k
Definition: ERF_TurbStruct.H:399
MYNNLevel25 pbl_mynn
Definition: ERF_TurbStruct.H:420
amrex::Real pbl_mrf_const_b
Definition: ERF_TurbStruct.H:456
amrex::Real Rt_min
Definition: ERF_TurbStruct.H:394
amrex::Real pbl_mrf_coriolis_freq
Definition: ERF_TurbStruct.H:453
amrex::Real Ri_crit
Definition: ERF_TurbStruct.H:410
RANSType rans_type
Definition: ERF_TurbStruct.H:413
amrex::Real pbl_mrf_Ribcr
Definition: ERF_TurbStruct.H:454
amrex::Real Ck
Definition: ERF_TurbStruct.H:388
amrex::Real Cmu0
Definition: ERF_TurbStruct.H:391
bool pbl_ysu_use_consistent_coriolis
Definition: ERF_TurbStruct.H:442
amrex::Real Cb
Definition: ERF_TurbStruct.H:392
amrex::Real pbl_mrf_sf
Definition: ERF_TurbStruct.H:458
bool use_Ri_correction
Definition: ERF_TurbStruct.H:409
amrex::Real pbl_ysu_land_Ribcr
Definition: ERF_TurbStruct.H:448
bool mrf_moistvars
Definition: ERF_TurbStruct.H:460
amrex::Real Cs
Definition: ERF_TurbStruct.H:382
amrex::Real Ce
Definition: ERF_TurbStruct.H:386
amrex::Real pbl_ysu_coriolis_freq
Definition: ERF_TurbStruct.H:439
bool pbl_ysu_force_over_water
Definition: ERF_TurbStruct.H:445
amrex::Real Pr_t
Definition: ERF_TurbStruct.H:374
amrex::Real theta_ref
Definition: ERF_TurbStruct.H:402
amrex::Real Sc_t
Definition: ERF_TurbStruct.H:378
amrex::Real Rt_crit
Definition: ERF_TurbStruct.H:393
amrex::Real pbl_ysu_unst_Ribcr
Definition: ERF_TurbStruct.H:450
bool mix_isotropic
Definition: ERF_TurbStruct.H:407
amrex::Real l_g_max
Definition: ERF_TurbStruct.H:395

◆ init_params()

void TurbChoice::init_params ( int  lev,
int  max_level,
std::string  pp_prefix 
)
inline
45  {
46  amrex::ParmParse pp(pp_prefix);
47 
48  // Which LES closure?
49  std::string les_type_string = "None";
50  query_one_or_per_level(pp, "les_type", les_type, lev, max_level);
51 
52  // Handle 2-D Smag
53  if (les_type == LESType::Smagorinsky2D) {
54  les_type = LESType::Smagorinsky;
55  smag2d = true;
56  }
57 
58  // Which RANS closure?
59  std::string rans_type_string = "None";
60  query_one_or_per_level(pp, "rans_type", rans_type, lev, max_level);
61 
62  if ((rans_type != RANSType::None) && (les_type != LESType::None)) {
63  amrex::Error("Hybrid RANS-LES not implemented");
64  }
65 
66  // Which PBL Closure
67  static std::string pbl_type_string = "None";
68  query_one_or_per_level(pp, "pbl_type", pbl_type, lev, max_level);
69 
70  // Do some more stuff for PBL Modeling
71  if (pbl_type != PBLType::None) {
72  // Check for compatibility between PBL, LES, Molec Transport
73  if (les_type != LESType::None) {
74  amrex::Print() << "Selected a PBL model and an LES model: "
75  << "Using PBL for vertical transport, LES for horizontal"
76  << std::endl;
77  }
78  if (les_type == LESType::Smagorinsky) {
79  if (!smag2d)
80  amrex::Error("If using Smagorinsky with a PBL model, the 2-D "
81  "formulation should be used");
82  } else if (les_type == LESType::Deardorff) {
83  amrex::Error(
84  "It is not recommended to use Deardorff LES and a PBL model");
85  }
86 
87  if (pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
88  query_one_or_per_level(pp, "pbl_mynn_A1", pbl_mynn.A1, lev, max_level);
89  query_one_or_per_level(pp, "pbl_mynn_A2", pbl_mynn.A2, lev, max_level);
90  query_one_or_per_level(pp, "pbl_mynn_B1", pbl_mynn.B1, lev, max_level);
91  query_one_or_per_level(pp, "pbl_mynn_B2", pbl_mynn.B2, lev, max_level);
92  query_one_or_per_level(pp, "pbl_mynn_C1", pbl_mynn.C1, lev, max_level);
93  query_one_or_per_level(pp, "pbl_mynn_C2", pbl_mynn.C2, lev, max_level);
94  query_one_or_per_level(pp, "pbl_mynn_C3", pbl_mynn.C3, lev, max_level);
95  query_one_or_per_level(pp, "pbl_mynn_C4", pbl_mynn.C4, lev, max_level);
96  query_one_or_per_level(pp, "pbl_mynn_C5", pbl_mynn.C5, lev, max_level);
101  pp, "pbl_mynn_diffuse_moistvars", pbl_mynn.diffuse_moistvars, lev,
102  max_level);
104  pp, "pbl_mynn_SMmin", pbl_mynn.SMmin, lev, max_level);
106  pp, "pbl_mynn_SMmax", pbl_mynn.SMmax, lev, max_level);
108  pp, "pbl_mynn_SHmin", pbl_mynn.SHmin, lev, max_level);
110  pp, "pbl_mynn_SHmax", pbl_mynn.SHmax, lev, max_level);
112  pp, "pbl_mynn_SQfactor", pbl_mynn.SQfac, lev, max_level);
113  } else if (pbl_type == PBLType::YSU) {
115  pp, "pbl_ysu_coriolis_freq", pbl_ysu_coriolis_freq, lev, max_level);
117  pp, "pbl_ysu_use_consistent_coriolis",
118  pbl_ysu_use_consistent_coriolis, lev, max_level);
120  pp, "pbl_ysu_force_over_water", pbl_ysu_force_over_water, lev,
121  max_level);
123  pp, "pbl_ysu_land_Ribcr", pbl_ysu_land_Ribcr, lev, max_level);
125  pp, "pbl_ysu_unst_Ribcr", pbl_ysu_unst_Ribcr, lev, max_level);
126  } else if (pbl_type == PBLType::MRF) {
128  pp, "pbl_mrf_coriolis_freq", pbl_mrf_coriolis_freq, lev, max_level);
130  pp, "pbl_mrf_Ribcr", pbl_mrf_Ribcr, lev, max_level);
132  pp, "pbl_mrf_const_b", pbl_mrf_const_b, lev, max_level);
133  query_one_or_per_level(pp, "pbl_mrf_sf", pbl_mrf_sf, lev, max_level);
135  pp, "mrf_moistvars", mrf_moistvars, lev, max_level);
136  } else if (pbl_type == PBLType::SHOC) {
137 #ifndef ERF_USE_SHOC
138  amrex::Abort("You set use_shoc to true but didn't build with SHOC; you must rebuild the executable");
139 #endif
140  std::string zlo_bc = "none";
141  amrex::ParmParse pp_bc("zlo");
142  pp_bc.get("type",zlo_bc);
143  if (amrex::toLower(zlo_bc) != "surface_layer") {
144  amrex::Abort("You must use the surface_layer BC at zlo with SHOC.");
145  }
146  }
147  }
148 
149  // Flags for QKE/TKE euation
150  if (pbl_type == PBLType::MYJ || pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
151  // Add sources/sinks to QKE/TKE? (MYJ does this inline)
152  if (pbl_type == PBLType::MYNN25 || pbl_type == PBLType::MYNNEDMF) {
153  use_pbl_tke = true;
154  }
155  // Advect QKE/TKE?
156  query_one_or_per_level(pp, "advect_tke" , advect_tke , lev, max_level);
157  // Apply numerical diffusion to QKE/TKE?
158  query_one_or_per_level(pp, "diffuse_tke_3D", diffuse_tke_3D, lev, max_level);
159  }
160 
161  // LES constants...
162  query_one_or_per_level(pp, "Cs", Cs, lev, max_level);
163 
164  query_one_or_per_level(pp, "Pr_t", Pr_t, lev, max_level);
165  query_one_or_per_level(pp, "Sc_t", Sc_t, lev, max_level);
166 
167  // Compute relevant forms of diffusion parameters
168  Pr_t_inv = amrex::Real(1.0) / Pr_t;
169  Sc_t_inv = amrex::Real(1.0) / Sc_t;
170 
171  if (les_type == LESType::Deardorff) {
172  query_one_or_per_level(pp, "Ck", Ck, lev, max_level);
173  query_one_or_per_level(pp, "Ce", Ce, lev, max_level);
174  query_one_or_per_level(pp, "Ce_wall", Ce_wall, lev, max_level);
175  }
176 
177  // To quantify atmospheric stability for subgrid modeling
178  query_one_or_per_level(pp, "thermal_stratification", strat_type, lev, max_level);
179  if (strat_type == StratType::theta) {
180  amrex::Print() << "Thermal stratification based on gradient of potential temperature" << std::endl;
181  } else if (strat_type == StratType::thetav) {
182  amrex::Print() << "Thermal stratification based on gradient of virtual potential temperature" << std::endl;
183  } else if (strat_type == StratType::thetal) {
184  amrex::Print() << "Thermal stratification based on gradient of linearized liquid-water potential temperature" << std::endl;
185  }
186 
187  // k-eqn constants
188  query_one_or_per_level(pp, "Cmu0", Cmu0, lev, max_level);
189  query_one_or_per_level(pp, "Cb", Cb, lev, max_level);
190  query_one_or_per_level(pp, "Rt_crit", Rt_crit, lev, max_level);
191  query_one_or_per_level(pp, "Rt_min", Rt_min, lev, max_level);
192  query_one_or_per_level(pp, "max_geom_lscale", l_g_max, lev, max_level);
193  query_one_or_per_level(pp, "dirichlet_k", dirichlet_k, lev, max_level);
194 
195  // Common inputs (LES or RANS)
196  if (!query_one_or_per_level(pp, "sigma_k", sigma_k, lev, max_level) && rans_type == RANSType::kEqn) {
197  amrex::Print() << "Overriding default sigma_k for k-eqn RANS" << std::endl;
198  sigma_k = 1.0;
199  };
200  query_one_or_per_level(pp, "theta_ref", theta_ref, lev, max_level);
201 
202  query_one_or_per_level(pp, "mix_isotropic", mix_isotropic, lev, max_level);
203  query_one_or_per_level(pp, "use_Ri_correction", use_Ri_correction, lev, max_level);
204  query_one_or_per_level(pp, "Ri_crit", Ri_crit, lev, max_level);
205 
206  // Set common flags
207  use_kturb =
208  ((les_type != LESType::None) || (rans_type != RANSType::None) ||
209  (pbl_type != PBLType::None));
210  use_keqn =
211  ((les_type == LESType::Deardorff) || (rans_type == RANSType::kEqn));
212  use_tke =
213  ((les_type == LESType::Deardorff) || (rans_type == RANSType::kEqn) ||
214  (pbl_type == PBLType::MYJ) || (pbl_type == PBLType::MYNN25) ||
215  (pbl_type == PBLType::MYNNEDMF) || (pbl_type == PBLType::SHOC));
216 
217  if (use_tke) {
218  query_one_or_per_level(pp, "init_tke_from_ustar", init_tke_from_ustar, lev, max_level);
219  }
220 
221  // Validate inputs
222  if (les_type == LESType::Smagorinsky) {
223  if (Cs == 0) {
224  amrex::Error("Need to specify Cs for Smagorsinky LES");
225  }
226  if (smag2d && mix_isotropic) {
227  amrex::Print() << "Turning off mix_isotropic for 2-D Smagorinsky" << std::endl;
228  mix_isotropic = false;
229  }
230  }
231  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:233
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:56
amrex::Real SHmax
Definition: ERF_MYNNStruct.H:58
amrex::Real SQfac
Definition: ERF_MYNNStruct.H:52
amrex::Real SHmin
Definition: ERF_MYNNStruct.H:57
amrex::Real SMmin
Definition: ERF_MYNNStruct.H:55
bool diffuse_moistvars
Definition: ERF_MYNNStruct.H:63
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:67
bool use_pbl_tke
Definition: ERF_TurbStruct.H:428
amrex::Real Sc_t_inv
Definition: ERF_TurbStruct.H:379
bool use_keqn
Definition: ERF_TurbStruct.H:425
MYNNLevel2 pbl_mynn_level2
Definition: ERF_TurbStruct.H:421
StratType strat_type
Definition: ERF_TurbStruct.H:404
bool dirichlet_k
Definition: ERF_TurbStruct.H:415
bool advect_tke
Definition: ERF_TurbStruct.H:463
bool init_tke_from_ustar
Definition: ERF_TurbStruct.H:435
bool use_tke
Definition: ERF_TurbStruct.H:430
bool diffuse_tke_3D
Definition: ERF_TurbStruct.H:465
amrex::Real Pr_t_inv
Definition: ERF_TurbStruct.H:375
bool use_kturb
Definition: ERF_TurbStruct.H:424
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 = 0.35

◆ Ce

amrex::Real TurbChoice::Ce = 0.93

◆ Ce_wall

amrex::Real TurbChoice::Ce_wall = 0.0

◆ Ck

amrex::Real TurbChoice::Ck = 0.1

◆ Cmu0

◆ Cs

amrex::Real TurbChoice::Cs = 0.0

◆ 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 = 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
Initial value:
=
7.8

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

◆ pbl_mrf_coriolis_freq

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

Referenced by display(), and init_params().

◆ pbl_mrf_Ribcr

amrex::Real TurbChoice::pbl_mrf_Ribcr
Initial value:
=
0.5

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

◆ pbl_mrf_sf

amrex::Real TurbChoice::pbl_mrf_sf
Initial value:
=
0.1

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

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

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

◆ pbl_ysu_unst_Ribcr

amrex::Real TurbChoice::pbl_ysu_unst_Ribcr
Initial value:
=
0.0

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 = amrex::Real(1.0) / amrex::Real(3.0)

Referenced by display(), and init_params().

◆ Pr_t_inv

◆ rans_type

◆ Ri_crit

amrex::Real TurbChoice::Ri_crit = 0.25

◆ Rt_crit

amrex::Real TurbChoice::Rt_crit = -1.0

◆ Rt_min

amrex::Real TurbChoice::Rt_min = -3.0

◆ Sc_t

amrex::Real TurbChoice::Sc_t = 1.0

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: