ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_LagrangianMicrophysics.H
Go to the documentation of this file.
1 /*! @file ERF_LagrangianMicrophysics.H
2  * \brief Contains the Lagrangian microphysics class
3  */
4 
5 #ifndef ERF_LAGRANGIANMICROPHYSICS_H
6 #define ERF_LAGRANGIANMICROPHYSICS_H
7 
8 #ifdef ERF_USE_PARTICLES
9 
10 #include <utility>
11 #include <string>
12 
14 #include "ERF_SuperDropletsMoist.H"
15 #include "ERF_Microphysics.H"
16 
17 /* forward declaration */
18 class ERFPC;
19 
20 /*! \brief Eulerian microphysics interface
21  *
22  * One key difference from #EulerianMicrophysics is that only the base AMR
23  * level has the moisture model. Thus, for higher AMR levels, the microphysics
24  * interface need not do anything. The number of conserved state variables for
25  * moisture (RhoQ1, RhoQ2, ...) are the same for all AMR levels; however, for
26  * levels other than the base, they just contain and evolve zeros. */
27 class LagrangianMicrophysics : public Microphysics {
28 
29 public:
30 
31  /*! \brief Null constructor */
32  LagrangianMicrophysics () { }
33 
34  /*! \brief default destructor */
35  ~LagrangianMicrophysics () = default;
36 
37  /*! \brief Constructor: create the moisture model */
38  LagrangianMicrophysics (const int& /* nlev */, /*!< Number of AMR levels */
39  const MoistureType& a_model_type /*!< moisture model */ )
40  {
41  AMREX_ASSERT( Microphysics::modelType(a_model_type) == MoistureModelType::Lagrangian );
42  if (a_model_type == MoistureType::SuperDroplets) {
43  SetModel<SuperDropletsMoist>();
44  amrex::Print() << "Super-Droplets moisture model!\n";
45  } else {
46  amrex::Abort("Unknown Lagrangian moisture model!") ;
47  }
48  }
49 
50  /*! \brief Define the moisture model */
51  void Define (const int& lev, /*!< AMR level */
52  SolverChoice& sc /*!< Solver choice object */) override
53  {
54  if (lev > 0) return;
55  m_moist_model->Define(sc);
56  }
57 
58  /*! \brief Initialize the moisture model */
59  void Init (const int& lev, /*!< AMR level */
60  const amrex::MultiFab& cons_in, /*!< Conserved state variables */
61  const amrex::BoxArray& grids, /*!< Grids */
62  const amrex::Geometry& geom, /*!< Geometry */
63  const amrex::Real& dt_advance, /*!< Time step */
64  std::unique_ptr<amrex::MultiFab>& z_phys_nd,/*< Nodal z heights */
65  std::unique_ptr<amrex::MultiFab>& detJ_cc /*< CC Jacobian determinants */) override
66  {
67  if (lev == 0) {
68  m_moist_model->Init(cons_in, grids, geom, dt_advance,
69  z_phys_nd, detJ_cc);
70  } else {
71  m_moist_model->InitLevel(lev, cons_in);
72  }
73  }
74 
75  /*! \brief finish initializations steps that require flow variables */
76  void FinishInit (const int& lev, /*!< AMR level */
77  amrex::MultiFab& cons_in, /*!< Conserved state variable */
78  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& z_phys_nd /*!< terrain */) override
79  {
80  if (lev > 0) return;
81  m_moist_model->FinishInit( lev, cons_in, z_phys_nd );
82  }
83 
84  /*! \brief Advance the moisture model for one time step */
85  void Advance (const int& lev, /*!< AMR level */
86  const amrex::Real& dt_advance, /*!< Time step */
87  const int& iter, /*!< iteration number */
88  const amrex::Real& time, /*!< current time */
89  const SolverChoice& /*solverChoice*/, /*!< Solver choice object */
90  amrex::Vector<amrex::Vector<amrex::MultiFab>>& a_vars, /*!< Dycore state variables */
91  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& a_z/*!< terrain */,
92  const amrex::GpuArray<ERF_BC, AMREX_SPACEDIM*2>& a_phys_bc_types ) override
93  {
94  m_moist_model->SetCurrentLevel(lev);
95  m_moist_model->Advance(dt_advance, iter, time, a_vars, a_z, a_phys_bc_types);
96  }
97 
98  /*! \brief update microphysics variables from ERF state variables */
99  void Update_Micro_Vars_Lev (const int& lev, /*! AMR level */
100  amrex::MultiFab& cons_in /*!< Conserved state variables */) override
101  {
102  m_moist_model->SetCurrentLevel(lev);
103  m_moist_model->Update_Micro_Vars(cons_in);
104  }
105 
106  /*! \brief update ERF state variables from microphysics variables */
107  void Update_State_Vars_Lev (const int& lev, /*!< AMR level */
108  amrex::MultiFab& cons_in, /*!< Conserved state variables */
109  const amrex::MultiFab& z_phys_nd /*!< Nodal terrain heights */) override
110  {
111  m_moist_model->SetCurrentLevel(lev);
112  m_moist_model->Update_State_Vars(cons_in, z_phys_nd);
113  }
114 
115  /*! \brief average down moisture multifabs from finest_level to 0 */
116  void AverageDownMicroVars (const int finest_level) override
117  {
118  m_moist_model->AverageDownMicroVars(finest_level);
119  }
120 
121  /*! \brief get pointer to a moisture variable */
122  amrex::MultiFab* Get_Qmoist_Ptr (const int& lev, /*!< AMR level */
123  const int& varIdx /*!< moisture variable index */) override
124  {
125  m_moist_model->SetCurrentLevel(lev);
126  return m_moist_model->Qmoist_Ptr(varIdx);
127  }
128 
129  /*! \brief get the number of moisture model variables */
130  int Get_Qmoist_Size (const int& /*a_lev*/ /*!< AMR level */) override
131  {
132  return m_moist_model->Qmoist_Size();
133  }
134 
135  /*! \brief get the number of microphysics conserved moist (water-related) state variables */
136  int Get_Qstate_Moist_Size () override
137  {
138  return m_moist_model->Qstate_Moist_Size();
139  }
140 
141  /*! \brief get the number of microphysics conserved moist (water-related) state variables that are number concentrations */
142  int Get_Qstate_Moist_NumConc_Size () override
143  {
144  return m_moist_model->Qstate_Moist_NumConc_Size();
145  }
146 
147  /*! \brief get the number of microphysics conserved non-moist (non-water, i.e., other vapor/condensed species) state variables */
148  int Get_Qstate_NonMoist_Size () override
149  {
150  return m_moist_model->Qstate_NonMoist_Size();
151  }
152 
153  /*! \brief initialize particles in particle container at level a_lev */
154  inline void initParticles ( const int a_lev, /*!< AMR level */
155  std::unique_ptr<amrex::MultiFab>& z_phys_nd /*!< Nodal z heights at level a_lev */ )
156  {
157  m_moist_model->InitParticles(a_lev, z_phys_nd);
158  }
159 
160  /*! \brief restart particles in particle container */
161  inline void restartParticles ( amrex::ParGDBBase* a_gdb, const std::string& a_fname)
162  {
163  m_moist_model->RestartParticles(a_gdb, a_fname);
164  }
165 
166  /*! \brief get the particle container from the moisture model */
167  inline ERFPC* getParticleContainer () const
168  {
169  return m_moist_model->getParticleContainer();
170  }
171 
172  /*! \brief get the name of the moisture model's particle container */
173  inline const std::string& getName () const
174  {
175  return m_moist_model->getName();
176  }
177 
178  /*! \brief get the indices and names of moisture model variables for restart
179  at a given level */
180  void Get_Qmoist_Restart_Vars ( const int a_lev, /*!< level */
181  const SolverChoice& a_sc, /*!< Solver choice object */
182  std::vector<int>& a_idx, /*!< indices */
183  std::vector<std::string>& a_names /*!< names */ ) const override
184  {
185  if (a_lev == 0) {
186  m_moist_model->Qmoist_Restart_Vars( a_sc, a_idx, a_names );
187  } else {
188  a_idx.clear();
189  a_names.clear();
190  }
191  }
192 
193  /*! \brief Returns a list of additional plot variable names */
194  virtual void GetPlotVarNames (amrex::Vector<std::string>& a_vec) const override
195  {
196  m_moist_model->GetPlotVarNames(a_vec);
197  }
198 
199  /*! \brief Fills in a MultiFab for plotting */
200  virtual void GetPlotVar (const std::string& a_name,
201  amrex::MultiFab& a_mf,
202  const int a_lev) const override
203  {
204  m_moist_model->SetCurrentLevel(a_lev);
205  m_moist_model->GetPlotVar(a_name, a_mf, a_lev);
206  }
207 
208  /*! \brief get the diagnostics interval */
209  inline int getDiagnosticsInterval () const
210  {
211  return m_moist_model->getDiagnosticsInterval();
212  }
213 
214  /*! \brief Populates dz_min in micro model for precipitation */
215  virtual void Set_dzmin (const int /*a_lev*/,
216  const amrex::Real dz_min) const override
217  {
218  m_moist_model->Set_dzmin(dz_min);
219  }
220 
221 protected:
222 
223  /*! \brief Create and set the specified moisture model */
224  template<class NewMoistModel>
225  void SetModel ()
226  {
227  m_moist_model = std::make_unique<NewMoistModel>();
228  }
229 
230  std::unique_ptr<NullMoistLagrangian> m_moist_model; /*!< moisture model */
231 };
232 
233 #endif
234 #endif
Contains the base class for microphysics.
Contains the Lagrangian moisture model base class.
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Base class for microphysics interface.
Definition: ERF_Microphysics.H:15
virtual int Get_Qmoist_Size(const int &)=0
get the number of moisture model variables
virtual int Get_Qstate_NonMoist_Size()=0
get the number of microphysics conserved non-moist (non-water, i.e., other vapor/condensed species) s...
virtual void Define(const int &, SolverChoice &)=0
define the microphysics object
virtual void AverageDownMicroVars(const int)
average down moisture multifabs from finest_level to 0 (AMR coverage); no-op for non-Lagrangian model...
Definition: ERF_Microphysics.H:59
virtual int Get_Qstate_Moist_NumConc_Size()=0
get the number of microphysics conserved moist (water-related) state variables that are number concen...
virtual void Init(const int &, const amrex::MultiFab &, const amrex::BoxArray &, const amrex::Geometry &, const amrex::Real &, std::unique_ptr< amrex::MultiFab > &, std::unique_ptr< amrex::MultiFab > &)=0
initialize the microphysics object
virtual int Get_Qstate_Moist_Size()=0
get the number of microphysics conserved moist (water-related) state variables
virtual void Update_Micro_Vars_Lev(const int &, amrex::MultiFab &)=0
update microphysics variables from ERF state variables
virtual void Update_State_Vars_Lev(const int &, amrex::MultiFab &, const amrex::MultiFab &)=0
update ERF state variables from microphysics variables
virtual amrex::MultiFab * Get_Qmoist_Ptr(const int &, const int &)=0
get pointer to a moisture variable
virtual void GetPlotVarNames(amrex::Vector< std::string > &a_vec) const =0
Returns a list of additional plot variable names.
virtual void Get_Qmoist_Restart_Vars(int, const SolverChoice &, std::vector< int > &, std::vector< std::string > &) const =0
get the indices and names of moisture model variables for restart at a given level
virtual void Set_dzmin(const int lev, const amrex::Real dz_min) const =0
Import minimum dz at this level.
static MoistureModelType modelType(const MoistureType a_moisture_type)
query if a specified moisture model is Eulerian or Lagrangian
Definition: ERF_Microphysics.H:102
virtual void FinishInit(const int &, amrex::MultiFab &, const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &)=0
finish initializations steps that require flow variables
virtual void GetPlotVar(const std::string &a_name, amrex::MultiFab &a_mf, const int a_lev) const =0
Fills in a MultiFab for plotting.
virtual void Advance(const int &, const amrex::Real &, const int &, const amrex::Real &, const SolverChoice &, amrex::Vector< amrex::Vector< amrex::MultiFab >> &, const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &, const amrex::GpuArray< ERF_BC, AMREX_SPACEDIM *2 > &)=0
advance microphysics for one time step
Definition: ERF_DataStruct.H:141