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