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  // Which PBL Closure
46  static std::string pbl_type_string = "None";
47  query_one_or_per_level(pp, "pbl_type", pbl_type, lev, max_level);
48 
49  // Do some more stuff for PBL Modeling
50  if (pbl_type != PBLType::None) {
51  // Check for compatibility between PBL, LES, Molec Transport
52  if (les_type != LESType::None) {
53  amrex::Print() << "Selected a PBL model and an LES model: " <<
54  "Using PBL for vertical transport, LES for horizontal" << std::endl;
55  } else if (les_type == LESType::Deardorff) {
56  amrex::Error("It is not recommended to use Deardorff LES and a PBL model");
57  }
58 
59  if (pbl_type == PBLType::MYNN25) {
60  query_one_or_per_level(pp, "pbl_mynn_A1", pbl_mynn.A1, lev, max_level);
61  query_one_or_per_level(pp, "pbl_mynn_A2", pbl_mynn.A2, lev, max_level);
62  query_one_or_per_level(pp, "pbl_mynn_B1", pbl_mynn.B1, lev, max_level);
63  query_one_or_per_level(pp, "pbl_mynn_B2", pbl_mynn.B2, lev, max_level);
64  query_one_or_per_level(pp, "pbl_mynn_C1", pbl_mynn.C1, lev, max_level);
65  query_one_or_per_level(pp, "pbl_mynn_C2", pbl_mynn.C2, lev, max_level);
66  query_one_or_per_level(pp, "pbl_mynn_C3", pbl_mynn.C3, lev, max_level);
67  query_one_or_per_level(pp, "pbl_mynn_C4", pbl_mynn.C4, lev, max_level);
68  query_one_or_per_level(pp, "pbl_mynn_C5", pbl_mynn.C5, lev, max_level);
71  query_one_or_per_level(pp, "pbl_mynn_diffuse_moistvars", pbl_mynn.diffuse_moistvars, lev, max_level);
72  query_one_or_per_level(pp, "pbl_mynn_SMmin", pbl_mynn.SMmin, lev, max_level);
73  query_one_or_per_level(pp, "pbl_mynn_SMmax", pbl_mynn.SMmax, lev, max_level);
74  query_one_or_per_level(pp, "pbl_mynn_SHmin", pbl_mynn.SHmin, lev, max_level);
75  query_one_or_per_level(pp, "pbl_mynn_SHmax", pbl_mynn.SHmax, lev, max_level);
76  } else if (pbl_type == PBLType::YSU) {
77  query_one_or_per_level(pp, "pbl_ysu_coriolis_freq", pbl_ysu_coriolis_freq, lev, max_level);
78  query_one_or_per_level(pp, "pbl_ysu_use_consistent_coriolis", pbl_ysu_use_consistent_coriolis, lev, max_level);
79  query_one_or_per_level(pp, "pbl_ysu_force_over_water", pbl_ysu_force_over_water, lev, max_level);
80  query_one_or_per_level(pp, "pbl_ysu_land_Ribcr", pbl_ysu_land_Ribcr, lev, max_level);
81  query_one_or_per_level(pp, "pbl_ysu_unst_Ribcr", pbl_ysu_unst_Ribcr, lev, max_level);
82  }
83  }
84 
85  // Right now, solving the QKE equation is only supported when MYNN PBL is turned on
86  if (pbl_type == PBLType::MYNN25) {
87  use_KE = true;
88  query_one_or_per_level(pp, "advect_KE" , advect_KE, lev, max_level);
89  query_one_or_per_level(pp, "diffuse_KE_3D", diffuse_KE_3D, lev, max_level);
90  }
91 
92  // LES constants...
93  query_one_or_per_level(pp, "Cs" ,Cs, lev, max_level);
94  query_one_or_per_level(pp, "CI" ,CI, lev, max_level);
95  query_one_or_per_level(pp, "Pr_t",Pr_t, lev, max_level);
96  query_one_or_per_level(pp, "Sc_t",Sc_t, lev, max_level);
97 
98  // Compute relevant forms of diffusion parameters
99  Pr_t_inv = amrex::Real(1.0) / Pr_t;
100  Sc_t_inv = amrex::Real(1.0) / Sc_t;
101 
102  query_one_or_per_level(pp, "Ce" , Ce, lev, max_level);
103  query_one_or_per_level(pp, "Ce_wall" , Ce_wall, lev, max_level);
104  query_one_or_per_level(pp, "sigma_k" , sigma_k, lev, max_level);
105 
106  if (les_type == LESType::Deardorff) {
107  query_one_or_per_level(pp, "Ck" , Ck, lev, max_level);
108  }
109 
110  query_one_or_per_level(pp, "theta_ref", theta_ref, lev, max_level);
111 
112  // Validate inputs
113  if (les_type == LESType::Smagorinsky) {
114  if (Cs == 0) {
115  amrex::Error("Need to specify Cs for Smagorsinky LES");
116  }
117  }
118  }
119 
120  void display(int lev)
121  {
122  amrex::Print() << "Turbulence Settings at level " << lev << std::endl;
123 
124  if (les_type == LESType::None && pbl_type == PBLType::None) {
125  amrex::Print() << " Using DNS model at level " << lev << std::endl;
126  } else if (les_type == LESType::Smagorinsky) {
127  amrex::Print() << " Using Smagorinsky LES model at level " << lev << std::endl;
128  } else if (les_type == LESType::Deardorff) {
129  amrex::Print() << " Using Deardorff LES model at level " << lev << std::endl;
130  } else if (pbl_type == PBLType::MYNN25) {
131  amrex::Print() << " Using MYNN2.5 PBL model at level " << lev << std::endl;
132  } else if (pbl_type == PBLType::YSU) {
133  amrex::Print() << " Using YSU PBL model at level " << lev << std::endl;
134  }
135 
136  if (les_type != LESType::None) {
137  if (les_type == LESType::Smagorinsky) {
138  amrex::Print() << "Cs : " << Cs << std::endl;
139  amrex::Print() << "CI : " << CI << std::endl;
140  amrex::Print() << "Pr_t : " << Pr_t << std::endl;
141  amrex::Print() << "Sc_t : " << Sc_t << std::endl;
142  }
143  if (les_type == LESType::Deardorff) {
144  amrex::Print() << "Ce : " << Ce << std::endl;
145  amrex::Print() << "Ce at wall : " << Ce_wall << std::endl;
146  amrex::Print() << "Ck : " << Ck << std::endl;
147  amrex::Print() << "sigma_k : " << sigma_k << std::endl;
148  }
149  amrex::Print() << "reference theta : " << theta_ref << std::endl;
150  }
151 
152  if (pbl_type == PBLType::MYNN25) {
153  amrex::Print() << "pbl_mynn_A1 : " << pbl_mynn.A1 << std::endl;
154  amrex::Print() << "pbl_mynn_A2 : " << pbl_mynn.A2 << std::endl;
155  amrex::Print() << "pbl_mynn_B1 : " << pbl_mynn.B1 << std::endl;
156  amrex::Print() << "pbl_mynn_B2 : " << pbl_mynn.B2 << std::endl;
157  amrex::Print() << "pbl_mynn_C1 : " << pbl_mynn.C1 << std::endl;
158  amrex::Print() << "pbl_mynn_C2 : " << pbl_mynn.C2 << std::endl;
159  amrex::Print() << "pbl_mynn_C3 : " << pbl_mynn.C3 << std::endl;
160  amrex::Print() << "pbl_mynn_C4 : " << pbl_mynn.C4 << std::endl;
161  amrex::Print() << "pbl_mynn_C5 : " << pbl_mynn.C5 << std::endl;
162  } else if (pbl_type == PBLType::YSU) {
163  amrex::Print() << "pbl_ysu_coriolis_freq : " << pbl_ysu_coriolis_freq << std::endl;
164  amrex::Print() << "pbl_ysu_use_consistent_coriolis : " << pbl_ysu_use_consistent_coriolis << std::endl;
165  amrex::Print() << "pbl_ysu_force_over_water : " << pbl_ysu_force_over_water << std::endl;
166  amrex::Print() << "pbl_ysu_land_Ribcr : " << pbl_ysu_land_Ribcr << std::endl;
167  amrex::Print() << "pbl_ysu_unst_Ribcr : " << pbl_ysu_unst_Ribcr << std::endl;
168  }
169  }
170 
171  // Default prefix
172  std::string pp_prefix {"erf"};
173 
174  // LES model
175  LESType les_type;
176  // Smagorinsky Cs coefficient
177  amrex::Real Cs = 0.0;
178  // Smagorinsky CI coefficient
179  amrex::Real CI = 0.0;
180  // Smagorinsky Turbulent Prandtl Number
181  amrex::Real Pr_t = amrex::Real(1.0) / amrex::Real(3.0);
182  amrex::Real Pr_t_inv = amrex::Real(3.0);
183  // Smagorinsky Turbulent Schmidt Number
184  amrex::Real Sc_t = 1.0;
185  amrex::Real Sc_t_inv = 1.0;
186 
187  // Deardorff Ce coefficient
188  amrex::Real Ce = 0.93;
189  amrex::Real Ce_wall = 0.0; // if > 0, then set Ce to this at k=0
190  // Deardorff Ck coefficient
191  amrex::Real Ck = 0.1;
192  // Deardorff sigma_k coefficient
193  amrex::Real sigma_k = 0.5;
194 
195  amrex::Real theta_ref = 300.0;
196 
197  // RANS type
198  RANSType rans_type;
199 
200  // PBL model
201  PBLType pbl_type;
202 
204  MYNNLevel2 pbl_mynn_level2; // for filtering
205 
206  // Model coefficients - YSU
207  // TODO: Add parmparse for all of these above
208  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
209  bool pbl_ysu_use_consistent_coriolis = false; // ignore input pbl_ysu_coriolis_freq, take value from ERF coriolis forcing instead
210  bool pbl_ysu_force_over_water = false; // Force YSU to act as if it is over water regardless of other inputs (for testing)
211  amrex::Real pbl_ysu_land_Ribcr = 0.25; // Critical Bulk Richardson number of Land for stable conditions
212  amrex::Real pbl_ysu_unst_Ribcr = 0.0; // Critical Bulk Richardson number for unstable conditions
213 
214  // QKE stuff - default is to use it, if MYNN2.5 PBL is used default is turb transport in Z-direction only
215  bool use_KE = true;
216  bool diffuse_KE_3D = true;
217  bool advect_KE = true;
218 };
219 #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:193
MYNNLevel25 pbl_mynn
Definition: ERF_TurbStruct.H:203
PBLType pbl_type
Definition: ERF_TurbStruct.H:201
bool use_KE
Definition: ERF_TurbStruct.H:215
amrex::Real Sc_t_inv
Definition: ERF_TurbStruct.H:185
void init_params(int lev, int max_level)
Definition: ERF_TurbStruct.H:33
MYNNLevel2 pbl_mynn_level2
Definition: ERF_TurbStruct.H:204
RANSType rans_type
Definition: ERF_TurbStruct.H:198
bool diffuse_KE_3D
Definition: ERF_TurbStruct.H:216
amrex::Real Ck
Definition: ERF_TurbStruct.H:191
std::string pp_prefix
Definition: ERF_TurbStruct.H:172
bool pbl_ysu_use_consistent_coriolis
Definition: ERF_TurbStruct.H:209
amrex::Real pbl_ysu_land_Ribcr
Definition: ERF_TurbStruct.H:211
amrex::Real Cs
Definition: ERF_TurbStruct.H:177
amrex::Real Pr_t_inv
Definition: ERF_TurbStruct.H:182
amrex::Real Ce_wall
Definition: ERF_TurbStruct.H:189
amrex::Real Ce
Definition: ERF_TurbStruct.H:188
LESType les_type
Definition: ERF_TurbStruct.H:175
amrex::Real pbl_ysu_coriolis_freq
Definition: ERF_TurbStruct.H:208
bool pbl_ysu_force_over_water
Definition: ERF_TurbStruct.H:210
void display(int lev)
Definition: ERF_TurbStruct.H:120
amrex::Real Pr_t
Definition: ERF_TurbStruct.H:181
amrex::Real theta_ref
Definition: ERF_TurbStruct.H:195
amrex::Real Sc_t
Definition: ERF_TurbStruct.H:184
bool advect_KE
Definition: ERF_TurbStruct.H:217
amrex::Real pbl_ysu_unst_Ribcr
Definition: ERF_TurbStruct.H:212
amrex::Real CI
Definition: ERF_TurbStruct.H:179