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, Deardorff);
7 
8 AMREX_ENUM(RANSType, None, kEqn);
9 
10 AMREX_ENUM(PBLType, None, MYNN25, YSU);
11 
12 template <typename T>
13 void query_one_or_per_level(const amrex::ParmParse& pp, const char* query_string, T& query_var, const int lev, const int maxlev)
14 {
15  int count = pp.countval(query_string);
16  if (count == 0) {
17  return; // nothing to do
18  } else if (count == 1) {
19  pp.query(query_string, query_var);
20  } else if (count == maxlev + 1) {
21  pp.query(query_string, query_var, lev);
22  } else {
23  amrex::Error("For parmparse variable " + pp.prefixedName(query_string)
24  + ": if specified, specify once total or once for each level");
25  }
26 }
27 
28 /**
29  * Container holding quantities related to turbulence parametrizations
30  */
31 struct TurbChoice {
32  public:
33  void init_params(int lev, int max_level)
34  {
35  amrex::ParmParse pp(pp_prefix);
36 
37  // Which LES closure?
38  std::string les_type_string = "None";
39  query_one_or_per_level(pp, "les_type", les_type, lev, max_level);
40 
41  // Which RANS closure?
42  std::string rans_type_string = "None";
43  query_one_or_per_level(pp, "rans_type", rans_type, lev, max_level);
44 
45  if ((rans_type != RANSType::None) && (les_type != LESType::None)) {
46  amrex::Error("Hybrid RANS-LES not implemented");
47  }
48 
49  // Which PBL Closure
50  static std::string pbl_type_string = "None";
51  query_one_or_per_level(pp, "pbl_type", pbl_type, lev, max_level);
52 
53  // Do some more stuff for PBL Modeling
54  if (pbl_type != PBLType::None) {
55  // Check for compatibility between PBL, LES, Molec Transport
56  if (les_type != LESType::None) {
57  amrex::Print() << "Selected a PBL model and an LES model: " <<
58  "Using PBL for vertical transport, LES for horizontal" << std::endl;
59  } else if (les_type == LESType::Deardorff) {
60  amrex::Error("It is not recommended to use Deardorff LES and a PBL model");
61  }
62 
63  if (pbl_type == PBLType::MYNN25) {
64  query_one_or_per_level(pp, "pbl_mynn_A1", pbl_mynn.A1, lev, max_level);
65  query_one_or_per_level(pp, "pbl_mynn_A2", pbl_mynn.A2, lev, max_level);
66  query_one_or_per_level(pp, "pbl_mynn_B1", pbl_mynn.B1, lev, max_level);
67  query_one_or_per_level(pp, "pbl_mynn_B2", pbl_mynn.B2, lev, max_level);
68  query_one_or_per_level(pp, "pbl_mynn_C1", pbl_mynn.C1, lev, max_level);
69  query_one_or_per_level(pp, "pbl_mynn_C2", pbl_mynn.C2, lev, max_level);
70  query_one_or_per_level(pp, "pbl_mynn_C3", pbl_mynn.C3, lev, max_level);
71  query_one_or_per_level(pp, "pbl_mynn_C4", pbl_mynn.C4, lev, max_level);
72  query_one_or_per_level(pp, "pbl_mynn_C5", pbl_mynn.C5, lev, max_level);
75  query_one_or_per_level(pp, "pbl_mynn_diffuse_moistvars", pbl_mynn.diffuse_moistvars, lev, max_level);
76  query_one_or_per_level(pp, "pbl_mynn_SMmin", pbl_mynn.SMmin, lev, max_level);
77  query_one_or_per_level(pp, "pbl_mynn_SMmax", pbl_mynn.SMmax, lev, max_level);
78  query_one_or_per_level(pp, "pbl_mynn_SHmin", pbl_mynn.SHmin, lev, max_level);
79  query_one_or_per_level(pp, "pbl_mynn_SHmax", pbl_mynn.SHmax, lev, max_level);
80  } else if (pbl_type == PBLType::YSU) {
81  query_one_or_per_level(pp, "pbl_ysu_coriolis_freq", pbl_ysu_coriolis_freq, lev, max_level);
82  query_one_or_per_level(pp, "pbl_ysu_use_consistent_coriolis", pbl_ysu_use_consistent_coriolis, lev, max_level);
83  query_one_or_per_level(pp, "pbl_ysu_force_over_water", pbl_ysu_force_over_water, lev, max_level);
84  query_one_or_per_level(pp, "pbl_ysu_land_Ribcr", pbl_ysu_land_Ribcr, lev, max_level);
85  query_one_or_per_level(pp, "pbl_ysu_unst_Ribcr", pbl_ysu_unst_Ribcr, lev, max_level);
86  }
87  }
88 
89  // Right now, solving the QKE equation is only supported when MYNN PBL is turned on
90  if (pbl_type == PBLType::MYNN25) {
91  use_KE = true;
92  query_one_or_per_level(pp, "advect_KE" , advect_KE, lev, max_level);
93  query_one_or_per_level(pp, "diffuse_KE_3D", diffuse_KE_3D, lev, max_level);
94  }
95 
96  // LES constants...
97  query_one_or_per_level(pp, "Cs" ,Cs, lev, max_level);
98  query_one_or_per_level(pp, "Pr_t",Pr_t, lev, max_level);
99  query_one_or_per_level(pp, "Sc_t",Sc_t, lev, max_level);
100 
101  // Compute relevant forms of diffusion parameters
102  Pr_t_inv = amrex::Real(1.0) / Pr_t;
103  Sc_t_inv = amrex::Real(1.0) / Sc_t;
104 
105  query_one_or_per_level(pp, "Ce" , Ce, lev, max_level);
106  query_one_or_per_level(pp, "Ce_wall" , Ce_wall, lev, max_level);
107 
108  if (les_type == LESType::Deardorff) {
109  query_one_or_per_level(pp, "Ck" , Ck, lev, max_level);
110  }
111 
112  // k-eqn constants
113  query_one_or_per_level(pp, "Cmu0" , Cmu0 , lev, max_level);
114  query_one_or_per_level(pp, "Cb" , Cb , lev, max_level);
115  query_one_or_per_level(pp, "Rt_crit" , Rt_crit, lev, max_level);
116  query_one_or_per_level(pp, "Rt_min" , Rt_min , lev, max_level);
117  query_one_or_per_level(pp, "max_geom_lscale", l_g_max, lev, max_level);
118 
119  // Common inputs (LES or RANS)
120  query_one_or_per_level(pp, "sigma_k" , sigma_k, lev, max_level);
121  query_one_or_per_level(pp, "theta_ref", theta_ref, lev, max_level);
122 
123  // Validate inputs
124  if (les_type == LESType::Smagorinsky) {
125  if (Cs == 0) {
126  amrex::Error("Need to specify Cs for Smagorsinky LES");
127  }
128  }
129  }
130 
131  void display(int lev)
132  {
133  amrex::Print() << "Turbulence Settings at level " << lev << std::endl;
134 
135  if (les_type == LESType::None && rans_type == RANSType::None && pbl_type == PBLType::None) {
136  amrex::Print() << " Using DNS model at level " << lev << std::endl;
137  } else if (les_type == LESType::Smagorinsky) {
138  amrex::Print() << " Using Smagorinsky LES model at level " << lev << std::endl;
139  } else if (les_type == LESType::Deardorff) {
140  amrex::Print() << " Using Deardorff LES model at level " << lev << std::endl;
141  } else if (rans_type == RANSType::kEqn) {
142  amrex::Print() << " Using Axell & Liungman one-equation RANS k model at level " << lev << std::endl;
143  } else if (pbl_type == PBLType::MYNN25) {
144  amrex::Print() << " Using MYNN2.5 PBL model at level " << lev << std::endl;
145  } else if (pbl_type == PBLType::YSU) {
146  amrex::Print() << " Using YSU PBL model at level " << lev << std::endl;
147  }
148 
149  if (les_type != LESType::None) {
150  if (les_type == LESType::Smagorinsky) {
151  amrex::Print() << "Cs : " << Cs << std::endl;
152  }
153  if (les_type == LESType::Deardorff) {
154  amrex::Print() << "Ce : " << Ce << std::endl;
155  amrex::Print() << "Ce at wall : " << Ce_wall << std::endl;
156  amrex::Print() << "Ck : " << Ck << std::endl;
157  amrex::Print() << "sigma_k : " << sigma_k << std::endl;
158  amrex::Print() << "reference theta : " << theta_ref << std::endl;
159  // Sullivan et al 1994, Eqn 14
160  amrex::Real Cs_equiv = std::sqrt(Ck * std::sqrt(Ck/Ce));
161  amrex::Print() << "equivalent Cs : " << Cs_equiv << std::endl;
162  }
163  }
164 
165  if (rans_type != RANSType::None) {
166  if (rans_type == RANSType::kEqn) {
167  amrex::Print() << "Cmu0 : " << Cmu0 << std::endl;
168  amrex::Print() << "sigma_k : " << sigma_k << std::endl;
169  amrex::Print() << "Cb : " << Cb << std::endl;
170  amrex::Print() << "Rt_crit : " << Rt_crit << std::endl;
171  amrex::Print() << "Rt_min : " << Rt_min << std::endl;
172  amrex::Print() << "max_geom_lscale : " << l_g_max << std::endl;
173  amrex::Print() << "reference theta : " << theta_ref << std::endl;
174  }
175  }
176 
177  if ((les_type != LESType::None) || (rans_type != RANSType::None)) {
178  amrex::Print() << "Pr_t : " << Pr_t << std::endl;
179  amrex::Print() << "Sc_t : " << Sc_t << std::endl;
180  }
181 
182  if (pbl_type == PBLType::MYNN25) {
183  amrex::Print() << " pbl_mynn_A1 : " << pbl_mynn.A1 << std::endl;
184  amrex::Print() << " pbl_mynn_A2 : " << pbl_mynn.A2 << std::endl;
185  amrex::Print() << " pbl_mynn_B1 : " << pbl_mynn.B1 << std::endl;
186  amrex::Print() << " pbl_mynn_B2 : " << pbl_mynn.B2 << std::endl;
187  amrex::Print() << " pbl_mynn_C1 : " << pbl_mynn.C1 << std::endl;
188  amrex::Print() << " pbl_mynn_C2 : " << pbl_mynn.C2 << std::endl;
189  amrex::Print() << " pbl_mynn_C3 : " << pbl_mynn.C3 << std::endl;
190  amrex::Print() << " pbl_mynn_C4 : " << pbl_mynn.C4 << std::endl;
191  amrex::Print() << " pbl_mynn_C5 : " << pbl_mynn.C5 << std::endl;
192  } else if (pbl_type == PBLType::YSU) {
193  amrex::Print() << " pbl_ysu_coriolis_freq : " << pbl_ysu_coriolis_freq << std::endl;
194  amrex::Print() << " pbl_ysu_use_consistent_coriolis : " << pbl_ysu_use_consistent_coriolis << std::endl;
195  amrex::Print() << " pbl_ysu_force_over_water : " << pbl_ysu_force_over_water << std::endl;
196  amrex::Print() << " pbl_ysu_land_Ribcr : " << pbl_ysu_land_Ribcr << std::endl;
197  amrex::Print() << " pbl_ysu_unst_Ribcr : " << pbl_ysu_unst_Ribcr << std::endl;
198  }
199  }
200 
201  // Default prefix
202  std::string pp_prefix {"erf"};
203 
204  // LES model
205  LESType les_type;
206 
207  // Turbulent Prandtl number
208  amrex::Real Pr_t = amrex::Real(1.0) / amrex::Real(3.0);
209  amrex::Real Pr_t_inv = amrex::Real(3.0);
210 
211  // Turbulent Schmidt number
212  amrex::Real Sc_t = 1.0;
213  amrex::Real Sc_t_inv = 1.0;
214 
215  // Smagorinsky Cs coefficient
216  amrex::Real Cs = 0.0;
217 
218  // Deardorff
219  amrex::Real Ce = 0.93;
220  amrex::Real Ce_wall = 0.0; // if > 0, then set Ce to this at k=0
221  amrex::Real Ck = 0.1;
222 
223  // k-eqn RANS coefficients (Axell & Liungman 2001)
224  amrex::Real Cmu0 = 0.5562;
225  amrex::Real Cb = 0.35;
226  amrex::Real Rt_crit = -1.0;
227  amrex::Real Rt_min = -3.0;
228  amrex::Real l_g_max = 1e99;
229 
230  // Deardorff or k-eqn RANS
231  // - diffusivity of tke is 1/sigma_k times the eddy viscosity
232  amrex::Real sigma_k = 0.5;
233  // - reference potential temperature used to quantify the stratification in
234  // a stable region
235  amrex::Real theta_ref = 300.0;
236 
237  // RANS type
238  RANSType rans_type;
239 
240  // PBL model
241  PBLType pbl_type;
242 
244  MYNNLevel2 pbl_mynn_level2; // for filtering
245 
246  // Model coefficients - YSU
247  // TODO: Add parmparse for all of these above
248  amrex::Real pbl_ysu_coriolis_freq = 1.0e-4; // 1e-4 is hardcoded in WRF, we let the user specify or tske the value from ERF coriolis forcing
249  bool pbl_ysu_use_consistent_coriolis = false; // ignore input pbl_ysu_coriolis_freq, take value from ERF coriolis forcing instead
250  bool pbl_ysu_force_over_water = false; // Force YSU to act as if it is over water regardless of other inputs (for testing)
251  amrex::Real pbl_ysu_land_Ribcr = 0.25; // Critical Bulk Richardson number of Land for stable conditions
252  amrex::Real pbl_ysu_unst_Ribcr = 0.0; // Critical Bulk Richardson number for unstable conditions
253 
254  // QKE stuff - default is to use it, if MYNN2.5 PBL is used default is turb transport in Z-direction only
255  bool use_KE = true;
256  bool diffuse_KE_3D = true;
257  bool advect_KE = true;
258 };
259 #endif
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:219
AMREX_ENUM(LESType, None, Smagorinsky, Deardorff)
void 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:13
@ T
Definition: ERF_IndexDefines.H:99
Definition: ERF_MYNNStruct.H:9
amrex::Real SMmax
Definition: ERF_MYNNStruct.H:53
amrex::Real SHmax
Definition: ERF_MYNNStruct.H:55
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 SHmin
Definition: ERF_MYNNStruct.H:54
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 SMmin
Definition: ERF_MYNNStruct.H:52
amrex::Real A1
Definition: ERF_MYNNStruct.H:41
bool diffuse_moistvars
Definition: ERF_MYNNStruct.H:60
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:64
Definition: ERF_TurbStruct.H:31
amrex::Real sigma_k
Definition: ERF_TurbStruct.H:232
MYNNLevel25 pbl_mynn
Definition: ERF_TurbStruct.H:243
PBLType pbl_type
Definition: ERF_TurbStruct.H:241
bool use_KE
Definition: ERF_TurbStruct.H:255
amrex::Real Sc_t_inv
Definition: ERF_TurbStruct.H:213
amrex::Real Rt_min
Definition: ERF_TurbStruct.H:227
void init_params(int lev, int max_level)
Definition: ERF_TurbStruct.H:33
MYNNLevel2 pbl_mynn_level2
Definition: ERF_TurbStruct.H:244
RANSType rans_type
Definition: ERF_TurbStruct.H:238
bool diffuse_KE_3D
Definition: ERF_TurbStruct.H:256
amrex::Real Ck
Definition: ERF_TurbStruct.H:221
amrex::Real Cmu0
Definition: ERF_TurbStruct.H:224
std::string pp_prefix
Definition: ERF_TurbStruct.H:202
bool pbl_ysu_use_consistent_coriolis
Definition: ERF_TurbStruct.H:249
amrex::Real Cb
Definition: ERF_TurbStruct.H:225
amrex::Real pbl_ysu_land_Ribcr
Definition: ERF_TurbStruct.H:251
amrex::Real Cs
Definition: ERF_TurbStruct.H:216
amrex::Real Pr_t_inv
Definition: ERF_TurbStruct.H:209
amrex::Real Ce_wall
Definition: ERF_TurbStruct.H:220
amrex::Real Ce
Definition: ERF_TurbStruct.H:219
LESType les_type
Definition: ERF_TurbStruct.H:205
amrex::Real pbl_ysu_coriolis_freq
Definition: ERF_TurbStruct.H:248
bool pbl_ysu_force_over_water
Definition: ERF_TurbStruct.H:250
void display(int lev)
Definition: ERF_TurbStruct.H:131
amrex::Real Pr_t
Definition: ERF_TurbStruct.H:208
amrex::Real theta_ref
Definition: ERF_TurbStruct.H:235
amrex::Real Sc_t
Definition: ERF_TurbStruct.H:212
bool advect_KE
Definition: ERF_TurbStruct.H:257
amrex::Real Rt_crit
Definition: ERF_TurbStruct.H:226
amrex::Real pbl_ysu_unst_Ribcr
Definition: ERF_TurbStruct.H:252
amrex::Real l_g_max
Definition: ERF_TurbStruct.H:228