ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_EulerianMicrophysics.H
Go to the documentation of this file.
1 /*! @file ERF_EulerianMicrophysics.H
2  * \brief Contains the Eulerian microphysics class
3  */
4 
5 #ifndef ERF_EULERIANMICROPHYSICS_H
6 #define ERF_EULERIANMICROPHYSICS_H
7 
8 #include "ERF_NullMoist.H"
9 #include "ERF_SAM.H"
10 #include "ERF_Kessler.H"
11 #include "ERF_Microphysics.H"
12 
13 /*! \brief Eulerian microphysics interface */
15 
16 public:
17 
18  /*! \brief Null constructor */
20 
21  /*! \brief default destructor */
22  ~EulerianMicrophysics () = default;
23 
24  /*! \brief Constructor: create the moisture model */
25  EulerianMicrophysics (const int& nlev, /*!< Number of AMR levels */
26  const MoistureType& a_model_type /*!< moisture model */)
27  {
28  AMREX_ASSERT( Microphysics::modelType(a_model_type) == MoistureModelType::Eulerian );
29  m_moist_model.resize(nlev);
30  if (a_model_type == MoistureType::SAM ||
31  a_model_type == MoistureType::SAM_NoIce ||
32  a_model_type == MoistureType::SAM_NoPrecip_NoIce) {
33  SetModel<SAM>();
34  } else if (a_model_type == MoistureType::Kessler ||
35  a_model_type == MoistureType::Kessler_NoRain) {
36  SetModel<Kessler>();
37  } else if (a_model_type == MoistureType::None) {
38  SetModel<NullMoist>();
39  } else {
40  amrex::Abort("EulerianMicrophysics: Dont know this moisture_type!") ;
41  }
42  }
43 
44  /*! \brief Define the moisture model */
45  void Define (const int& lev, /*!< AMR level */
46  SolverChoice& sc /*!< Solver choice object */) override
47  {
48  m_moist_model[lev]->Define(sc);
49  }
50 
51  /*! \brief Initialize the moisture model */
52  void Init (const int& lev, /*!< AMR level */
53  const amrex::MultiFab& cons_in, /*!< Conserved state variables */
54  const amrex::BoxArray& grids, /*!< Grids */
55  const amrex::Geometry& geom, /*!< Geometry */
56  const amrex::Real& dt_advance, /*!< Time step */
57  std::unique_ptr<amrex::MultiFab>& z_phys_nd, /*< Nodal z heights */
58  std::unique_ptr<amrex::MultiFab>& detJ_cc /*< CC Jacobian determinants */) override
59  {
60  m_moist_model[lev]->Init(cons_in, grids, geom, dt_advance,
61  z_phys_nd, detJ_cc);
62  }
63 
64  /*! \brief Advance the moisture model for one time step */
65  void Advance (const int& lev, /*!< AMR level */
66  const amrex::Real& dt_advance, /*!< Time step */
67  const int&, /*!< iteration number */
68  const amrex::Real&, /*!< current time */
69  const SolverChoice &solverChoice, /*!< Solver choice object */
70  amrex::Vector<amrex::Vector<amrex::MultiFab>>&, /*!< Dycore state variables */
71  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& /*!< terrain */) override
72  {
73  m_moist_model[lev]->Advance(dt_advance, solverChoice);
74  }
75 
76  /*! \brief update microphysics variables from ERF state variables */
77  void Update_Micro_Vars_Lev (const int& lev, /*! AMR level */
78  amrex::MultiFab& cons_in /*!< Conserved state variables */) override
79  {
80  m_moist_model[lev]->Update_Micro_Vars(cons_in);
81  }
82 
83  /*! \brief update ERF state variables from microphysics variables */
84  void Update_State_Vars_Lev (const int& lev, /*!< AMR level */
85  amrex::MultiFab& cons_in /*!< Conserved state variables */) override
86  {
87  m_moist_model[lev]->Update_State_Vars(cons_in);
88  }
89 
90  /*! \brief get pointer to a moisture variable */
91  amrex::MultiFab* Get_Qmoist_Ptr (const int& lev, /*!< AMR level */
92  const int& varIdx /*!< moisture variable index */) override
93  {
94  return m_moist_model[lev]->Qmoist_Ptr(varIdx);
95  }
96 
97  /*! \brief get the number of moisture model variables */
98  int Get_Qmoist_Size (const int& /* lev */) override
99  {
100  return m_moist_model[0]->Qmoist_Size();
101  }
102 
103  /*! \brief get the number of moisture-model-related conserved state variables */
104  int Get_Qstate_Size () override
105  {
106  return m_moist_model[0]->Qstate_Size();
107  }
108 
109  /*! \brief get the indices and names of moisture model variables for restart
110  at a given level */
111  void Get_Qmoist_Restart_Vars ( const int a_lev, /*!< level */
112  const SolverChoice& a_sc, /*!< Solver choice object */
113  std::vector<int>& a_idx, /*!< indices */
114  std::vector<std::string>& a_names /*!< names */ ) const override
115  {
116  m_moist_model[a_lev]->Qmoist_Restart_Vars( a_sc, a_idx, a_names );
117  }
118 
119 protected:
120 
121  /*! \brief Create and set the specified moisture model */
122  template<class NewMoistModel>
123  void SetModel ()
124  {
125  for (int lev(0); lev<m_moist_model.size(); ++lev) {
126  m_moist_model[lev] = std::make_unique<NewMoistModel>();
127  }
128  }
129 
130 private:
131  amrex::Vector<std::unique_ptr<NullMoist>> m_moist_model; /*!< moisture model */
132 };
133 #endif
Contains the base class for microphysics.
Eulerian microphysics interface.
Definition: ERF_EulerianMicrophysics.H:14
void Get_Qmoist_Restart_Vars(const int a_lev, const SolverChoice &a_sc, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
get the indices and names of moisture model variables for restart at a given level
Definition: ERF_EulerianMicrophysics.H:111
void Update_Micro_Vars_Lev(const int &lev, amrex::MultiFab &cons_in) override
update microphysics variables from ERF state variables
Definition: ERF_EulerianMicrophysics.H:77
EulerianMicrophysics(const int &nlev, const MoistureType &a_model_type)
Constructor: create the moisture model.
Definition: ERF_EulerianMicrophysics.H:25
~EulerianMicrophysics()=default
default destructor
amrex::Vector< std::unique_ptr< NullMoist > > m_moist_model
Definition: ERF_EulerianMicrophysics.H:131
int Get_Qstate_Size() override
get the number of moisture-model-related conserved state variables
Definition: ERF_EulerianMicrophysics.H:104
amrex::MultiFab * Get_Qmoist_Ptr(const int &lev, const int &varIdx) override
get pointer to a moisture variable
Definition: ERF_EulerianMicrophysics.H:91
void Define(const int &lev, SolverChoice &sc) override
Define the moisture model.
Definition: ERF_EulerianMicrophysics.H:45
void Advance(const int &lev, const amrex::Real &dt_advance, const int &, const amrex::Real &, const SolverChoice &solverChoice, amrex::Vector< amrex::Vector< amrex::MultiFab >> &, const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &) override
Advance the moisture model for one time step.
Definition: ERF_EulerianMicrophysics.H:65
void Update_State_Vars_Lev(const int &lev, amrex::MultiFab &cons_in) override
update ERF state variables from microphysics variables
Definition: ERF_EulerianMicrophysics.H:84
void SetModel()
Create and set the specified moisture model.
Definition: ERF_EulerianMicrophysics.H:123
void Init(const int &lev, 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
Initialize the moisture model.
Definition: ERF_EulerianMicrophysics.H:52
int Get_Qmoist_Size(const int &) override
get the number of moisture model variables
Definition: ERF_EulerianMicrophysics.H:98
EulerianMicrophysics()
Null constructor.
Definition: ERF_EulerianMicrophysics.H:19
Base class for microphysics interface.
Definition: ERF_Microphysics.H:13
static MoistureModelType modelType(const MoistureType a_moisture_type)
query if a specified moisture model is Eulerian or Lagrangian
Definition: ERF_Microphysics.H:64
Definition: ERF_DataStruct.H:78