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