ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
TurbChoice Struct Reference

#include <TurbStruct.H>

Collaboration diagram for TurbChoice:

Public Member Functions

void init_params (int lev, int max_level)
 
void display (int lev)
 

Public Attributes

std::string pp_prefix {"erf"}
 
LESType les_type
 
amrex::Real Cs = 0.0
 
amrex::Real CI = 0.0
 
amrex::Real Pr_t = 1. / 3.
 
amrex::Real Pr_t_inv = 3.0
 
amrex::Real Sc_t = 1.0
 
amrex::Real Sc_t_inv = 1.0
 
amrex::Real Ce = 0.93
 
amrex::Real Ce_wall = 0.0
 
amrex::Real Ck = 0.1
 
amrex::Real sigma_k = 0.5
 
amrex::Real theta_ref = 300.0
 
PBLType pbl_type
 
amrex::Real pbl_mynn_A1 = 1.18
 
amrex::Real pbl_mynn_A2 = 0.665
 
amrex::Real pbl_mynn_B1 = 24.0
 
amrex::Real pbl_mynn_B2 = 15.0
 
amrex::Real pbl_mynn_C1 = 0.137
 
amrex::Real pbl_mynn_C2 = 0.75
 
amrex::Real pbl_mynn_C3 = 0.352
 
amrex::Real pbl_mynn_C4 = 0.0
 
amrex::Real pbl_mynn_C5 = 0.2
 
amrex::Real pbl_ysu_coriolis_freq = 1.0e-4
 
bool pbl_ysu_over_land = false
 
amrex::Real pbl_ysu_land_Ribcr = 0.25
 
amrex::Real pbl_ysu_unst_Ribcr = 0.0
 
bool use_QKE = false
 
bool diffuse_QKE_3D = false
 
bool advect_QKE = true
 

Detailed Description

Container holding quantities related to turbulence parametrizations

Member Function Documentation

◆ display()

void TurbChoice::display ( int  lev)
inline
208  {
209  amrex::Print() << "Turbulence Settings at level " << lev << std::endl;
210 
212  amrex::Print() << "Using DNS model at level " << lev << std::endl;
213  } else if (les_type == LESType::Smagorinsky) {
214  amrex::Print() << "Using Smagorinsky LES model at level " << lev << std::endl;
215  } else if (les_type == LESType::Deardorff) {
216  amrex::Print() << "Using Deardorff LES model at level " << lev << std::endl;
217  } else if (pbl_type == PBLType::MYNN25) {
218  amrex::Print() << "Using MYNN2.5 PBL model at level " << lev << std::endl;
219  } else if (pbl_type == PBLType::YSU) {
220  amrex::Print() << "Using YSU PBL model at level " << lev << std::endl;
221  }
222 
223  if (les_type != LESType::None) {
225  amrex::Print() << "Cs : " << Cs << std::endl;
226  amrex::Print() << "CI : " << CI << std::endl;
227  amrex::Print() << "Pr_t : " << Pr_t << std::endl;
228  amrex::Print() << "Sc_t : " << Sc_t << std::endl;
229  }
230  if (les_type == LESType::Deardorff) {
231  amrex::Print() << "Ce : " << Ce << std::endl;
232  amrex::Print() << "Ce at wall : " << Ce_wall << std::endl;
233  amrex::Print() << "Ck : " << Ck << std::endl;
234  amrex::Print() << "sigma_k : " << sigma_k << std::endl;
235  }
236  amrex::Print() << "reference theta : " << theta_ref << std::endl;
237  }
238 
239  if (pbl_type != PBLType::None) {
240  amrex::Print() << "pbl_mynn_A1 : " << pbl_mynn_A1 << std::endl;
241  amrex::Print() << "pbl_mynn_A2 : " << pbl_mynn_A2 << std::endl;
242  amrex::Print() << "pbl_mynn_B1 : " << pbl_mynn_B1 << std::endl;
243  amrex::Print() << "pbl_mynn_B2 : " << pbl_mynn_B2 << std::endl;
244  amrex::Print() << "pbl_mynn_C1 : " << pbl_mynn_C1 << std::endl;
245  amrex::Print() << "pbl_mynn_C2 : " << pbl_mynn_C2 << std::endl;
246  amrex::Print() << "pbl_mynn_C3 : " << pbl_mynn_C3 << std::endl;
247  amrex::Print() << "pbl_mynn_C4 : " << pbl_mynn_C4 << std::endl;
248  amrex::Print() << "pbl_mynn_C5 : " << pbl_mynn_C5 << std::endl;
249  }
250  }
amrex::Real pbl_mynn_A2
Definition: TurbStruct.H:282
amrex::Real pbl_mynn_C2
Definition: TurbStruct.H:286
amrex::Real sigma_k
Definition: TurbStruct.H:274
PBLType pbl_type
Definition: TurbStruct.H:279
amrex::Real pbl_mynn_A1
Definition: TurbStruct.H:281
amrex::Real pbl_mynn_C3
Definition: TurbStruct.H:287
amrex::Real Ck
Definition: TurbStruct.H:272
amrex::Real pbl_mynn_B2
Definition: TurbStruct.H:284
amrex::Real pbl_mynn_C1
Definition: TurbStruct.H:285
amrex::Real pbl_mynn_C5
Definition: TurbStruct.H:289
amrex::Real Cs
Definition: TurbStruct.H:258
amrex::Real Ce_wall
Definition: TurbStruct.H:270
amrex::Real Ce
Definition: TurbStruct.H:269
LESType les_type
Definition: TurbStruct.H:256
amrex::Real pbl_mynn_B1
Definition: TurbStruct.H:283
amrex::Real Pr_t
Definition: TurbStruct.H:262
amrex::Real theta_ref
Definition: TurbStruct.H:276
amrex::Real Sc_t
Definition: TurbStruct.H:265
amrex::Real pbl_mynn_C4
Definition: TurbStruct.H:288
amrex::Real CI
Definition: TurbStruct.H:260

◆ init_params()

void TurbChoice::init_params ( int  lev,
int  max_level 
)
inline
19  {
20  amrex::ParmParse pp(pp_prefix);
21 
22  int nvals_les = pp.countval("les_type");
23  int nvals_pbl = pp.countval("pbl_type");
24 
25  // If nvals_les and nvals_pbl are both zero then we skip everything below
26  if ( (nvals_les == 0) && (nvals_pbl == 0) ) {
29  return;
30  }
31 
32  // If nvals_les and nvals_pbl are both > 0 then they must be the same, and either 1 or max_levels
33  if ( (nvals_les > 0) && (nvals_pbl > 0) ) {
34  if (nvals_les != nvals_pbl) {
35  amrex::Error("If specifying both, we must specify the same number of values for les_type and pbl_type");
36  }
37  if (nvals_les != 1 && nvals_les != max_level+1) {
38  amrex::Error("If specifying both, we must either input one value for all levels, or one value per level");
39  }
40  }
41 
42  // Here we cover the case where either one is 1 and the other is 0, or they both are = 1
43  if (nvals_les == 1 || nvals_pbl == 1) {
44 
45  // Which LES closure?
46  std::string les_type_string = "None";
47  pp.query("les_type",les_type_string);
48 
49  if (!les_type_string.compare("Smagorinsky")) {
51  } else if (!les_type_string.compare("Deardorff")) {
53  } else if (!les_type_string.compare("None")) {
54  les_type = LESType::None; // Means DNS
55  } else {
56  amrex::Error("Don't know this les_type");
57  }
58 
59  // Which PBL Closure
60  static std::string pbl_type_string = "None";
61  pp.query("pbl_type",pbl_type_string);
62  if (pbl_type_string == "MYNN2.5") {
64  } else if (pbl_type_string == "YSU") {
66  } else if (pbl_type_string == "None") {
68  } else {
69  amrex::Error("Don't know this pbl_type");
70  }
71 
72  // Do some more stuff for PBL Modeling
73  if (pbl_type != PBLType::None) {
74  // Check for compatibility between PBL, LES, Molec Transport
76  amrex::Error("It is not recommended to use Deardorff LES and a PBL model");
77  } else if (les_type != LESType::None) {
78  amrex::Print() << "Selected a PBL model and an LES model: " <<
79  "Using PBL for vertical transport, LES for horizontal" << std::endl;
80  }
81  pp.query("pbl_mynn_A1", pbl_mynn_A1);
82  pp.query("pbl_mynn_A2", pbl_mynn_A2);
83  pp.query("pbl_mynn_B1", pbl_mynn_B1);
84  pp.query("pbl_mynn_B2", pbl_mynn_B2);
85  pp.query("pbl_mynn_C1", pbl_mynn_C1);
86  pp.query("pbl_mynn_C2", pbl_mynn_C2);
87  pp.query("pbl_mynn_C3", pbl_mynn_C3);
88  pp.query("pbl_mynn_C4", pbl_mynn_C4);
89  pp.query("pbl_mynn_C5", pbl_mynn_C5);
90  }
91 
92  // Right now, solving the QKE equation is only supported when MYNN PBL is turned on
93  if (pbl_type == PBLType::MYNN25 ||
94  pbl_type == PBLType::YSU ) {
95  use_QKE = true;
96  pp.query("diffuse_QKE_3D", diffuse_QKE_3D);
97  pp.query("advect_QKE", advect_QKE);
98  }
99 
100  // LES constants...
101  pp.query("Cs" , Cs);
102  pp.query("CI" , CI);
103  pp.query("Pr_t", Pr_t);
104  pp.query("Sc_t", Sc_t);
105 
106  // Compute relevant forms of diffusion parameters
107  Pr_t_inv = 1.0 / Pr_t;
108  Sc_t_inv = 1.0 / Sc_t;
109 
110  pp.query("Ce" , Ce);
111  pp.query("Ce_wall" , Ce_wall);
112  pp.query("sigma_k" , sigma_k);
113 
114  if (les_type == LESType::Deardorff) {
115  pp.query("Ck" , Ck);
116  }
117 
118  pp.query("theta_ref", theta_ref);
119 
120  // Validate inputs
122  if (Cs == 0) {
123  amrex::Error("Need to specify Cs for Smagorsinky LES");
124  }
125  }
126 
127  // Here we cover the case where either one is > 1 and the other is 0, or they both are the same and > 1
128  } else {
129 
130  // Which LES closure?
131  std::string les_type_string = "None";
132  pp.get("les_type", les_type_string, lev);
133 
134  if (!les_type_string.compare("Smagorinsky")) {
136  } else if (!les_type_string.compare("Deardorff")) {
138  } else if (!les_type_string.compare("None")) {
139  les_type = LESType::None; // Means DNS
140  } else {
141  amrex::Error("Don't know this les_type");
142  }
143 
144  // Which PBL Closure
145  static std::string pbl_type_string = "None";
146  pp.get("pbl_type", pbl_type_string, lev);
147  if (pbl_type_string == "MYNN2.5") {
149  } else if (pbl_type_string == "YSU") {
151  } else if (pbl_type_string == "None") {
153  } else {
154  amrex::Error("Don't know this pbl_type");
155  }
156 
157  // Do some more stuff for PBL Modeling
158  if (pbl_type != PBLType::None) {
159  // Check for compatibility between PBL, LES, Molec Transport
160  if (les_type != LESType::None) {
161  amrex::Print() << "Selected a PBL model and an LES model: " <<
162  "Using PBL for vertical transport, LES for horizontal" << std::endl;
163  } else if (les_type == LESType::Deardorff) {
164  amrex::Error("It is not recommended to use Deardorff LES and a PBL model");
165  }
166  pp.query("pbl_mynn_A1", pbl_mynn_A1, lev);
167  pp.query("pbl_mynn_A2", pbl_mynn_A2, lev);
168  pp.query("pbl_mynn_B1", pbl_mynn_B1, lev);
169  pp.query("pbl_mynn_B2", pbl_mynn_B2, lev);
170  pp.query("pbl_mynn_C1", pbl_mynn_C1, lev);
171  pp.query("pbl_mynn_C2", pbl_mynn_C2, lev);
172  pp.query("pbl_mynn_C3", pbl_mynn_C3, lev);
173  pp.query("pbl_mynn_C4", pbl_mynn_C4, lev);
174  pp.query("pbl_mynn_C5", pbl_mynn_C5, lev);
175  }
176 
177  // Right now, solving the QKE equation is only supported when MYNN PBL is turned on
178  if (pbl_type == PBLType::MYNN25 ||
179  pbl_type == PBLType::YSU ) {
180  use_QKE = true;
181  pp.query("diffuse_QKE_3D", diffuse_QKE_3D, lev);
182  pp.query("advect_QKE" , advect_QKE, lev);
183  }
184 
185  // LES constants...
186  pp.query("Cs" ,Cs, lev);
187  pp.query("CI" ,CI, lev);
188  pp.query("Pr_t",Pr_t, lev);
189  pp.query("Sc_t",Sc_t, lev);
190 
191  // Compute relevant forms of diffusion parameters
192  Pr_t_inv = 1.0 / Pr_t;
193  Sc_t_inv = 1.0 / Sc_t;
194 
195  pp.query("Ce" , Ce, lev);
196  pp.query("Ce_wall" , Ce_wall, lev);
197  pp.query("sigma_k" , sigma_k, lev);
198 
199  if (les_type == LESType::Deardorff) {
200  pp.query("Ck" , Ck, lev);
201  }
202 
203  pp.query("theta_ref", theta_ref, lev);
204  }
205  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: Microphysics_Utils.H:183
bool advect_QKE
Definition: TurbStruct.H:301
amrex::Real Sc_t_inv
Definition: TurbStruct.H:266
std::string pp_prefix
Definition: TurbStruct.H:253
bool use_QKE
Definition: TurbStruct.H:299
amrex::Real Pr_t_inv
Definition: TurbStruct.H:263
bool diffuse_QKE_3D
Definition: TurbStruct.H:300
Here is the call graph for this function:

Member Data Documentation

◆ advect_QKE

bool TurbChoice::advect_QKE = true

◆ Ce

amrex::Real TurbChoice::Ce = 0.93

◆ Ce_wall

amrex::Real TurbChoice::Ce_wall = 0.0

◆ CI

amrex::Real TurbChoice::CI = 0.0

Referenced by display(), and init_params().

◆ Ck

amrex::Real TurbChoice::Ck = 0.1

◆ Cs

amrex::Real TurbChoice::Cs = 0.0

◆ diffuse_QKE_3D

bool TurbChoice::diffuse_QKE_3D = false

◆ les_type

◆ pbl_mynn_A1

amrex::Real TurbChoice::pbl_mynn_A1 = 1.18

◆ pbl_mynn_A2

amrex::Real TurbChoice::pbl_mynn_A2 = 0.665

◆ pbl_mynn_B1

amrex::Real TurbChoice::pbl_mynn_B1 = 24.0

◆ pbl_mynn_B2

amrex::Real TurbChoice::pbl_mynn_B2 = 15.0

◆ pbl_mynn_C1

amrex::Real TurbChoice::pbl_mynn_C1 = 0.137

◆ pbl_mynn_C2

amrex::Real TurbChoice::pbl_mynn_C2 = 0.75

◆ pbl_mynn_C3

amrex::Real TurbChoice::pbl_mynn_C3 = 0.352

◆ pbl_mynn_C4

amrex::Real TurbChoice::pbl_mynn_C4 = 0.0

Referenced by display(), and init_params().

◆ pbl_mynn_C5

amrex::Real TurbChoice::pbl_mynn_C5 = 0.2

◆ pbl_type

◆ pbl_ysu_coriolis_freq

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

◆ pbl_ysu_land_Ribcr

amrex::Real TurbChoice::pbl_ysu_land_Ribcr = 0.25

◆ pbl_ysu_over_land

bool TurbChoice::pbl_ysu_over_land = false

◆ pbl_ysu_unst_Ribcr

amrex::Real TurbChoice::pbl_ysu_unst_Ribcr = 0.0

◆ pp_prefix

std::string TurbChoice::pp_prefix {"erf"}

Referenced by init_params().

◆ Pr_t

amrex::Real TurbChoice::Pr_t = 1. / 3.

Referenced by display(), and init_params().

◆ Pr_t_inv

amrex::Real TurbChoice::Pr_t_inv = 3.0

◆ Sc_t

amrex::Real TurbChoice::Sc_t = 1.0

Referenced by display(), and init_params().

◆ Sc_t_inv

amrex::Real TurbChoice::Sc_t_inv = 1.0

◆ sigma_k

amrex::Real TurbChoice::sigma_k = 0.5

◆ theta_ref

amrex::Real TurbChoice::theta_ref = 300.0

◆ use_QKE


The documentation for this struct was generated from the following file: