ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_SuperDropletsMoist.H
Go to the documentation of this file.
1 #ifndef SUPERDROPLETSMOIST_H
2 #define SUPERDROPLETSMOIST_H
3 
4 #ifdef ERF_USE_PARTICLES
5 
6 #include <limits>
7 #include <string>
8 #include <AMReX_Enum.H>
9 #include <AMReX_Geometry.H>
10 #include <AMReX_iMultiFab.H>
11 #include <AMReX_MultiFabUtil.H>
13 #include "ERF_SuperDropletPC.H"
14 
15 namespace MicVar_SD {
16  enum {
17  rho = 0, /*!< Density */
18  theta, /*!< Potential temperature */
19  temperature, /*!< Temperature */
20  pressure, /*!< Pressure */
21  q_t, /*!< Total (vapour + cloud) */
22  q_v, /*!< Water vapour */
23  q_c, /*!< Cloud condensate */
24  dqcdt, /*!< Condensation rate */
25  q_r, /*!< Rain */
26  rh, /*!< relative humidity */
27  rain_accum, /*!< rain accumulation (mm) */
28  NumVars
29  };
30 }
31 
32 namespace MicVar_SD_Species {
33  enum {
34  q_t = 0, /*!< Total density ratio */
35  q_v, /*!< Vapour density ratio */
36  q_c, /*!< Condensate density ratio */
37  sr, /*!< Saturation ratio */
38  accum, /*!< Ground accumulation */
39  NumVars
40  };
41 }
42 
43 namespace MicVar_SD_Aerosols {
44  enum {
45  accum = 0, /*!< Ground accumulation */
46  NumVars
47  };
48 }
49 
50 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
51 static int s_qt_idx (const int a_i /*!< species index */ )
52 {
53  AMREX_ALWAYS_ASSERT(a_i > 0);
54  return MicVar_SD::NumVars
57 }
58 
59 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
60 static int s_qv_idx (const int a_i /*!< species index */ )
61 {
62  AMREX_ALWAYS_ASSERT(a_i > 0);
63  return MicVar_SD::NumVars
65  + MicVar_SD_Species::q_v;
66 }
67 
68 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
69 static int s_qc_idx (const int a_i /*!< species index */ )
70 {
71  AMREX_ALWAYS_ASSERT(a_i > 0);
72  return MicVar_SD::NumVars
74  + MicVar_SD_Species::q_c;
75 }
76 
77 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
78 static int s_sr_idx (const int a_i /*!< species index */ )
79 {
80  AMREX_ALWAYS_ASSERT(a_i > 0);
81  return MicVar_SD::NumVars
83  + MicVar_SD_Species::sr;
84 }
85 
86 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
87 static int s_accum_idx (const int a_i /*!< species index */ )
88 {
89  AMREX_ALWAYS_ASSERT(a_i > 0);
90  return MicVar_SD::NumVars
92  + MicVar_SD_Species::accum;
93 }
94 
95 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
96 static int a_accum_idx (const int a_ns, /*!< number of species */
97  const int a_i /*!< aerosol index */ )
98 {
99  AMREX_ALWAYS_ASSERT(a_i >= 0);
100  return MicVar_SD::NumVars
101  + (a_ns-1)*MicVar_SD_Species::NumVars
103  + MicVar_SD_Aerosols::accum;
104 }
105 
106 /*! \brief List of super-droplet moisture model initializations */
107 AMREX_ENUM(SDMoistInit,
108  uniform,
109  condensate_density
110 );
111 
112 AMREX_ENUM(SDMSimulationDim,
113  one_d_z, two_d_xz, two_d_yz, three_d
114 );
115 
116 class SuperDropletsMoist : public NullMoistLagrangian {
117 
118  using FabPtr = std::shared_ptr<amrex::MultiFab>;
119  using MFPtr = std::unique_ptr<amrex::MultiFab>;
120  using BCTypeArr = amrex::GpuArray<ERF_BC, AMREX_SPACEDIM*2>;
121 
122  public:
123 
124  /*! \brief Constructor */
125  SuperDropletsMoist()
126  {
127  m_qmoist_size = 7;
128  m_qstate_moist_size = 3; // qv, qc, qr
129  m_qstate_nonmoist_size = -1; // intentional
130  m_name = "super_droplets_moisture";
131  m_idx_w = -1;
132  m_num_species = 0;
133  m_num_aerosols = 0;
134  readInputs();
135  }
136 
137  /*! \brief Destructor */
138  virtual ~SuperDropletsMoist() = default;
139 
140  /*! \brief Define super-droplet moisture model parameters */
141  virtual void Define (SolverChoice&) override;
142 
143  /*! \brief Initialize super-droplet moisture model data structures */
144  virtual void Init ( const amrex::MultiFab&,
145  const amrex::BoxArray&,
146  const amrex::Geometry&,
147  const amrex::Real&,
148  MFPtr&,
149  MFPtr& ) override;
150 
151  /*! \brief Initialize particles at AMR level */
152  virtual void InitParticles ( const int, MFPtr& ) override;
153 
154  /*! \brief Restart particles */
155  virtual void RestartParticles ( amrex::ParGDBBase*, const std::string& ) override;
156 
157  /*! \brief finish initializations */
158  virtual void FinishInit(const int&,
159  amrex::MultiFab&,
160  const amrex::Vector<MFPtr>&) override;
161 
162  void Advance ( const amrex::Real&, const SolverChoice& ) override
163  {
164  amrex::Abort("Do not use this advance() for Lagrangian microphysics.");
165  }
166 
167  /*! \brief Set the current AMR level being processed */
168  void SetCurrentLevel (const int& a_lev) override { m_current_lev = a_lev; }
169 
170  /*! \brief Get the diagnostics interval */
171  int getDiagnosticsInterval () const override { return m_diagnostics_iter; }
172 
173  /*! \brief Ensure m_mic_fab_vars[lev] is allocated for the given level */
174  void EnsureMicFabVars (const int a_lev, const amrex::MultiFab& a_cons_vars);
175 
176  /*! \brief Build a coarse-cell mask: 1 where exposed, 0 where covered by a finer level */
177  amrex::iMultiFab buildFineMask (const int a_lev) const;
178 
179  /*! \brief Initialize m_mic_fab_vars for a non-zero AMR level */
180  void InitLevel (const int a_lev, const amrex::MultiFab& a_cons_vars) override;
181 
182  /*! \brief Advance for one timestep
183  * \param[in] dt Timestep
184  * \param[in] iteration Current iteration number
185  * \param[in] time Current simulation time
186  * \param[in,out] cons_vars Array of conserved variables
187  * \param[in] mf_array Array of additional multifabs
188  * \param[in] bc_arr Array of boundary conditions
189  */
190  virtual void Advance (const amrex::Real& dt,
191  const int& iteration,
192  const amrex::Real& time,
193  amrex::Vector<amrex::Vector<amrex::MultiFab>>& cons_vars,
194  const amrex::Vector<MFPtr>& mf_array,
195  const BCTypeArr& bc_arr) override;
196 
197  /*! \brief Update microphysics variables
198  * \param[in] a_cons_vars Conservative variables multifab
199  */
200  virtual void Update_Micro_Vars (amrex::MultiFab& a_cons_vars) override;
201 
202  /*! \brief Update state variables from microphysics variables
203  * \param[in,out] a_cons_vars Conservative variables multifab to update
204  * \param[in] a_z_phys_nd Terrain heights (nodal)
205  */
206  virtual void Update_State_Vars (amrex::MultiFab& a_cons_vars,
207  const amrex::MultiFab& a_z_phys_nd) override;
208 
209  /*! \brief Average down moisture multifabs from finest_level down to 0 */
210  virtual void AverageDownMicroVars (const int finest_level) override;
211 
212  /*! \brief Copy moisture model vars from state to member multifabs
213  * \param[in] a_state State multifab to copy from
214  */
215  virtual void Copy_State_to_Micro (const amrex::MultiFab& a_state) override;
216 
217  /*! \brief Copy moisture model vars to state from member multifabs
218  * \param[out] a_state State multifab to copy to
219  */
220  virtual void Copy_Micro_to_State (amrex::MultiFab& a_state) override;
221 
222  /*! \brief returns the moisture variable asked for */
223  inline virtual amrex::MultiFab* Qmoist_Ptr (const int& a_idx) override
224  {
225  AMREX_ALWAYS_ASSERT( a_idx < m_qmoist_size );
226  const int lev = m_current_lev;
227  if (lev < 0 || lev >= static_cast<int>(m_mic_fab_vars.size()) ||
228  m_mic_var_map[a_idx] < 0 || m_mic_var_map[a_idx] >= static_cast<int>(m_mic_fab_vars[lev].size())) {
229  return nullptr;
230  }
231  return m_mic_fab_vars[lev][m_mic_var_map[a_idx]].get();
232  }
233 
234  /*! \brief returns number of variables in this moisture model */
235  inline virtual int Qmoist_Size () override
236  {
237  return m_qmoist_size;
238  }
239 
240  /*! \brief returns number of moist conserved variables */
241  inline virtual int Qstate_Moist_Size () override
242  {
243  return m_qstate_moist_size;
244  }
245 
246  /*! \brief returns number of non-water species conserved variables */
247  inline virtual int Qstate_NonMoist_Size () override
248  {
249  AMREX_ALWAYS_ASSERT(m_qstate_nonmoist_size >= 0);
250  return m_qstate_nonmoist_size;
251  }
252 
253  /*! \brief returns pointer to super-droplets particle container */
254  inline virtual ERFPC* getParticleContainer() override
255  {
256  return m_super_droplets;
257  }
258 
259  inline virtual const std::string& getName() const override
260  {
261  return m_name;
262  }
263 
264  /*! \brief compute cloud/rain mixing ratio */
265  virtual void computeQcQrWater (const amrex::MultiFab& a_z_phys_nd);
266 
267  /*! \brief compute condensate for all non-water species */
268  virtual void computeQcSpecies (const amrex::MultiFab& a_z_phys_nd)
269  {
270  for (int is = 1; is < m_num_species; is++) { computeQcSpecies(is, a_z_phys_nd); }
271  }
272 
273  /*! \brief compute condensate mixing ratio for a non-water species */
274  virtual void computeQcSpecies (const int a_i, const amrex::MultiFab& a_z_phys_nd);
275 
276  /*! \brief compute condensate mixing ratio for a species */
277  virtual void computeQc (const int a_i, const amrex::MultiFab& a_z_phys_nd)
278  {
279  if (a_i == m_idx_w) { computeQcQrWater(a_z_phys_nd); }
280  else { computeQcSpecies(a_i, a_z_phys_nd); }
281  }
282 
283  /*! \brief compute total water mixing ratio */
284  virtual void computeQtWater ();
285 
286  /*! \brief compute qt (total) for all non-water species */
287  virtual void computeQtSpecies ()
288  {
289  for (int is = 1; is < m_num_species; is++) { computeQtSpecies(is); }
290  }
291 
292  /*! \brief compute qt (total) for a non-water species */
293  virtual void computeQtSpecies (const int a_i);
294 
295  /*! \brief compute qt (total) for a species */
296  virtual void computeQt (const int a_i)
297  {
298  if (a_i == m_idx_w) { computeQtWater(); }
299  else { computeQtSpecies(a_i); }
300  }
301 
302  /*! \brief Compute rain accumulation */
303  virtual void rainAccumulation (const amrex::MultiFab& a_z_phys_nd);
304 
305  /*! \brief Compute non-water species accumulation */
306  virtual void speciesAccumulation (const amrex::MultiFab& a_z_phys_nd);
307 
308  /*! \brief Compute aerosol accumulation */
309  virtual void aerosolAccumulation (const amrex::MultiFab& a_z_phys_nd);
310 
311  /*! \brief Convert a multifab containing density of something to its mixing ratio */
312  void densityToRatio ( amrex::MultiFab&, const int a_comp = 0 );
313  /*! \brief Convert a multifab containing the mixing ratio of something to its density */
314  void ratioToDensity ( amrex::MultiFab&, const int a_comp = 0 );
315 
316  /*! \brief Compute phase changes (vapour <--> cloud water)
317  * \param[in] a_dt Timestep for phase change calculation
318  * \param[in] a_z Array of terrain heights
319  * \param[in] a_lev AMR level
320  */
321  virtual void phaseChange ( const amrex::Real& a_dt,
322  const amrex::Vector<MFPtr>& a_z,
323  const int a_lev);
324 
325  virtual void GetPlotVarNames (amrex::Vector<std::string>& a_names) const override
326  {
327  for (int v = 1; v < m_num_species; v++) {
328  a_names.push_back("qv_"+amrex::getEnumNameString(m_species[v]));
329  a_names.push_back("qc_"+amrex::getEnumNameString(m_species[v]));
330  a_names.push_back("qt_"+amrex::getEnumNameString(m_species[v]));
331  a_names.push_back("sat_ratio_"+amrex::getEnumNameString(m_species[v]));
332  a_names.push_back("accum_"+amrex::getEnumNameString(m_species[v]));
333  }
334  for (int v = 0; v < m_num_aerosols; v++) {
335  a_names.push_back("accum_"+amrex::getEnumNameString(m_aerosols[v]));
336  }
337  }
338 
339  virtual void GetPlotVar (const std::string& /*a_name*/,
340  amrex::MultiFab& /*a_mf*/) const override
341  {
342  amrex::Abort("SuperDropletsMoist::GetPlotVar() requires a level argument");
343  }
344 
345  virtual void GetPlotVar (const std::string& a_name,
346  amrex::MultiFab& a_mf,
347  const int a_lev) const override
348  {
349  a_mf.setVal(0.0);
350  AMREX_ASSERT(a_mf.nComp() >= 1);
351 
352  const int lev = a_lev;
353  if (lev < 0 || lev >= static_cast<int>(m_mic_fab_vars.size()) ||
354  m_mic_fab_vars[lev].empty()) {
355  return;
356  }
357 
358  const auto& lev_vec = m_mic_fab_vars[lev];
359 
360  // Helper to copy MultiFab if name matches
361  auto try_copy = [&](const std::string& prefix, int idx) -> bool {
362  if (a_name == prefix) {
363  if (idx >= 0 && idx < static_cast<int>(lev_vec.size()) && lev_vec[idx] &&
364  lev_vec[idx]->boxArray() == a_mf.boxArray() &&
365  lev_vec[idx]->DistributionMap() == a_mf.DistributionMap()) {
366  amrex::MultiFab::Copy(a_mf, *lev_vec[idx],
367  0, 0, 1, amrex::IntVect::TheZeroVector());
368  }
369  return true;
370  }
371  return false;
372  };
373 
374  // Species variables
375  for (int v = 1; v < m_num_species; v++) {
376  std::string sp_name = amrex::getEnumNameString(m_species[v]);
377  if (try_copy("qv_" + sp_name, s_qv_idx(v)) ||
378  try_copy("qc_" + sp_name, s_qc_idx(v)) ||
379  try_copy("qt_" + sp_name, s_qt_idx(v)) ||
380  try_copy("sat_ratio_" + sp_name, s_sr_idx(v)) ||
381  try_copy("accum_" + sp_name, s_accum_idx(v))) {
382  return;
383  }
384  }
385 
386  // Aerosol variables
387  for (int v = 0; v < m_num_aerosols; v++) {
388  std::string ae_name = amrex::getEnumNameString(m_aerosols[v]);
389  if (try_copy("accum_" + ae_name, a_accum_idx(m_num_species, v))) { return; }
390  }
391 
392  amrex::Abort("SuperDropletsMoist::GetPlotVar() called with invalid name");
393  }
394 
395  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
396  int q_qv_idx (const int a_i /*!< species index */ )
397  {
398  AMREX_ALWAYS_ASSERT(a_i > 0);
399  return RhoQ1_comp + m_qstate_moist_size + 2*(a_i-1) + 0;
400  }
401 
402  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
403  int q_qc_idx (const int a_i /*!< species index */ )
404  {
405  AMREX_ALWAYS_ASSERT(a_i > 0);
406  return RhoQ1_comp + m_qstate_moist_size + 2*(a_i-1) + 1;
407  }
408 
409  protected:
410 
411  bool m_flag_phase_change; /*!< Enable/disable phase change */
412  bool m_flag_advection; /*!< Enable/disable advection */
413  bool m_flag_coalescence; /*!< Enable/disable coalescence */
414 
415  amrex::Real m_Cp; /*!< specific heat at constant pressure for dry air */
416 
417  amrex::Real m_r_rain; /*!< minimum radius to be considered rain */
418 
419  /*! let initial superdroplets "relax" to a physically-appropriate size? */
420  bool m_init_phase_change;
421  /*! time (in seconds) of initial relaxation */
422  amrex::Real m_init_phase_change_time;
423 
424  int m_diagnostics_iter; /*!< number of iterations between computing diagnostics */
425 
426  /*! initialization type */
427  SDMoistInit m_init_type;
428 
429  /*! name of model and its particle container */
430  std::string m_name;
431  /*! Geometry object */
432  amrex::Geometry m_geom;
433 
434  /*! number of microphysics variables */
435  int m_qmoist_size;
436  /*! number of water-related state variables for this moisture model */
437  int m_qstate_moist_size;
438  /*! number of non-water-related state variables for this moisture model */
439  int m_qstate_nonmoist_size;
440 
441  amrex::Real m_dt; /*!< timestep */
442  amrex::Vector<int> m_mic_var_map; /*!< moisture model variables map */
443 
444  /*! moisture model variables - per level for AMR support */
445  amrex::Vector<amrex::Vector<FabPtr>> m_mic_fab_vars;
446 
447  /*! current AMR level being processed */
448  mutable int m_current_lev = 0;
449 
450  /*! names of vapour/condensate species */
451  std::vector<Species::Name> m_species;
452  int m_num_species;
453 
454  /*! names of aerosols */
455  std::vector<Species::Name> m_aerosols;
456  int m_num_aerosols;
457 
458  /*! number of substeps for phase change */
459  int m_num_substeps_phase_change;
460 
461  /*! kinematic mode? */
462  bool m_kinematic_mode;
463 
464  /*! Dimensionality of simulation*/
465  SDMSimulationDim m_dimensionality;
466 
467  /*! recycle particles */
468  bool m_recycle_particles;
469 
470  /*! species index of water */
471  int m_idx_w;
472 
473  /*! particle container for super-droplets
474  * Owned by ERF::particleData which deletes it during teardown;
475  * raw pointer here to avoid double-free. */
476  SuperDropletPC* m_super_droplets;
477 
478  /*! \brief read inputs */
479  virtual void readInputs();
480 
481  private:
482 
483 };
484 
485 #endif
486 #endif
AMREX_ENUM(InitType, None, Input_Sounding, NCFile, WRFInput, Metgrid, Uniform, ConstantDensity, ConstantDensityLinearTheta, Isentropic, MoistBaseState, HindCast)
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
rho
Definition: ERF_InitCustomPert_Bubble.H:106
int m_num_species
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:28
std::vector< std::string > m_species
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:27
Real q_t
Definition: ERF_InitCustomPert_SquallLine.H:24
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int idx(int i, int j, int k, int nx, int ny)
Definition: ERF_InitForEnsemble.cpp:287
Contains the Lagrangian moisture model base class.
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ NumVars
Definition: ERF_MM5.H:21
@ theta
Definition: ERF_MM5.H:20
@ rain_accum
Definition: ERF_Kessler.H:33
Definition: ERF_DataStruct.H:130