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
50  NumVars
51  };
52 }
53 
54 class SAM : public NullMoist {
55 
56  using FabPtr = std::shared_ptr<amrex::MultiFab>;
57 
58 public:
59  // constructor
60  SAM () {}
61 
62  // destructor
63  virtual ~SAM () = default;
64 
65  // cloud physics
66  void Cloud (const SolverChoice& sc);
67 
68  // ice physics
69  void IceFall (const SolverChoice& sc);
70 
71  // precip
72  void Precip (const SolverChoice& sc);
73 
74  // precip fall
75  void PrecipFall (const SolverChoice& sc);
76 
77  // Set up for first time
78  void
79  Define (SolverChoice& sc) override
80  {
81  docloud = sc.do_cloud;
82  doprecip = sc.do_precip;
83  m_fac_cond = lcond / sc.c_p;
84  m_fac_fus = lfus / sc.c_p;
85  m_fac_sub = lsub / sc.c_p;
86  m_gOcp = CONST_GRAV / sc.c_p;
87  m_axis = sc.ave_plane;
88  m_rdOcp = sc.rdOcp;
89  }
90 
91  // init
92  void
93  Init (const amrex::MultiFab& cons_in,
94  const amrex::BoxArray& grids,
95  const amrex::Geometry& geom,
96  const amrex::Real& dt_advance,
97  std::unique_ptr<amrex::MultiFab>& z_phys_nd,
98  std::unique_ptr<amrex::MultiFab>& detJ_cc) override;
99 
100  // Copy state into micro vars
101  void
102  Copy_State_to_Micro (const amrex::MultiFab& cons_in) override;
103 
104  // Copy state into micro vars
105  void
106  Copy_Micro_to_State (amrex::MultiFab& cons_in) override;
107 
108  void
109  Update_Micro_Vars (amrex::MultiFab& cons_in) override
110  {
111  this->Copy_State_to_Micro(cons_in);
112  this->Compute_Coefficients();
113  }
114 
115  void
116  Update_State_Vars (amrex::MultiFab& cons_in) override
117  {
118  this->Copy_Micro_to_State(cons_in);
119  }
120 
121  // wrapper to do all the updating
122  void
123  Advance (const amrex::Real& dt_advance,
124  const SolverChoice& sc) override
125  {
126  dt = dt_advance;
127 
128  this->Cloud(sc);
129  this->IceFall(sc);
130  this->Precip(sc);
131  this->PrecipFall(sc);
132  }
133 
134  amrex::MultiFab*
135  Qmoist_Ptr (const int& varIdx) override
136  {
137  AMREX_ALWAYS_ASSERT(varIdx < m_qmoist_size);
138  return mic_fab_vars[MicVarMap[varIdx]].get();
139  }
140 
141  void
143 
144  int
145  Qmoist_Size () override { return SAM::m_qmoist_size; }
146 
147  int
148  Qstate_Size () override { return SAM::m_qstate_size; }
149 
150  AMREX_GPU_HOST_DEVICE
151  AMREX_FORCE_INLINE
152  static amrex::Real
153  NewtonIterSat (int& i, int& j, int& k,
154  const int& SAM_moisture_type,
155  const amrex::Real& fac_cond,
156  const amrex::Real& fac_fus,
157  const amrex::Real& fac_sub,
158  const amrex::Real& an,
159  const amrex::Real& bn,
160  const amrex::Array4<amrex::Real>& tabs_array,
161  const amrex::Array4<amrex::Real>& pres_array,
162  const amrex::Array4<amrex::Real>& qv_array,
163  const amrex::Array4<amrex::Real>& qc_array,
164  const amrex::Array4<amrex::Real>& qi_array,
165  const amrex::Array4<amrex::Real>& qn_array,
166  const amrex::Array4<amrex::Real>& qt_array)
167  {
168  // Solution tolerance
169  amrex::Real tol = 1.0e-4;
170 
171  // Saturation moisture fractions
172  amrex::Real omn, domn;
173  amrex::Real qsat, dqsat;
174  amrex::Real qsatw, dqsatw;
175  amrex::Real qsati, dqsati;
176 
177  // Newton iteration vars
178  int niter;
179  amrex::Real fff, dfff, dtabs;
180  amrex::Real lstar, dlstar;
181  amrex::Real lstarw, lstari;
182  amrex::Real delta_qv, delta_qc, delta_qi;
183 
184  // Initial guess for temperature & pressure
185  amrex::Real tabs = tabs_array(i,j,k);
186  amrex::Real pres = pres_array(i,j,k);
187 
188  niter = 0;
189  dtabs = 1;
190  //==================================================
191  // Newton iteration to qv=qsat (cloud phase only)
192  //==================================================
193  do {
194  // Latent heats and their derivatives wrt to T
195  lstarw = fac_cond;
196  lstari = fac_fus;
197  domn = 0.0;
198 
199  // Saturation moisture fractions
200  erf_qsatw(tabs, pres, qsatw);
201  erf_qsati(tabs, pres, qsati);
202  erf_dtqsatw(tabs, pres, dqsatw);
203  erf_dtqsati(tabs, pres, dqsati);
204 
205  if (SAM_moisture_type == 1) {
206  // Cloud ice not permitted (condensation & fusion)
207  if(tabs >= tbgmax) {
208  omn = 1.0;
209  }
210  // Cloud water not permitted (sublimation & fusion)
211  else if(tabs <= tbgmin) {
212  omn = 0.0;
213  lstarw = fac_sub;
214  }
215  // Mixed cloud phase (condensation & fusion)
216  else {
217  omn = an*tabs-bn;
218  domn = an;
219  }
220  } else if (SAM_moisture_type == 2) {
221  omn = 1.0;
222  domn = 0.0;
223  }
224 
225  // Linear combination of each component
226  qsat = omn * qsatw + (1.0-omn) * qsati;
227  dqsat = omn * dqsatw + (1.0-omn) * dqsati
228  + domn * qsatw - domn * qsati;
229  lstar = omn * lstarw + (1.0-omn) * lstari;
230  dlstar = domn * lstarw - domn * lstari;
231 
232  // Function for root finding:
233  // 0 = -T_new + T_old + L_eff/C_p * (qv - qsat)
234  fff = -tabs + tabs_array(i,j,k) + lstar*(qv_array(i,j,k) - qsat);
235 
236  // Derivative of function (T_new iterated on)
237  dfff = -1.0 + dlstar*(qv_array(i,j,k) - qsat) - lstar*dqsat;
238 
239  // Update the temperature
240  dtabs = -fff/dfff;
241  tabs += dtabs;
242 
243  // Update iteration
244  niter = niter+1;
245  } while(std::abs(dtabs) > tol && niter < 20);
246 
247  // Update qsat from last iteration (dq = dq/dt * dt)
248  qsat += dqsat*dtabs;
249 
250  // Changes in each component
251  delta_qv = qv_array(i,j,k) - qsat;
252  delta_qc = std::max(-qc_array(i,j,k), delta_qv * omn);
253  delta_qi = std::max(-qi_array(i,j,k), delta_qv * (1.0-omn));
254 
255  // Partition the change in non-precipitating q
256  qv_array(i,j,k) = qsat;
257  qc_array(i,j,k) += delta_qc;
258  qi_array(i,j,k) += delta_qi;
259  qn_array(i,j,k) = qc_array(i,j,k) + qi_array(i,j,k);
260  qt_array(i,j,k) = qv_array(i,j,k) + qn_array(i,j,k);
261 
262  // Return to temperature
263  return tabs;
264  }
265 
266  void
268  std::vector<int>& a_idx,
269  std::vector<std::string>& a_names) const override
270  {
271  a_idx.clear();
272  a_names.clear();
273  if (a_sc.moisture_type == MoistureType::SAM) {
274  a_idx.push_back( 8); a_names.push_back("RainAccum");
275  a_idx.push_back( 9); a_names.push_back("SnowAccum");
276  a_idx.push_back(10); a_names.push_back("GraupAccum");
277  } else if (a_sc.moisture_type == MoistureType::SAM_NoIce) {
278  a_idx.push_back( 8); a_names.push_back("RainAccum");
279  }
280  }
281 
282 private:
283  // Number of qmoist variables (qt, qv, qcl, qci, qp, qpr, qps, qpg)
284  int m_qmoist_size = 11;
285 
286  // Number of qstate variables
287  int m_qstate_size = 6;
288 
289  // CFL MAX for vertical advection
290  static constexpr amrex::Real CFL_MAX = 0.5;
291 
292  // MicVar map (Qmoist indices -> MicVar enum)
293  amrex::Vector<int> MicVarMap;
294 
295  // geometry
296  amrex::Geometry m_geom;
297 
298  // valid boxes on which to evolve the solution
299  amrex::BoxArray m_gtoe;
300 
301  // timestep
302  amrex::Real dt;
303 
304  // number of vertical levels
305  int nlev, zlo, zhi;
306 
307  // plane average axis
308  int m_axis;
309 
310  // model options
312 
313  // constants
314  amrex::Real m_fac_cond;
315  amrex::Real m_fac_fus;
316  amrex::Real m_fac_sub;
317  amrex::Real m_gOcp;
318  amrex::Real m_rdOcp;
319 
320  // Pointer to terrain data
321  amrex::MultiFab* m_z_phys_nd;
322  amrex::MultiFab* m_detJ_cc;
323 
324  // independent variables
325  amrex::Array<FabPtr, MicVar::NumVars> mic_fab_vars;
326 
327  // microphysics parameters/coefficients
328  amrex::TableData<amrex::Real, 1> accrrc;
329  amrex::TableData<amrex::Real, 1> accrsi;
330  amrex::TableData<amrex::Real, 1> accrsc;
331  amrex::TableData<amrex::Real, 1> coefice;
332  amrex::TableData<amrex::Real, 1> evaps1;
333  amrex::TableData<amrex::Real, 1> evaps2;
334  amrex::TableData<amrex::Real, 1> accrgi;
335  amrex::TableData<amrex::Real, 1> accrgc;
336  amrex::TableData<amrex::Real, 1> evapg1;
337  amrex::TableData<amrex::Real, 1> evapg2;
338  amrex::TableData<amrex::Real, 1> evapr1;
339  amrex::TableData<amrex::Real, 1> evapr2;
340 
341  // vertical plane average data
342  amrex::TableData<amrex::Real, 1> rho1d;
343  amrex::TableData<amrex::Real, 1> pres1d;
344  amrex::TableData<amrex::Real, 1> tabs1d;
345  amrex::TableData<amrex::Real, 1> qt1d;
346  amrex::TableData<amrex::Real, 1> qv1d;
347  amrex::TableData<amrex::Real, 1> qn1d;
348 
349  amrex::TableData<amrex::Real, 1> gamaz;
350  amrex::TableData<amrex::Real, 1> zmid; // mid value of vertical coordinate in physical domain
351 };
352 #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
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:21
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsatw(amrex::Real t, amrex::Real p, amrex::Real &dtqsatw)
Definition: ERF_MicrophysicsUtils.H:170
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsati(amrex::Real t, amrex::Real p, amrex::Real &dtqsati)
Definition: ERF_MicrophysicsUtils.H:165
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:158
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsati(amrex::Real t, amrex::Real p, amrex::Real &qsati)
Definition: ERF_MicrophysicsUtils.H:151
Definition: ERF_NullMoist.H:8
Definition: ERF_SAM.H:54
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:123
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 &fac_fus, 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:153
int zlo
Definition: ERF_SAM.H:305
void Precip(const SolverChoice &sc)
Definition: ERF_Precip.cpp:10
amrex::TableData< amrex::Real, 1 > rho1d
Definition: ERF_SAM.H:342
amrex::TableData< amrex::Real, 1 > evapr2
Definition: ERF_SAM.H:339
SAM()
Definition: ERF_SAM.H:60
amrex::MultiFab * m_detJ_cc
Definition: ERF_SAM.H:322
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitSAM.cpp:88
amrex::TableData< amrex::Real, 1 > evapg1
Definition: ERF_SAM.H:336
void Update_Micro_Vars(amrex::MultiFab &cons_in) override
Definition: ERF_SAM.H:109
amrex::TableData< amrex::Real, 1 > accrgi
Definition: ERF_SAM.H:334
std::shared_ptr< amrex::MultiFab > FabPtr
Definition: ERF_SAM.H:56
virtual ~SAM()=default
amrex::TableData< amrex::Real, 1 > evapg2
Definition: ERF_SAM.H:337
void IceFall(const SolverChoice &sc)
Definition: ERF_IceFall.cpp:11
amrex::Array< FabPtr, MicVar::NumVars > mic_fab_vars
Definition: ERF_SAM.H:325
int m_axis
Definition: ERF_SAM.H:308
amrex::Real m_gOcp
Definition: ERF_SAM.H:317
amrex::Real m_rdOcp
Definition: ERF_SAM.H:318
void Cloud(const SolverChoice &sc)
Definition: ERF_CloudSAM.cpp:12
int m_qmoist_size
Definition: ERF_SAM.H:284
amrex::MultiFab * m_z_phys_nd
Definition: ERF_SAM.H:321
amrex::Real m_fac_fus
Definition: ERF_SAM.H:315
bool doprecip
Definition: ERF_SAM.H:311
amrex::TableData< amrex::Real, 1 > evaps1
Definition: ERF_SAM.H:332
amrex::BoxArray m_gtoe
Definition: ERF_SAM.H:299
amrex::TableData< amrex::Real, 1 > zmid
Definition: ERF_SAM.H:350
amrex::Vector< int > MicVarMap
Definition: ERF_SAM.H:293
amrex::TableData< amrex::Real, 1 > qv1d
Definition: ERF_SAM.H:346
amrex::Real m_fac_cond
Definition: ERF_SAM.H:314
amrex::TableData< amrex::Real, 1 > accrsi
Definition: ERF_SAM.H:329
amrex::TableData< amrex::Real, 1 > accrrc
Definition: ERF_SAM.H:328
bool docloud
Definition: ERF_SAM.H:311
void Compute_Coefficients()
Definition: ERF_InitSAM.cpp:140
amrex::TableData< amrex::Real, 1 > pres1d
Definition: ERF_SAM.H:343
int nlev
Definition: ERF_SAM.H:305
int m_qstate_size
Definition: ERF_SAM.H:287
amrex::Geometry m_geom
Definition: ERF_SAM.H:296
amrex::TableData< amrex::Real, 1 > evapr1
Definition: ERF_SAM.H:338
amrex::Real m_fac_sub
Definition: ERF_SAM.H:316
amrex::TableData< amrex::Real, 1 > evaps2
Definition: ERF_SAM.H:333
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
amrex::TableData< amrex::Real, 1 > qn1d
Definition: ERF_SAM.H:347
amrex::TableData< amrex::Real, 1 > accrgc
Definition: ERF_SAM.H:335
amrex::TableData< amrex::Real, 1 > accrsc
Definition: ERF_SAM.H:330
amrex::TableData< amrex::Real, 1 > qt1d
Definition: ERF_SAM.H:345
amrex::TableData< amrex::Real, 1 > coefice
Definition: ERF_SAM.H:331
int zhi
Definition: ERF_SAM.H:305
static constexpr amrex::Real CFL_MAX
Definition: ERF_SAM.H:290
void Define(SolverChoice &sc) override
Definition: ERF_SAM.H:79
int Qstate_Size() override
Definition: ERF_SAM.H:148
void Qmoist_Restart_Vars(const SolverChoice &a_sc, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
Definition: ERF_SAM.H:267
int Qmoist_Size() override
Definition: ERF_SAM.H:145
void PrecipFall(const SolverChoice &sc)
Definition: ERF_PrecipFall.cpp:16
void Update_State_Vars(amrex::MultiFab &cons_in) override
Definition: ERF_SAM.H:116
amrex::TableData< amrex::Real, 1 > tabs1d
Definition: ERF_SAM.H:344
amrex::TableData< amrex::Real, 1 > gamaz
Definition: ERF_SAM.H:349
amrex::MultiFab * Qmoist_Ptr(const int &varIdx) override
Definition: ERF_SAM.H:135
amrex::Real dt
Definition: ERF_SAM.H:302
Definition: ERF_SAM.H:27
@ NumVars
Definition: ERF_SAM.H:50
@ 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
@ omega
Definition: ERF_SAM.H:49
@ 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:82
amrex::Real rdOcp
Definition: ERF_DataStruct.H:619
amrex::Real c_p
Definition: ERF_DataStruct.H:618
bool do_precip
Definition: ERF_DataStruct.H:678
MoistureType moisture_type
Definition: ERF_DataStruct.H:664
bool do_cloud
Definition: ERF_DataStruct.H:677
int ave_plane
Definition: ERF_DataStruct.H:675