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 condensate (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_axis = sc.ave_plane;
74  m_rdOcp = sc.rdOcp;
75  m_do_cond = (!sc.use_shoc);
76  }
77 
78  // init
79  void
80  Init (const amrex::MultiFab& cons_in,
81  const amrex::BoxArray& grids,
82  const amrex::Geometry& geom,
83  const amrex::Real& dt_advance,
84  std::unique_ptr<amrex::MultiFab>& z_phys_nd,
85  std::unique_ptr<amrex::MultiFab>& detJ_cc) override;
86 
87  // Import minimum dz at this level
88  void
89  Set_dzmin (const amrex::Real dz_min) override
90  {
91  m_dzmin = dz_min;
92  }
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,
111  const amrex::MultiFab& /*z_phys_nd*/) override
112  {
113  this->Copy_Micro_to_State(cons_in);
114  }
115 
116  // wrapper to do all the updating
117  void
118  Advance (const amrex::Real& dt_advance,
119  const SolverChoice& sc) override;
120 
121  amrex::MultiFab*
122  Qmoist_Ptr (const int& varIdx) override
123  {
125  return mic_fab_vars[MicVarMap[varIdx]].get();
126  }
127 
128  void
130 
131  int
132  Qmoist_Size () override { return Morrison::m_qmoist_size; }
133 
134  int
136 
137  void
138  Set_RealWidth (const int real_width) override { m_real_width = real_width; }
139 
140  AMREX_GPU_HOST_DEVICE
141  AMREX_FORCE_INLINE
142  static amrex::Real
143  NewtonIterSat (int& i, int& j, int& k,
144  const amrex::Real& fac_cond,
145  const amrex::Real& fac_fus,
146  const amrex::Real& fac_sub,
147  const amrex::Real& an,
148  const amrex::Real& bn,
149  const amrex::Array4<amrex::Real>& tabs_array,
150  const amrex::Array4<amrex::Real>& pres_array,
151  const amrex::Array4<amrex::Real>& qv_array,
152  const amrex::Array4<amrex::Real>& qc_array,
153  const amrex::Array4<amrex::Real>& qi_array,
154  const amrex::Array4<amrex::Real>& qn_array,
155  const amrex::Array4<amrex::Real>& qt_array)
156  {
157  // Solution tolerance
158  amrex::Real tol = amrex::Real(1.0e-4);
159 
160  // Saturation moisture fractions
161  amrex::Real omn, domn;
162  amrex::Real qsat, dqsat;
163  amrex::Real qsatw, dqsatw;
164  amrex::Real qsati, dqsati;
165 
166  // Newton iteration vars
167  int niter;
168  amrex::Real fff, dfff, dtabs;
169  amrex::Real lstar, dlstar;
170  amrex::Real lstarw, lstari;
171  amrex::Real delta_qv, delta_qc, delta_qi;
172 
173  // Initial guess for temperature & pressure
174  amrex::Real tabs = tabs_array(i,j,k);
175  amrex::Real pres = pres_array(i,j,k);
176 
177  niter = 0;
178  dtabs = 1;
179  //==================================================
180  // Newton iteration to qv=qsat (cloud phase only)
181  //==================================================
182  do {
183  // Latent heats and their derivatives wrt to T
184  lstarw = fac_cond;
185  lstari = fac_fus;
186  domn = zero;
187 
188  // Saturation moisture fractions
189  erf_qsatw(tabs, pres, qsatw);
190  erf_qsati(tabs, pres, qsati);
191  erf_dtqsatw(tabs, pres, dqsatw);
192  erf_dtqsati(tabs, pres, dqsati);
193 
194  // Cloud ice not permitted (condensation & fusion)
195  if(tabs >= tbgmax) {
196  omn = one;
197  }
198  // Cloud water not permitted (sublimation & fusion)
199  else if(tabs <= tbgmin) {
200  omn = zero;
201  lstarw = fac_sub;
202  }
203  // Mixed cloud phase (condensation & fusion)
204  else {
205  omn = an*tabs-bn;
206  domn = an;
207  }
208 
209  // Linear combination of each component
210  qsat = omn * qsatw + (one-omn) * qsati;
211  dqsat = omn * dqsatw + (one-omn) * dqsati
212  + domn * qsatw - domn * qsati;
213  lstar = omn * lstarw + (one-omn) * lstari;
214  dlstar = domn * lstarw - domn * lstari;
215 
216  // Function for root finding:
217  // 0 = -T_new + T_old + L_eff/C_p * (qv - qsat)
218  fff = -tabs + tabs_array(i,j,k) + lstar*(qv_array(i,j,k) - qsat);
219 
220  // Derivative of function (T_new iterated on)
221  dfff = -one + dlstar*(qv_array(i,j,k) - qsat) - lstar*dqsat;
222 
223  // Update the temperature
224  dtabs = -fff/dfff;
225  tabs += dtabs;
226 
227  // Update iteration
228  niter = niter+1;
229  } while(std::abs(dtabs) > tol && niter < 20);
230 
231  // Update qsat from last iteration (dq = dq/dt * dt)
232  qsat += dqsat*dtabs;
233 
234  // Changes in each component
235  delta_qv = qv_array(i,j,k) - qsat;
236  delta_qc = std::max(-qc_array(i,j,k), delta_qv * omn);
237  delta_qi = std::max(-qi_array(i,j,k), delta_qv * (one-omn));
238 
239  // Partition the change in non-precipitating q
240  qv_array(i,j,k) = qsat;
241  qc_array(i,j,k) += delta_qc;
242  qi_array(i,j,k) += delta_qi;
243  qn_array(i,j,k) = qc_array(i,j,k) + qi_array(i,j,k);
244  qt_array(i,j,k) = qv_array(i,j,k) + qn_array(i,j,k);
245 
246  // Return to temperature
247  return tabs;
248  }
249 
250  void
252  std::vector<int>& a_idx,
253  std::vector<std::string>& a_names) const override
254  {
255  a_idx.clear();
256  a_names.clear();
257 
258  // The ordering here needs to match that in
259  // MicVarMap = {MicVar_Morr::nc, MicVar_Morr::nr, MicVar_Morr::ni, MicVar_Morr::ns, MicVar_Morr::ng,
260  // MicVar_Morr::rain_accum, MicVar_Morr::snow_accum, MicVar_Morr::graup_accum};
261  //
262  a_idx.push_back(0); a_names.push_back("Nc");
263  a_idx.push_back(1); a_names.push_back("Nr");
264  a_idx.push_back(2); a_names.push_back("Ni");
265  a_idx.push_back(3); a_names.push_back("Ns");
266  a_idx.push_back(4); a_names.push_back("Ng");
267  a_idx.push_back(5); a_names.push_back("RainAccum");
268  a_idx.push_back(6); a_names.push_back("SnowAccum");
269  a_idx.push_back(7); a_names.push_back("GraupAccum");
270  }
271 
272  void
273  GetPlotVarNames (amrex::Vector<std::string>& a_vec) const override;
274 
275  void
276  GetPlotVar (const std::string& a_name,
277  amrex::MultiFab& a_mf) const override;
278 
279 private:
280  // Number of qmoist variables (Ncloud, Nrain, Nice, Nsnow, Ngraup, rain_accum, snow_accum, graup_accum)
281  int m_qmoist_size = 8;
282 
283  // Number of qstate variables
285 
286  // CFL MAX for vertical advection
287  static constexpr amrex::Real CFL_MAX = myhalf;
288 
289  // MicVar map (Qmoist indices -> MicVar enum)
290  amrex::Vector<int> MicVarMap;
291 
292  // geometry
293  amrex::Geometry m_geom;
294 
295  // Nudging + set width
296  int m_real_width{0};
297 
298  // number of vertical levels
299  int nlev, zlo, zhi;
300 
301  // plane average axis
302  int m_axis;
303 
304  // constants
306  bool m_do_cond;
307 
308  // Minimum dz at this level
310 
311  // Pointer to terrain data
312  amrex::MultiFab* m_z_phys_nd;
313  amrex::MultiFab* m_detJ_cc;
314 
315  // independent variables
316  amrex::Array<FabPtr, MicVar_Morr::NumVars> mic_fab_vars;
317 
318  // microphysics parameters/coefficients
319  amrex::TableData<amrex::Real, 1> accrrc;
320  amrex::TableData<amrex::Real, 1> accrsi;
321  amrex::TableData<amrex::Real, 1> accrsc;
322  amrex::TableData<amrex::Real, 1> coefice;
323  amrex::TableData<amrex::Real, 1> evaps1;
324  amrex::TableData<amrex::Real, 1> evaps2;
325  amrex::TableData<amrex::Real, 1> accrgi;
326  amrex::TableData<amrex::Real, 1> accrgc;
327  amrex::TableData<amrex::Real, 1> evapg1;
328  amrex::TableData<amrex::Real, 1> evapg2;
329  amrex::TableData<amrex::Real, 1> evapr1;
330  amrex::TableData<amrex::Real, 1> evapr2;
331 
332  // vertical plane average data
333  amrex::TableData<amrex::Real, 1> rho1d;
334  amrex::TableData<amrex::Real, 1> pres1d;
335  amrex::TableData<amrex::Real, 1> tabs1d;
336  amrex::TableData<amrex::Real, 1> qt1d;
337  amrex::TableData<amrex::Real, 1> qv1d;
338  amrex::TableData<amrex::Real, 1> qn1d;
339 };
340 #endif
constexpr amrex::Real tbgmax
Definition: ERF_Constants.H:43
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
constexpr amrex::Real tbgmin
Definition: ERF_Constants.H:42
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsatw(amrex::Real t, amrex::Real p, amrex::Real &dtqsatw)
Definition: ERF_MicrophysicsUtils.H:186
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsati(amrex::Real t, amrex::Real p, amrex::Real &dtqsati)
Definition: ERF_MicrophysicsUtils.H:178
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:171
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsati(amrex::Real t, amrex::Real p, amrex::Real &qsati)
Definition: ERF_MicrophysicsUtils.H:160
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF_Morrison.H:58
virtual ~Morrison()=default
amrex::Geometry m_geom
Definition: ERF_Morrison.H:293
Morrison()
Definition: ERF_Morrison.H:64
void Advance(const amrex::Real &dt_advance, const SolverChoice &sc) override
Definition: ERF_AdvanceMorrison.cpp:500
amrex::TableData< amrex::Real, 1 > evapg2
Definition: ERF_Morrison.H:328
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:251
amrex::TableData< amrex::Real, 1 > qn1d
Definition: ERF_Morrison.H:338
int Qstate_Moist_Size() override
Definition: ERF_Morrison.H:135
int n_qstate_moist_size
Definition: ERF_Morrison.H:284
void Update_Micro_Vars(amrex::MultiFab &cons_in) override
Definition: ERF_Morrison.H:103
amrex::TableData< amrex::Real, 1 > pres1d
Definition: ERF_Morrison.H:334
amrex::TableData< amrex::Real, 1 > accrgi
Definition: ERF_Morrison.H:325
void Update_State_Vars(amrex::MultiFab &cons_in, const amrex::MultiFab &) override
Definition: ERF_Morrison.H:110
static constexpr amrex::Real CFL_MAX
Definition: ERF_Morrison.H:287
amrex::TableData< amrex::Real, 1 > accrgc
Definition: ERF_Morrison.H:326
void Set_dzmin(const amrex::Real dz_min) override
Definition: ERF_Morrison.H:89
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitMorrison.cpp:90
amrex::Real m_dzmin
Definition: ERF_Morrison.H:309
amrex::TableData< amrex::Real, 1 > tabs1d
Definition: ERF_Morrison.H:335
amrex::Vector< int > MicVarMap
Definition: ERF_Morrison.H:290
amrex::MultiFab * Qmoist_Ptr(const int &varIdx) override
Definition: ERF_Morrison.H:122
amrex::MultiFab * m_z_phys_nd
Definition: ERF_Morrison.H:312
amrex::TableData< amrex::Real, 1 > qt1d
Definition: ERF_Morrison.H:336
int nlev
Definition: ERF_Morrison.H:299
void Compute_Coefficients()
Definition: ERF_InitMorrison.cpp:144
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:313
amrex::TableData< amrex::Real, 1 > qv1d
Definition: ERF_Morrison.H:337
void Set_RealWidth(const int real_width) override
Definition: ERF_Morrison.H:138
amrex::TableData< amrex::Real, 1 > coefice
Definition: ERF_Morrison.H:322
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:143
int m_axis
Definition: ERF_Morrison.H:302
amrex::TableData< amrex::Real, 1 > accrsi
Definition: ERF_Morrison.H:320
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:132
int m_real_width
Definition: ERF_Morrison.H:296
amrex::TableData< amrex::Real, 1 > accrrc
Definition: ERF_Morrison.H:319
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:327
int m_qmoist_size
Definition: ERF_Morrison.H:281
amrex::TableData< amrex::Real, 1 > evapr1
Definition: ERF_Morrison.H:329
int zlo
Definition: ERF_Morrison.H:299
amrex::TableData< amrex::Real, 1 > evaps1
Definition: ERF_Morrison.H:323
amrex::TableData< amrex::Real, 1 > evaps2
Definition: ERF_Morrison.H:324
amrex::Real m_rdOcp
Definition: ERF_Morrison.H:305
amrex::Array< FabPtr, MicVar_Morr::NumVars > mic_fab_vars
Definition: ERF_Morrison.H:316
amrex::TableData< amrex::Real, 1 > rho1d
Definition: ERF_Morrison.H:333
std::shared_ptr< amrex::MultiFab > FabPtr
Definition: ERF_Morrison.H:60
int zhi
Definition: ERF_Morrison.H:299
amrex::TableData< amrex::Real, 1 > evapr2
Definition: ERF_Morrison.H:330
amrex::TableData< amrex::Real, 1 > accrsc
Definition: ERF_Morrison.H:321
bool m_do_cond
Definition: ERF_Morrison.H:306
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:130
amrex::Real rdOcp
Definition: ERF_DataStruct.H:1169
bool use_shoc
Definition: ERF_DataStruct.H:1200
int ave_plane
Definition: ERF_DataStruct.H:1263