ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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_Microphysics.H"
15 
16 /* forward declaration */
17 class ERFPC;
18 
19 /*! \brief Eulerian microphysics interface
20  *
21  * One key difference from #EulerianMicrophysics is that only the base AMR
22  * level has the moisture model. Thus, for higher AMR levels, the microphysics
23  * interface need not do anything. The number of conserved state variables for
24  * moisture (RhoQ1, RhoQ2, ...) are the same for all AMR levels; however, for
25  * levels other than the base, they just contain and evolve zeros. */
26 class LagrangianMicrophysics : public Microphysics {
27 
28 public:
29 
30  /*! \brief Null constructor */
31  LagrangianMicrophysics () { }
32 
33  /*! \brief default destructor */
34  ~LagrangianMicrophysics () = default;
35 
36  /*! \brief Constructor: create the moisture model */
37  LagrangianMicrophysics (const int& /* nlev */, /*!< Number of AMR levels */
38  const MoistureType& a_model_type /*!< moisture model */ )
39  {
40  AMREX_ASSERT( Microphysics::modelType(a_model_type) == MoistureModelType::Lagrangian );
41  amrex::ignore_unused(a_model_type);
42  amrex::Abort("No Lagrangian moisture model implemented yet!") ;
43  }
44 
45  /*! \brief Define the moisture model */
46  void Define (const int& lev, /*!< AMR level */
47  SolverChoice& sc /*!< Solver choice object */) override
48  {
49  if (lev > 0) return;
50  m_moist_model->Define(sc);
51  }
52 
53  /*! \brief Initialize the moisture model */
54  void Init (const int& lev, /*!< AMR level */
55  const amrex::MultiFab& cons_in, /*!< Conserved state variables */
56  const amrex::BoxArray& grids, /*!< Grids */
57  const amrex::Geometry& geom, /*!< Geometry */
58  const amrex::Real& dt_advance, /*!< Time step */
59  std::unique_ptr<amrex::MultiFab>& z_phys_nd,/*< Nodal z heights */
60  std::unique_ptr<amrex::MultiFab>& detJ_cc /*< CC Jacobian determinants */) override
61  {
62  if (lev > 0) return;
63  m_moist_model->Init(cons_in, grids, geom, dt_advance,
64  z_phys_nd, detJ_cc);
65  }
66 
67  /*! \brief Advance the moisture model for one time step */
68  void Advance (const int& lev, /*!< AMR level */
69  const amrex::Real& dt_advance, /*!< Time step */
70  const int& iter, /*!< iteration number */
71  const amrex::Real& time, /*!< current time */
72  const SolverChoice& /*solverChoice*/, /*!< Solver choice object */
73  amrex::Vector<amrex::Vector<amrex::MultiFab>>& a_vars, /*!< Dycore state variables */
74  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& a_z/*!< terrain */,
75  const amrex::GpuArray<ERF_BC, AMREX_SPACEDIM*2>& a_phys_bc_types ) override
76  {
77  if (lev > 0) return;
78  m_moist_model->Advance(dt_advance, iter, time, a_vars, a_z, a_phys_bc_types);
79  }
80 
81  /*! \brief update microphysics variables from ERF state variables */
82  void Update_Micro_Vars_Lev (const int& lev, /*! AMR level */
83  amrex::MultiFab& cons_in /*!< Conserved state variables */) override
84  {
85  if (lev > 0) return;
86  m_moist_model->Update_Micro_Vars(cons_in);
87  }
88 
89  /*! \brief update ERF state variables from microphysics variables */
90  void Update_State_Vars_Lev (const int& lev, /*!< AMR level */
91  amrex::MultiFab& cons_in /*!< Conserved state variables */) override
92  {
93  if (lev > 0) return;
94  m_moist_model->Update_State_Vars(cons_in);
95  }
96 
97  /*! \brief get pointer to a moisture variable */
98  amrex::MultiFab* Get_Qmoist_Ptr (const int& lev, /*!< AMR level */
99  const int& varIdx /*!< moisture variable index */) override
100  {
101  return (lev > 0 ? nullptr : m_moist_model->Qmoist_Ptr(varIdx));
102  }
103 
104  /*! \brief get the number of moisture model variables */
105  int Get_Qmoist_Size (const int& a_lev /*!< AMR level */) override
106  {
107  return (a_lev > 0 ? 0 : m_moist_model->Qmoist_Size());
108  }
109 
110  /*! \brief get the number of microphysics conserved moist (water-related) state variables */
111  int Get_Qstate_Moist_Size () override
112  {
113  return m_moist_model->Qstate_Moist_Size();
114  }
115 
116  /*! \brief get the number of microphysics conserved non-moist (non-water, i.e., other vapor/condensed species) state variables */
117  int Get_Qstate_NonMoist_Size () override
118  {
119  return m_moist_model->Qstate_NonMoist_Size();
120  }
121 
122  /*! \brief initialize particles in particle container */
123  inline void initParticles ( std::unique_ptr<amrex::MultiFab>& z_phys_nd /*!< Nodal z heights */ )
124  {
125  m_moist_model->InitParticles(z_phys_nd);
126  }
127 
128  /*! \brief restart particles in particle container */
129  inline void restartParticles ( amrex::ParGDBBase* a_gdb, const std::string& a_fname)
130  {
131  m_moist_model->RestartParticles(a_gdb, a_fname);
132  }
133 
134  /*! \brief get the particle container from the moisture model */
135  inline ERFPC* getParticleContainer () const
136  {
137  return m_moist_model->getParticleContainer();
138  }
139 
140  /*! \brief get the name of the moisture model's particle container */
141  inline const std::string& getName () const
142  {
143  return m_moist_model->getName();
144  }
145 
146  /*! \brief get the indices and names of moisture model variables for restart
147  at a given level */
148  void Get_Qmoist_Restart_Vars ( const int a_lev, /*!< level */
149  const SolverChoice& a_sc, /*!< Solver choice object */
150  std::vector<int>& a_idx, /*!< indices */
151  std::vector<std::string>& a_names /*!< names */ ) const override
152  {
153  if (a_lev == 0) {
154  m_moist_model->Qmoist_Restart_Vars( a_sc, a_idx, a_names );
155  } else {
156  a_idx.clear();
157  a_names.clear();
158  }
159  }
160 
161  /*! \brief Returns a list of additional plot variable names */
162  virtual void GetPlotVarNames (amrex::Vector<std::string>& a_vec) const override
163  {
164  m_moist_model->GetPlotVarNames(a_vec);
165  }
166 
167  /*! \brief Fills in a MultiFab for plotting */
168  virtual void GetPlotVar (const std::string& a_name,
169  amrex::MultiFab& a_mf,
170  const int /* a_lev */) const override
171  {
172  m_moist_model->GetPlotVar(a_name, a_mf);
173  }
174 
175 protected:
176 
177  /*! \brief Create and set the specified moisture model */
178  template<class NewMoistModel>
179  void SetModel ()
180  {
181  m_moist_model = std::make_unique<NewMoistModel>();
182  }
183 
184  std::unique_ptr<NullMoistLagrangian> m_moist_model; /*!< moisture model */
185 };
186 
187 #endif
188 #endif
Contains the base class for microphysics.
Contains the Lagrangian moisture model base class.
Base class for microphysics interface.
Definition: ERF_Microphysics.H:14
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 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 void Update_State_Vars_Lev(const int &, amrex::MultiFab &)=0
update ERF state variables from microphysics variables
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 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
static MoistureModelType modelType(const MoistureType a_moisture_type)
query if a specified moisture model is Eulerian or Lagrangian
Definition: ERF_Microphysics.H:86
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:123