ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_SAM.H
Go to the documentation of this file.
1 /*
2  * Implementation 1-moment microphysics model
3  * NOTE: this model is based on the SAM code, and the Klemp's paper
4  * 1): Joseph, Klemp, the simulation of three-dimensional convective storm dynamics,
5  * Journal of the atmospheric sciences, vol35, p1070
6  * 2): Marat Khairoutdinov and David Randall, cloud resolving modeling of the ARM summer 1997 IOP:
7  * model formulation, results, unvertainties, and sensitivities, Journal of the atmospheric sciences, vol60, p607
8  */
9 #ifndef ERF_SAM_H
10 #define ERF_SAM_H
11 
12 #include <string>
13 #include <vector>
14 #include <memory>
15 
16 #include <AMReX_FArrayBox.H>
17 #include <AMReX_Geometry.H>
18 #include <AMReX_TableData.H>
19 #include <AMReX_MultiFabUtil.H>
20 
21 #include "ERF_Constants.H"
22 #include "ERF_MicrophysicsUtils.H"
23 #include "ERF_IndexDefines.H"
24 #include "ERF_DataStruct.H"
25 #include "ERF_NullMoist.H"
26 
27 namespace MicVar {
28  enum {
29  // independent variables
30  rho=0, // density
31  theta, // liquid/ice water potential temperature
32  tabs, // temperature
33  pres, // pressure
34  // non-precipitating vars
35  qt, // total cloud
36  qn, // cloud condesnate (liquid + ice)
37  qv, // cloud vapor
38  qcl, // cloud water
39  qci, // cloud ice
40  // precipitating vars
41  qp, // total precip
42  qpr, // precip rain
43  qps, // precip ice
44  qpg, // graupel
45  // derived vars
49  NumVars
50  };
51 }
52 
53 class SAM : public NullMoist {
54 
55  using FabPtr = std::shared_ptr<amrex::MultiFab>;
56 
57 public:
58  // constructor
59  SAM () {}
60 
61  // destructor
62  virtual ~SAM () = default;
63 
64  // cloud physics
65  void Cloud (const SolverChoice& sc);
66 
67  // ice physics
68  void IceFall (const SolverChoice& sc);
69 
70  // precip
71  void Precip (const SolverChoice& sc);
72 
73  // precip fall
74  void PrecipFall (const SolverChoice& sc);
75 
76  // Set up for first time
77  void
78  Define (SolverChoice& sc) override
79  {
80  m_fac_cond = lcond / sc.c_p;
81  m_fac_fus = lfus / sc.c_p;
82  m_fac_sub = lsub / sc.c_p;
83  m_axis = sc.ave_plane;
84  m_rdOcp = sc.rdOcp;
85  m_do_cond = (!sc.use_shoc);
86  }
87 
88  // init
89  void
90  Init (const amrex::MultiFab& cons_in,
91  const amrex::BoxArray& grids,
92  const amrex::Geometry& geom,
93  const amrex::Real& dt_advance,
94  std::unique_ptr<amrex::MultiFab>& z_phys_nd,
95  std::unique_ptr<amrex::MultiFab>& detJ_cc) override;
96 
97  // Import minimum dz at this level
98  void
99  Set_dzmin (const amrex::Real dz_min) override
100  {
101  m_dzmin = dz_min;
102  }
103 
104  // Copy state into micro vars
105  void
106  Copy_State_to_Micro (const amrex::MultiFab& cons_in) override;
107 
108  // Copy state into micro vars
109  void
110  Copy_Micro_to_State (amrex::MultiFab& cons_in) override;
111 
112  void
113  Update_Micro_Vars (amrex::MultiFab& cons_in) override
114  {
115  this->Copy_State_to_Micro(cons_in);
116  this->Compute_Coefficients();
117  }
118 
119  void
120  Update_State_Vars (amrex::MultiFab& cons_in) override
121  {
122  this->Copy_Micro_to_State(cons_in);
123  }
124 
125  // wrapper to do all the updating
126  void
127  Advance (const amrex::Real& dt_advance,
128  const SolverChoice& sc) override
129  {
130  dt = dt_advance;
131 
132  this->Cloud(sc);
133  this->IceFall(sc);
134  this->Precip(sc);
135  this->PrecipFall(sc);
136  }
137 
138  amrex::MultiFab*
139  Qmoist_Ptr (const int& varIdx) override
140  {
141  AMREX_ALWAYS_ASSERT(varIdx < m_qmoist_size);
142  return mic_fab_vars[MicVarMap[varIdx]].get();
143  }
144 
145  void
147 
148  int
149  Qmoist_Size () override { return SAM::m_qmoist_size; }
150 
151  int
153 
154  AMREX_GPU_HOST_DEVICE
155  AMREX_FORCE_INLINE
156  static amrex::Real
157  NewtonIterSat (int& i, int& j, int& k,
158  const int& SAM_moisture_type,
159  const amrex::Real& fac_cond,
160  const amrex::Real& /*fac_fus*/,
161  const amrex::Real& fac_sub,
162  const amrex::Real& an,
163  const amrex::Real& bn,
164  const amrex::Array4<amrex::Real>& tabs_array,
165  const amrex::Array4<amrex::Real>& pres_array,
166  const amrex::Array4<amrex::Real>& qv_array,
167  const amrex::Array4<amrex::Real>& qc_array,
168  const amrex::Array4<amrex::Real>& qi_array,
169  const amrex::Array4<amrex::Real>& qn_array,
170  const amrex::Array4<amrex::Real>& qt_array)
171  {
172  // Solution tolerance
173  amrex::Real tol = 1.0e-4;
174 
175  // Saturation moisture fractions
176  amrex::Real omn, domn;
177  amrex::Real qsat, dqsat;
178  amrex::Real qsatw, dqsatw;
179  amrex::Real qsati, dqsati;
180 
181  // Newton iteration vars
182  int niter;
183  amrex::Real fff, dfff, dtabs;
184  amrex::Real lstar, dlstar;
185  amrex::Real lstarw, lstari;
186  amrex::Real delta_qv, delta_qc, delta_qi;
187 
188  // Initial guess for temperature & pressure
189  amrex::Real tabs = tabs_array(i,j,k);
190  amrex::Real pres = pres_array(i,j,k);
191 
192  niter = 0;
193  dtabs = 1;
194  //==================================================
195  // Newton iteration to qv=qsat (cloud phase only)
196  //==================================================
197  do {
198  // Latent heats and their derivatives wrt to T
199  lstarw = fac_cond;
200  lstari = fac_sub;
201  domn = 0.0;
202 
203  // Saturation moisture fractions
204  erf_qsatw(tabs, pres, qsatw);
205  erf_qsati(tabs, pres, qsati);
206  erf_dtqsatw(tabs, pres, dqsatw);
207  erf_dtqsati(tabs, pres, dqsati);
208 
209  if (SAM_moisture_type == 1) {
210  // Cloud ice not permitted (condensation & fusion)
211  if(tabs >= tbgmax) {
212  omn = 1.0;
213  }
214  // Cloud water not permitted (sublimation & fusion)
215  else if(tabs <= tbgmin) {
216  omn = 0.0;
217  }
218  // Mixed cloud phase (condensation & fusion)
219  else {
220  omn = an*tabs-bn;
221  domn = an;
222  }
223  } else if (SAM_moisture_type == 2) {
224  omn = 1.0;
225  domn = 0.0;
226  }
227 
228  // Linear combination of each component
229  qsat = omn * qsatw + (1.0-omn) * qsati;
230  dqsat = omn * dqsatw + (1.0-omn) * dqsati
231  + domn * qsatw - domn * qsati;
232  lstar = omn * lstarw + (1.0-omn) * lstari;
233  dlstar = domn * lstarw - domn * lstari;
234 
235  // Function for root finding:
236  // 0 = -T_new + T_old + L_eff/C_p * (qv - qsat)
237  fff = -tabs + tabs_array(i,j,k) + lstar*(qv_array(i,j,k) - qsat);
238 
239  // Derivative of function (T_new iterated on)
240  dfff = -1.0 + dlstar*(qv_array(i,j,k) - qsat) - lstar*dqsat;
241 
242  // Update the temperature
243  dtabs = -fff/dfff;
244  tabs += dtabs;
245 
246  // Update iteration
247  niter = niter+1;
248  } while(std::abs(dtabs) > tol && niter < 20);
249 
250  // Update qsat from last iteration (dq = dq/dt * dt)
251  qsat += dqsat*dtabs;
252 
253  // Changes in each component
254  delta_qv = qv_array(i,j,k) - qsat;
255  delta_qc = std::max(-qc_array(i,j,k), delta_qv * omn);
256  delta_qi = std::max(-qi_array(i,j,k), delta_qv * (1.0-omn));
257 
258  // Partition the change in non-precipitating q
259  qv_array(i,j,k) = qsat;
260  qc_array(i,j,k) += delta_qc;
261  qi_array(i,j,k) += delta_qi;
262  qn_array(i,j,k) = qc_array(i,j,k) + qi_array(i,j,k);
263  qt_array(i,j,k) = qv_array(i,j,k) + qn_array(i,j,k);
264 
265  // Return to temperature
266  return tabs;
267  }
268 
269  void
271  std::vector<int>& a_idx,
272  std::vector<std::string>& a_names) const override
273  {
274  a_idx.clear();
275  a_names.clear();
276 
277  // NOTE: These are the indices to access into qmoist (not mic_fab_vars)
278  a_idx.push_back(0); a_names.push_back("RainAccum");
279  a_idx.push_back(1); a_names.push_back("SnowAccum");
280  a_idx.push_back(2); a_names.push_back("GraupAccum");
281  }
282 
283 private:
284  // Number of qmoist variables (rain_accum, snow_accum, graup_accum)
285  int m_qmoist_size = 3;
286 
287  // Number of qstate variables
289 
290  // CFL MAX for vertical advection
291  static constexpr amrex::Real CFL_MAX = 0.5;
292 
293  // MicVar map (Qmoist indices -> MicVar enum)
294  amrex::Vector<int> MicVarMap;
295 
296  // geometry
297  amrex::Geometry m_geom;
298 
299  // valid boxes on which to evolve the solution
300  amrex::BoxArray m_gtoe;
301 
302  // timestep
304 
305  // number of vertical levels
306  int nlev, zlo, zhi;
307 
308  // plane average axis
309  int m_axis;
310 
311  // constants
316  bool m_do_cond;
317 
318  // Minimum dz at this level
320 
321  // Pointer to terrain data
322  amrex::MultiFab* m_z_phys_nd;
323  amrex::MultiFab* m_detJ_cc;
324 
325  // independent variables
326  amrex::Array<FabPtr, MicVar::NumVars> mic_fab_vars;
327 
328  // microphysics parameters/coefficients
329  amrex::TableData<amrex::Real, 1> accrrc;
330  amrex::TableData<amrex::Real, 1> accrsi;
331  amrex::TableData<amrex::Real, 1> accrsc;
332  amrex::TableData<amrex::Real, 1> coefice;
333  amrex::TableData<amrex::Real, 1> evaps1;
334  amrex::TableData<amrex::Real, 1> evaps2;
335  amrex::TableData<amrex::Real, 1> accrgi;
336  amrex::TableData<amrex::Real, 1> accrgc;
337  amrex::TableData<amrex::Real, 1> evapg1;
338  amrex::TableData<amrex::Real, 1> evapg2;
339  amrex::TableData<amrex::Real, 1> evapr1;
340  amrex::TableData<amrex::Real, 1> evapr2;
341 
342  // vertical plane average data
343  amrex::TableData<amrex::Real, 1> rho1d;
344  amrex::TableData<amrex::Real, 1> pres1d;
345  amrex::TableData<amrex::Real, 1> tabs1d;
346  amrex::TableData<amrex::Real, 1> qt1d;
347  amrex::TableData<amrex::Real, 1> qv1d;
348  amrex::TableData<amrex::Real, 1> qn1d;
349 };
350 #endif
constexpr amrex::Real lfus
Definition: ERF_Constants.H:67
constexpr amrex::Real tbgmax
Definition: ERF_Constants.H:32
constexpr amrex::Real lsub
Definition: ERF_Constants.H:68
constexpr amrex::Real lcond
Definition: ERF_Constants.H:66
constexpr amrex::Real tbgmin
Definition: ERF_Constants.H:31
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsatw(amrex::Real t, amrex::Real p, amrex::Real &dtqsatw)
Definition: ERF_MicrophysicsUtils.H:181
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsati(amrex::Real t, amrex::Real p, amrex::Real &dtqsati)
Definition: ERF_MicrophysicsUtils.H:173
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:166
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsati(amrex::Real t, amrex::Real p, amrex::Real &qsati)
Definition: ERF_MicrophysicsUtils.H:155
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF_NullMoist.H:8
Definition: ERF_SAM.H:53
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateSAM.cpp:15
void Advance(const amrex::Real &dt_advance, const SolverChoice &sc) override
Definition: ERF_SAM.H:127
int zlo
Definition: ERF_SAM.H:306
void Precip(const SolverChoice &sc)
Definition: ERF_Precip.cpp:10
amrex::TableData< amrex::Real, 1 > rho1d
Definition: ERF_SAM.H:343
amrex::TableData< amrex::Real, 1 > evapr2
Definition: ERF_SAM.H:340
SAM()
Definition: ERF_SAM.H:59
amrex::MultiFab * m_detJ_cc
Definition: ERF_SAM.H:323
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitSAM.cpp:85
amrex::TableData< amrex::Real, 1 > evapg1
Definition: ERF_SAM.H:337
void Update_Micro_Vars(amrex::MultiFab &cons_in) override
Definition: ERF_SAM.H:113
amrex::TableData< amrex::Real, 1 > accrgi
Definition: ERF_SAM.H:335
std::shared_ptr< amrex::MultiFab > FabPtr
Definition: ERF_SAM.H:55
virtual ~SAM()=default
amrex::TableData< amrex::Real, 1 > evapg2
Definition: ERF_SAM.H:338
void IceFall(const SolverChoice &sc)
Definition: ERF_IceFall.cpp:11
amrex::Array< FabPtr, MicVar::NumVars > mic_fab_vars
Definition: ERF_SAM.H:326
int m_axis
Definition: ERF_SAM.H:309
amrex::Real m_rdOcp
Definition: ERF_SAM.H:315
void Cloud(const SolverChoice &sc)
Definition: ERF_CloudSAM.cpp:12
int m_qmoist_size
Definition: ERF_SAM.H:285
amrex::MultiFab * m_z_phys_nd
Definition: ERF_SAM.H:322
amrex::Real m_fac_fus
Definition: ERF_SAM.H:313
amrex::TableData< amrex::Real, 1 > evaps1
Definition: ERF_SAM.H:333
amrex::BoxArray m_gtoe
Definition: ERF_SAM.H:300
amrex::Vector< int > MicVarMap
Definition: ERF_SAM.H:294
bool m_do_cond
Definition: ERF_SAM.H:316
amrex::TableData< amrex::Real, 1 > qv1d
Definition: ERF_SAM.H:347
amrex::Real m_fac_cond
Definition: ERF_SAM.H:312
amrex::TableData< amrex::Real, 1 > accrsi
Definition: ERF_SAM.H:330
amrex::TableData< amrex::Real, 1 > accrrc
Definition: ERF_SAM.H:329
void Compute_Coefficients()
Definition: ERF_InitSAM.cpp:137
amrex::TableData< amrex::Real, 1 > pres1d
Definition: ERF_SAM.H:344
int nlev
Definition: ERF_SAM.H:306
int Qstate_Moist_Size() override
Definition: ERF_SAM.H:152
amrex::Geometry m_geom
Definition: ERF_SAM.H:297
amrex::TableData< amrex::Real, 1 > evapr1
Definition: ERF_SAM.H:339
amrex::Real m_fac_sub
Definition: ERF_SAM.H:314
amrex::TableData< amrex::Real, 1 > evaps2
Definition: ERF_SAM.H:334
void Init(const amrex::MultiFab &cons_in, const amrex::BoxArray &grids, const amrex::Geometry &geom, const amrex::Real &dt_advance, std::unique_ptr< amrex::MultiFab > &z_phys_nd, std::unique_ptr< amrex::MultiFab > &detJ_cc) override
Definition: ERF_InitSAM.cpp:22
void Set_dzmin(const amrex::Real dz_min) override
Definition: ERF_SAM.H:99
amrex::TableData< amrex::Real, 1 > qn1d
Definition: ERF_SAM.H:348
amrex::TableData< amrex::Real, 1 > accrgc
Definition: ERF_SAM.H:336
amrex::TableData< amrex::Real, 1 > accrsc
Definition: ERF_SAM.H:331
int n_qstate_moist_size
Definition: ERF_SAM.H:288
amrex::TableData< amrex::Real, 1 > qt1d
Definition: ERF_SAM.H:346
amrex::TableData< amrex::Real, 1 > coefice
Definition: ERF_SAM.H:332
int zhi
Definition: ERF_SAM.H:306
static constexpr amrex::Real CFL_MAX
Definition: ERF_SAM.H:291
void Define(SolverChoice &sc) override
Definition: ERF_SAM.H:78
void Qmoist_Restart_Vars(const SolverChoice &, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
Definition: ERF_SAM.H:270
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real NewtonIterSat(int &i, int &j, int &k, const int &SAM_moisture_type, const amrex::Real &fac_cond, const amrex::Real &, const amrex::Real &fac_sub, const amrex::Real &an, const amrex::Real &bn, const amrex::Array4< amrex::Real > &tabs_array, const amrex::Array4< amrex::Real > &pres_array, const amrex::Array4< amrex::Real > &qv_array, const amrex::Array4< amrex::Real > &qc_array, const amrex::Array4< amrex::Real > &qi_array, const amrex::Array4< amrex::Real > &qn_array, const amrex::Array4< amrex::Real > &qt_array)
Definition: ERF_SAM.H:157
int Qmoist_Size() override
Definition: ERF_SAM.H:149
void PrecipFall(const SolverChoice &sc)
Definition: ERF_PrecipFall.cpp:16
void Update_State_Vars(amrex::MultiFab &cons_in) override
Definition: ERF_SAM.H:120
amrex::TableData< amrex::Real, 1 > tabs1d
Definition: ERF_SAM.H:345
amrex::MultiFab * Qmoist_Ptr(const int &varIdx) override
Definition: ERF_SAM.H:139
amrex::Real dt
Definition: ERF_SAM.H:303
amrex::Real m_dzmin
Definition: ERF_SAM.H:319
Definition: ERF_SAM.H:27
@ NumVars
Definition: ERF_SAM.H:49
@ pres
Definition: ERF_SAM.H:33
@ qpr
Definition: ERF_SAM.H:42
@ qci
Definition: ERF_SAM.H:39
@ qv
Definition: ERF_SAM.H:37
@ rho
Definition: ERF_SAM.H:30
@ snow_accum
Definition: ERF_SAM.H:47
@ qpg
Definition: ERF_SAM.H:44
@ qt
Definition: ERF_SAM.H:35
@ qn
Definition: ERF_SAM.H:36
@ rain_accum
Definition: ERF_SAM.H:46
@ theta
Definition: ERF_SAM.H:31
@ qcl
Definition: ERF_SAM.H:38
@ graup_accum
Definition: ERF_SAM.H:48
@ tabs
Definition: ERF_SAM.H:32
@ qps
Definition: ERF_SAM.H:43
@ qp
Definition: ERF_SAM.H:41
Definition: ERF_DataStruct.H:123
amrex::Real rdOcp
Definition: ERF_DataStruct.H:975
amrex::Real c_p
Definition: ERF_DataStruct.H:974
bool use_shoc
Definition: ERF_DataStruct.H:1004
int ave_plane
Definition: ERF_DataStruct.H:1035