5 #ifndef ERF_EULERIANMICROPHYSICS_H
6 #define ERF_EULERIANMICROPHYSICS_H
27 const MoistureType& a_model_type )
31 if (a_model_type == MoistureType::SAM ||
32 a_model_type == MoistureType::SAM_NoIce ||
33 a_model_type == MoistureType::SAM_NoPrecip_NoIce) {
35 }
else if (a_model_type == MoistureType::Kessler ||
36 a_model_type == MoistureType::Kessler_NoRain) {
38 }
else if (a_model_type == MoistureType::SatAdj) {
40 }
else if (a_model_type == MoistureType::None) {
41 SetModel<NullMoist>();
43 amrex::Abort(
"EulerianMicrophysics: Dont know this moisture_type!") ;
56 const amrex::MultiFab& cons_in,
57 const amrex::BoxArray& grids,
58 const amrex::Geometry& geom,
59 const amrex::Real& dt_advance,
60 std::unique_ptr<amrex::MultiFab>& z_phys_nd,
61 std::unique_ptr<amrex::MultiFab>& detJ_cc )
override
69 const amrex::Real& dt_advance,
73 amrex::Vector<amrex::Vector<amrex::MultiFab>>&,
74 const amrex::Vector<std::unique_ptr<amrex::MultiFab>>&,
75 const amrex::GpuArray<ERF_BC, AMREX_SPACEDIM*2>& )
override
82 amrex::MultiFab& cons_in )
override
89 amrex::MultiFab& cons_in )
override
96 const int& varIdx )
override
117 std::vector<int>& a_idx,
118 std::vector<std::string>& a_names )
const override
120 m_moist_model[a_lev]->Qmoist_Restart_Vars(a_sc, a_idx, a_names);
126 template<
class NewMoistModel>
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