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