ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
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"
23 
24 namespace MicVar_Morr {
25  enum {
26  // independent variables
27  rho=0, // density
28  theta, // liquid/ice water potential temperature
29  tabs, // temperature
30  pres, // pressure
31  // non-precipitating vars
32  qt, // total cloud
33  qn, // cloud condesnate (liquid + ice)
34  qv, // cloud vapor
35  qcl, // cloud water
36  qci, // cloud ice
37  // precipitating vars
38  qp, // total precip
39  qpr, // precip rain
40  qps, // precip ice
41  qpg, // graupel
42  // The "numbers" and the "accums" are written/read at checkpoint/restart
43  // number concentrations
44  nc, // cloud droplet number
45  nr, // rain number
46  ni, // cloud ice number
47  ns, // snow number
48  ng, // graupel number
49  // derived vars
54  NumVars
55  };
56 }
57 
58 class Morrison : public NullMoist {
59 
60  using FabPtr = std::shared_ptr<amrex::MultiFab>;
61 
62 public:
63  // constructor
64  Morrison () {}
65 
66  // destructor
67  virtual ~Morrison () = default;
68 
69  // Set up for first time
70  void
71  Define (SolverChoice& sc) override
72  {
73  m_fac_cond = lcond / sc.c_p;
74  m_fac_fus = lfus / sc.c_p;
75  m_fac_sub = lsub / sc.c_p;
76  m_axis = sc.ave_plane;
77  m_rdOcp = sc.rdOcp;
78  m_do_cond = (!sc.use_shoc);
79  }
80 
81  // init
82  void
83  Init (const amrex::MultiFab& cons_in,
84  const amrex::BoxArray& grids,
85  const amrex::Geometry& geom,
86  const amrex::Real& dt_advance,
87  std::unique_ptr<amrex::MultiFab>& z_phys_nd,
88  std::unique_ptr<amrex::MultiFab>& detJ_cc) override;
89 
90  // Import minimum dz at this level
91  void
92  Set_dzmin (const amrex::Real dz_min) override
93  {
94  m_dzmin = dz_min;
95  }
96 
97  // Copy state into micro vars
98  void
99  Copy_State_to_Micro (const amrex::MultiFab& cons_in) override;
100 
101  // Copy state into micro vars
102  void
103  Copy_Micro_to_State (amrex::MultiFab& cons_in) override;
104 
105  void
106  Update_Micro_Vars (amrex::MultiFab& cons_in) override
107  {
108  this->Copy_State_to_Micro(cons_in);
109  this->Compute_Coefficients();
110  }
111 
112  void
113  Update_State_Vars (amrex::MultiFab& cons_in) override
114  {
115  this->Copy_Micro_to_State(cons_in);
116  }
117 
118  // wrapper to do all the updating
119  void
120  Advance (const amrex::Real& dt_advance,
121  const SolverChoice& sc) override;
122 
123  amrex::MultiFab*
124  Qmoist_Ptr (const int& varIdx) override
125  {
126  AMREX_ALWAYS_ASSERT(varIdx < m_qmoist_size);
127  return mic_fab_vars[MicVarMap[varIdx]].get();
128  }
129 
130  void
132 
133  int
134  Qmoist_Size () override { return Morrison::m_qmoist_size; }
135 
136  int
138 
139  void
140  Set_RealWidth (const int real_width) override { m_real_width = real_width; }
141 
142  AMREX_GPU_HOST_DEVICE
143  AMREX_FORCE_INLINE
144  static amrex::Real
145  NewtonIterSat (int& i, int& j, int& k,
146  const amrex::Real& fac_cond,
147  const amrex::Real& fac_fus,
148  const amrex::Real& fac_sub,
149  const amrex::Real& an,
150  const amrex::Real& bn,
151  const amrex::Array4<amrex::Real>& tabs_array,
152  const amrex::Array4<amrex::Real>& pres_array,
153  const amrex::Array4<amrex::Real>& qv_array,
154  const amrex::Array4<amrex::Real>& qc_array,
155  const amrex::Array4<amrex::Real>& qi_array,
156  const amrex::Array4<amrex::Real>& qn_array,
157  const amrex::Array4<amrex::Real>& qt_array)
158  {
159  // Solution tolerance
160  amrex::Real tol = 1.0e-4;
161 
162  // Saturation moisture fractions
163  amrex::Real omn, domn;
164  amrex::Real qsat, dqsat;
165  amrex::Real qsatw, dqsatw;
166  amrex::Real qsati, dqsati;
167 
168  // Newton iteration vars
169  int niter;
170  amrex::Real fff, dfff, dtabs;
171  amrex::Real lstar, dlstar;
172  amrex::Real lstarw, lstari;
173  amrex::Real delta_qv, delta_qc, delta_qi;
174 
175  // Initial guess for temperature & pressure
176  amrex::Real tabs = tabs_array(i,j,k);
177  amrex::Real pres = pres_array(i,j,k);
178 
179  niter = 0;
180  dtabs = 1;
181  //==================================================
182  // Newton iteration to qv=qsat (cloud phase only)
183  //==================================================
184  do {
185  // Latent heats and their derivatives wrt to T
186  lstarw = fac_cond;
187  lstari = fac_fus;
188  domn = 0.0;
189 
190  // Saturation moisture fractions
191  erf_qsatw(tabs, pres, qsatw);
192  erf_qsati(tabs, pres, qsati);
193  erf_dtqsatw(tabs, pres, dqsatw);
194  erf_dtqsati(tabs, pres, dqsati);
195 
196  // Cloud ice not permitted (condensation & fusion)
197  if(tabs >= tbgmax) {
198  omn = 1.0;
199  }
200  // Cloud water not permitted (sublimation & fusion)
201  else if(tabs <= tbgmin) {
202  omn = 0.0;
203  lstarw = fac_sub;
204  }
205  // Mixed cloud phase (condensation & fusion)
206  else {
207  omn = an*tabs-bn;
208  domn = an;
209  }
210 
211  // Linear combination of each component
212  qsat = omn * qsatw + (1.0-omn) * qsati;
213  dqsat = omn * dqsatw + (1.0-omn) * dqsati
214  + domn * qsatw - domn * qsati;
215  lstar = omn * lstarw + (1.0-omn) * lstari;
216  dlstar = domn * lstarw - domn * lstari;
217 
218  // Function for root finding:
219  // 0 = -T_new + T_old + L_eff/C_p * (qv - qsat)
220  fff = -tabs + tabs_array(i,j,k) + lstar*(qv_array(i,j,k) - qsat);
221 
222  // Derivative of function (T_new iterated on)
223  dfff = -1.0 + dlstar*(qv_array(i,j,k) - qsat) - lstar*dqsat;
224 
225  // Update the temperature
226  dtabs = -fff/dfff;
227  tabs += dtabs;
228 
229  // Update iteration
230  niter = niter+1;
231  } while(std::abs(dtabs) > tol && niter < 20);
232 
233  // Update qsat from last iteration (dq = dq/dt * dt)
234  qsat += dqsat*dtabs;
235 
236  // Changes in each component
237  delta_qv = qv_array(i,j,k) - qsat;
238  delta_qc = std::max(-qc_array(i,j,k), delta_qv * omn);
239  delta_qi = std::max(-qi_array(i,j,k), delta_qv * (1.0-omn));
240 
241  // Partition the change in non-precipitating q
242  qv_array(i,j,k) = qsat;
243  qc_array(i,j,k) += delta_qc;
244  qi_array(i,j,k) += delta_qi;
245  qn_array(i,j,k) = qc_array(i,j,k) + qi_array(i,j,k);
246  qt_array(i,j,k) = qv_array(i,j,k) + qn_array(i,j,k);
247 
248  // Return to temperature
249  return tabs;
250  }
251 
252  void
254  std::vector<int>& a_idx,
255  std::vector<std::string>& a_names) const override
256  {
257  a_idx.clear();
258  a_names.clear();
259 
260  // The ordering here needs to match that in
261  // MicVarMap = {MicVar_Morr::nc, MicVar_Morr::nr, MicVar_Morr::ni, MicVar_Morr::ns, MicVar_Morr::ng,
262  // MicVar_Morr::rain_accum, MicVar_Morr::snow_accum, MicVar_Morr::graup_accum};
263  //
264  a_idx.push_back(0); a_names.push_back("Nc");
265  a_idx.push_back(1); a_names.push_back("Nr");
266  a_idx.push_back(2); a_names.push_back("Ni");
267  a_idx.push_back(3); a_names.push_back("Ns");
268  a_idx.push_back(4); a_names.push_back("Ng");
269  a_idx.push_back(5); a_names.push_back("RainAccum");
270  a_idx.push_back(6); a_names.push_back("SnowAccum");
271  a_idx.push_back(7); a_names.push_back("GraupAccum");
272  }
273 
274  void
275  GetPlotVarNames (amrex::Vector<std::string>& a_vec) const override;
276 
277  void
278  GetPlotVar (const std::string& a_name,
279  amrex::MultiFab& a_mf) const override;
280 
281 private:
282  // Number of qmoist variables (Ncloud, Nrain, Nice, Nsnow, Ngraup, rain_accum, snow_accum, graup_accum)
283  int m_qmoist_size = 8;
284 
285  // Number of qstate variables
287 
288  // CFL MAX for vertical advection
289  static constexpr amrex::Real CFL_MAX = 0.5;
290 
291  // MicVar map (Qmoist indices -> MicVar enum)
292  amrex::Vector<int> MicVarMap;
293 
294  // geometry
295  amrex::Geometry m_geom;
296 
297  // Nudging + set width
298  int m_real_width{0};
299 
300  // number of vertical levels
301  int nlev, zlo, zhi;
302 
303  // plane average axis
304  int m_axis;
305 
306  // constants
311  bool m_do_cond;
312 
313  // Minimum dz at this level
315 
316  // Pointer to terrain data
317  amrex::MultiFab* m_z_phys_nd;
318  amrex::MultiFab* m_detJ_cc;
319 
320  // independent variables
321  amrex::Array<FabPtr, MicVar_Morr::NumVars> mic_fab_vars;
322 
323  // microphysics parameters/coefficients
324  amrex::TableData<amrex::Real, 1> accrrc;
325  amrex::TableData<amrex::Real, 1> accrsi;
326  amrex::TableData<amrex::Real, 1> accrsc;
327  amrex::TableData<amrex::Real, 1> coefice;
328  amrex::TableData<amrex::Real, 1> evaps1;
329  amrex::TableData<amrex::Real, 1> evaps2;
330  amrex::TableData<amrex::Real, 1> accrgi;
331  amrex::TableData<amrex::Real, 1> accrgc;
332  amrex::TableData<amrex::Real, 1> evapg1;
333  amrex::TableData<amrex::Real, 1> evapg2;
334  amrex::TableData<amrex::Real, 1> evapr1;
335  amrex::TableData<amrex::Real, 1> evapr2;
336 
337  // vertical plane average data
338  amrex::TableData<amrex::Real, 1> rho1d;
339  amrex::TableData<amrex::Real, 1> pres1d;
340  amrex::TableData<amrex::Real, 1> tabs1d;
341  amrex::TableData<amrex::Real, 1> qt1d;
342  amrex::TableData<amrex::Real, 1> qv1d;
343  amrex::TableData<amrex::Real, 1> qn1d;
344 };
345 #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_Morrison.H:58
virtual ~Morrison()=default
amrex::Geometry m_geom
Definition: ERF_Morrison.H:295
Morrison()
Definition: ERF_Morrison.H:64
void Advance(const amrex::Real &dt_advance, const SolverChoice &sc) override
Definition: ERF_AdvanceMorrison.cpp:503
amrex::TableData< amrex::Real, 1 > evapg2
Definition: ERF_Morrison.H:333
void GetPlotVarNames(amrex::Vector< std::string > &a_vec) const override
Populate a vector with names of all available Morrison plot variables.
Definition: ERF_Morrison_Plot.cpp:101
void Qmoist_Restart_Vars(const SolverChoice &, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
Definition: ERF_Morrison.H:253
amrex::TableData< amrex::Real, 1 > qn1d
Definition: ERF_Morrison.H:343
int Qstate_Moist_Size() override
Definition: ERF_Morrison.H:137
int n_qstate_moist_size
Definition: ERF_Morrison.H:286
void Update_Micro_Vars(amrex::MultiFab &cons_in) override
Definition: ERF_Morrison.H:106
amrex::TableData< amrex::Real, 1 > pres1d
Definition: ERF_Morrison.H:339
amrex::TableData< amrex::Real, 1 > accrgi
Definition: ERF_Morrison.H:330
static constexpr amrex::Real CFL_MAX
Definition: ERF_Morrison.H:289
amrex::TableData< amrex::Real, 1 > accrgc
Definition: ERF_Morrison.H:331
void Set_dzmin(const amrex::Real dz_min) override
Definition: ERF_Morrison.H:92
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitMorrison.cpp:95
amrex::Real m_dzmin
Definition: ERF_Morrison.H:314
amrex::TableData< amrex::Real, 1 > tabs1d
Definition: ERF_Morrison.H:340
amrex::Vector< int > MicVarMap
Definition: ERF_Morrison.H:292
amrex::MultiFab * Qmoist_Ptr(const int &varIdx) override
Definition: ERF_Morrison.H:124
amrex::MultiFab * m_z_phys_nd
Definition: ERF_Morrison.H:317
amrex::TableData< amrex::Real, 1 > qt1d
Definition: ERF_Morrison.H:341
int nlev
Definition: ERF_Morrison.H:301
void Compute_Coefficients()
Definition: ERF_InitMorrison.cpp:149
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateMorrison.cpp:17
void Define(SolverChoice &sc) override
Definition: ERF_Morrison.H:71
amrex::MultiFab * m_detJ_cc
Definition: ERF_Morrison.H:318
amrex::TableData< amrex::Real, 1 > qv1d
Definition: ERF_Morrison.H:342
amrex::Real m_fac_cond
Definition: ERF_Morrison.H:307
amrex::Real m_fac_sub
Definition: ERF_Morrison.H:309
void Set_RealWidth(const int real_width) override
Definition: ERF_Morrison.H:140
amrex::TableData< amrex::Real, 1 > coefice
Definition: ERF_Morrison.H:327
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:145
int m_axis
Definition: ERF_Morrison.H:304
amrex::TableData< amrex::Real, 1 > accrsi
Definition: ERF_Morrison.H:325
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:134
int m_real_width
Definition: ERF_Morrison.H:298
amrex::TableData< amrex::Real, 1 > accrrc
Definition: ERF_Morrison.H:324
void GetPlotVar(const std::string &a_name, amrex::MultiFab &a_mf) const override
Extract a Morrison microphysics variable for plotting.
Definition: ERF_Morrison_Plot.cpp:130
amrex::TableData< amrex::Real, 1 > evapg1
Definition: ERF_Morrison.H:332
int m_qmoist_size
Definition: ERF_Morrison.H:283
amrex::TableData< amrex::Real, 1 > evapr1
Definition: ERF_Morrison.H:334
int zlo
Definition: ERF_Morrison.H:301
amrex::TableData< amrex::Real, 1 > evaps1
Definition: ERF_Morrison.H:328
amrex::TableData< amrex::Real, 1 > evaps2
Definition: ERF_Morrison.H:329
amrex::Real m_rdOcp
Definition: ERF_Morrison.H:310
amrex::Array< FabPtr, MicVar_Morr::NumVars > mic_fab_vars
Definition: ERF_Morrison.H:321
amrex::TableData< amrex::Real, 1 > rho1d
Definition: ERF_Morrison.H:338
std::shared_ptr< amrex::MultiFab > FabPtr
Definition: ERF_Morrison.H:60
int zhi
Definition: ERF_Morrison.H:301
amrex::TableData< amrex::Real, 1 > evapr2
Definition: ERF_Morrison.H:335
void Update_State_Vars(amrex::MultiFab &cons_in) override
Definition: ERF_Morrison.H:113
amrex::TableData< amrex::Real, 1 > accrsc
Definition: ERF_Morrison.H:326
bool m_do_cond
Definition: ERF_Morrison.H:311
amrex::Real m_fac_fus
Definition: ERF_Morrison.H:308
Definition: ERF_NullMoist.H:8
Definition: ERF_Morrison.H:24
@ qv
Definition: ERF_Morrison.H:34
@ ng
Definition: ERF_Morrison.H:48
@ nc
Definition: ERF_Morrison.H:44
@ qpg
Definition: ERF_Morrison.H:41
@ pres
Definition: ERF_Morrison.H:30
@ nr
Definition: ERF_Morrison.H:45
@ qcl
Definition: ERF_Morrison.H:35
@ tabs
Definition: ERF_Morrison.H:29
@ theta
Definition: ERF_Morrison.H:28
@ qp
Definition: ERF_Morrison.H:38
@ ni
Definition: ERF_Morrison.H:46
@ ns
Definition: ERF_Morrison.H:47
@ omega
Definition: ERF_Morrison.H:53
@ qps
Definition: ERF_Morrison.H:40
@ qn
Definition: ERF_Morrison.H:33
@ graup_accum
Definition: ERF_Morrison.H:52
@ rho
Definition: ERF_Morrison.H:27
@ qpr
Definition: ERF_Morrison.H:39
@ qci
Definition: ERF_Morrison.H:36
@ NumVars
Definition: ERF_Morrison.H:54
@ rain_accum
Definition: ERF_Morrison.H:50
@ snow_accum
Definition: ERF_Morrison.H:51
@ qt
Definition: ERF_Morrison.H:32
Definition: ERF_DataStruct.H:129
amrex::Real rdOcp
Definition: ERF_DataStruct.H:1124
amrex::Real c_p
Definition: ERF_DataStruct.H:1123
bool use_shoc
Definition: ERF_DataStruct.H:1155
int ave_plane
Definition: ERF_DataStruct.H:1186