ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_TurbStruct.H
Go to the documentation of this file.
1 #ifndef ERF_TURB_STRUCT_H_
2 #define ERF_TURB_STRUCT_H_
3 
4 #include <ERF_MYNNStruct.H>
5 
6 AMREX_ENUM(LESType, None, Smagorinsky, Smagorinsky2D, Deardorff);
7 
8 AMREX_ENUM(RANSType, None, kEqn);
9 
10 AMREX_ENUM(PBLType, None, MYJ, MYNN25, MYNNEDMF, YSU, MRF, SHOC);
11 
12 AMREX_ENUM(StratType, theta, thetav, thetal);
13 
14 template <typename T>
15 int
17  const amrex::ParmParse& pp,
18  const char* query_string,
19  T& query_var,
20  const int lev,
21  const int maxlev)
22 {
23  int count = pp.countval(query_string);
24  if (count == 0) {
25  return 0; // nothing to do
26  } else if (count == 1) {
27  // In this case, we assume the same value is used at every level
28  return pp.query(query_string, query_var);
29  } else if (count >= maxlev + 1) {
30  // In this case, there may be more values than levels so we
31  // will just read the number of values we need and ignore the
32  // extra levels
33  return pp.query(query_string, query_var, lev);
34  } else {
35  // In this case, there is more than one value AND there
36  // are more levels than values; as an example, if max_level = 2
37  // but count = 2 so we have three levels but only two values --
38  // it is not clear how to interpret this so we abort
39  amrex::Error(
40  "For parmparse variable " + pp.prefixedName(query_string) +
41  ": if specified, specify once total or at least once for each level");
42  return 0; // avoid compiler warning
43  }
44 }
45 
46 template <typename T>
47 int
49  const amrex::ParmParse& pp,
50  const char* query_string,
51  T& query_var,
52  const int lev,
53  const int maxlev)
54 {
55  int count = pp.countval(query_string);
56  if (count == 0) {
57  return 0; // nothing to do
58  } else if (count == 1) {
59  // In this case, we assume the same value is used at every level
60  return pp.query_enum_case_insensitive(query_string, query_var);
61  } else if (count >= maxlev + 1) {
62  // In this case, there may be more values than levels so we
63  // will just read the number of values we need and ignore the
64  // extra levels
65  return pp.query_enum_case_insensitive(query_string, query_var, lev);
66  } else {
67  // In this case, there is more than one value AND there
68  // are more levels than values; as an example, if max_level = 2
69  // but count = 2 so we have three levels but only two values --
70  // it is not clear how to interpret this so we abort
71  amrex::Error(
72  "For parmparse variable " + pp.prefixedName(query_string) +
73  ": if specified, specify once total or at least once for each level");
74  return 0; // avoid compiler warning
75  }
76 }
77 
78 /**
79  * Container holding quantities related to turbulence parametrizations
80  */
81 struct TurbChoice
82 {
83 public:
84  void init_params (int lev, int max_level, std::string pp_prefix)
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  }
273 
274  void check_params (amrex::GpuArray<ERF_BC, AMREX_SPACEDIM*2>& phys_bc_type)
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  }
293 
294  void display (int lev)
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  }
410 
411  // LES model
412  LESType les_type = LESType::None;
413 
414  // Turbulent Prandtl number
417 
418  // Turbulent Schmidt number
421 
422  // Smagorinsky
424  bool smag2d = false;
425 
426  // Deardorff
428  amrex::Real Ce_wall = zero; // if > 0, then set Ce to this at k=0
430 
431  // k-eqn RANS coefficients (Axell & Liungman 2001)
436  amrex::Real l_g_max = amrex::Real(30.0); // ~ kappa * (amrex::Real(0.1) * zi)
437 
438  // Deardorff or k-eqn RANS
439  // - diffusivity of tke is 1/sigma_k times the eddy viscosity
441  // - reference potential temperature used to quantify the stratification in
442  // a stable region
444  // - how to quantify stratification effects
446 
447  // Anisotropic length scales
448  bool mix_isotropic = true;
449 
450  bool use_Ri_correction = true;
452 
453  // RANS type
454  RANSType rans_type = RANSType::None;
455 
456  bool dirichlet_k = false;
457 
458  // PBL model
459  PBLType pbl_type = PBLType::None;
460 
462  MYNNLevel2 pbl_mynn_level2; // limiting in the decaying turbulence regime
463 
464  // Common Flags
465  bool use_kturb = false; // Any turbulence modeling?
466  bool use_keqn =
467  false; // Any microscale turbulence modeling (LES, RANS) with TKE closure?
468  // Then need to populate SmnSmn_lev for production term.
469  bool use_pbl_tke =
470  false; // Any mesoscale turbulence modeling (PBL) with TKE closure?
471  bool use_tke = false; // Any TKE closure (meso or microscale)?
472 
473  // Initialize TKE/QKE with linear profiles whose surface values are a
474  // function of the friction velocity calculated by the surface layer scheme
475  // (e.g., as done in MYNN-EDMF)
476  bool init_tke_from_ustar = false;
477 
478  // This is the value of tke_min in WRF 4.5
480 
481  // Model coefficients - YSU
482  // TODO: Add parmparse for all of these above
484  amrex::Real(1.0e-4); // 1e-4 is hardcoded in WRF, we let the user specify or take the
485  // value from ERF coriolis forcing
487  false; // ignore input pbl_ysu_coriolis_freq, take value from ERF coriolis
488  // forcing instead
490  false; // Force YSU to act as if it is over water regardless of other inputs
491  // (for testing)
493  fourth; // Critical Bulk Richardson number of Land for stable conditions
495  zero; // Critical Bulk Richardson number for unstable conditions
496  // Model coefficients - MRF
498  amrex::Real pbl_mrf_Ribcr = amrex::Real(0.5); // Critical Bulk Richardson number for MRF PBL model
499  amrex::Real pbl_mrf_const_b = amrex::Real(7.8); // Constant b in MRF PBL model, used to compute the PBL height
500  amrex::Real pbl_mrf_sf = amrex::Real(0.1); // MRF surface flux, used to compute the PBL height
501  bool mrf_moistvars = false; // if true, adds turbulence to moisture
502  // MYNN2.5 PBL model
503  // TKE/QKE stuff
504  bool advect_tke = true; // if MYNN2.5 PBL is used default is turb transport in
505  // Z-direction only
506  bool diffuse_tke_3D = true; // if numerical diffusion is turned on
507 };
508 #endif
constexpr amrex::Real three
Definition: ERF_Constants.H:9
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
@ no_slip_wall
@ surface_layer
ParmParse pp("prob")
Real T
Definition: ERF_InitCustomPert_Bubble.H:105
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_ENUM(LESType, None, Smagorinsky, Smagorinsky2D, Deardorff)
int query_one_or_per_level_enum_case_insensitive(const amrex::ParmParse &pp, const char *query_string, T &query_var, const int lev, const int maxlev)
Definition: ERF_TurbStruct.H:48
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
Definition: ERF_MYNNStruct.H:11
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 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 SHmin
Definition: ERF_MYNNStruct.H:59
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 SMmin
Definition: ERF_MYNNStruct.H:57
amrex::Real A1
Definition: ERF_MYNNStruct.H:43
bool diffuse_moistvars
Definition: ERF_MYNNStruct.H:65
Definition: ERF_MYNNStruct.H:68
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
Definition: ERF_TurbStruct.H:82
amrex::Real tke_min
Definition: ERF_TurbStruct.H:479
bool smag2d
Definition: ERF_TurbStruct.H:424
amrex::Real sigma_k
Definition: ERF_TurbStruct.H:440
MYNNLevel25 pbl_mynn
Definition: ERF_TurbStruct.H:461
PBLType pbl_type
Definition: ERF_TurbStruct.H:459
amrex::Real pbl_mrf_const_b
Definition: ERF_TurbStruct.H:499
bool use_pbl_tke
Definition: ERF_TurbStruct.H:469
void check_params(amrex::GpuArray< ERF_BC, AMREX_SPACEDIM *2 > &phys_bc_type)
Definition: ERF_TurbStruct.H:274
amrex::Real Sc_t_inv
Definition: ERF_TurbStruct.H:420
amrex::Real Rt_min
Definition: ERF_TurbStruct.H:435
bool use_keqn
Definition: ERF_TurbStruct.H:466
amrex::Real pbl_mrf_coriolis_freq
Definition: ERF_TurbStruct.H:497
MYNNLevel2 pbl_mynn_level2
Definition: ERF_TurbStruct.H:462
amrex::Real Ri_crit
Definition: ERF_TurbStruct.H:451
RANSType rans_type
Definition: ERF_TurbStruct.H:454
StratType strat_type
Definition: ERF_TurbStruct.H:445
bool dirichlet_k
Definition: ERF_TurbStruct.H:456
amrex::Real pbl_mrf_Ribcr
Definition: ERF_TurbStruct.H:498
bool advect_tke
Definition: ERF_TurbStruct.H:504
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
void init_params(int lev, int max_level, std::string pp_prefix)
Definition: ERF_TurbStruct.H:84
bool use_Ri_correction
Definition: ERF_TurbStruct.H:450
bool init_tke_from_ustar
Definition: ERF_TurbStruct.H:476
amrex::Real pbl_ysu_land_Ribcr
Definition: ERF_TurbStruct.H:492
bool mrf_moistvars
Definition: ERF_TurbStruct.H:501
bool use_tke
Definition: ERF_TurbStruct.H:471
bool diffuse_tke_3D
Definition: ERF_TurbStruct.H:506
amrex::Real Cs
Definition: ERF_TurbStruct.H:423
amrex::Real Pr_t_inv
Definition: ERF_TurbStruct.H:416
amrex::Real Ce_wall
Definition: ERF_TurbStruct.H:428
amrex::Real Ce
Definition: ERF_TurbStruct.H:427
LESType les_type
Definition: ERF_TurbStruct.H:412
amrex::Real pbl_ysu_coriolis_freq
Definition: ERF_TurbStruct.H:483
bool pbl_ysu_force_over_water
Definition: ERF_TurbStruct.H:489
void display(int lev)
Definition: ERF_TurbStruct.H:294
amrex::Real Pr_t
Definition: ERF_TurbStruct.H:415
bool use_kturb
Definition: ERF_TurbStruct.H:465
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