ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
EulerianMicrophysics Class Reference

Eulerian microphysics interface. More...

#include <ERF_EulerianMicrophysics.H>

Inheritance diagram for EulerianMicrophysics:
Collaboration diagram for EulerianMicrophysics:

Public Member Functions

 EulerianMicrophysics ()
 Null constructor. More...
 
 ~EulerianMicrophysics ()=default
 default destructor More...
 
 EulerianMicrophysics (const int &nlev, const MoistureType &a_model_type)
 Constructor: create the moisture model. More...
 
void Define (const int &lev, SolverChoice &sc) override
 Define the moisture model. More...
 
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. More...
 
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. More...
 
void Update_Micro_Vars_Lev (const int &lev, amrex::MultiFab &cons_in) override
 update microphysics variables from ERF state variables More...
 
void Update_State_Vars_Lev (const int &lev, amrex::MultiFab &cons_in) override
 update ERF state variables from microphysics variables More...
 
amrex::MultiFab * Get_Qmoist_Ptr (const int &lev, const int &varIdx) override
 get pointer to a moisture variable More...
 
int Get_Qmoist_Size (const int &) override
 get the number of moisture model variables More...
 
int Get_Qstate_Size () override
 get the number of moisture-model-related conserved state variables More...
 
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 More...
 
- Public Member Functions inherited from Microphysics
 Microphysics ()
 Null constructor. More...
 
virtual ~Microphysics ()=default
 default destructor More...
 

Protected Member Functions

template<class NewMoistModel >
void SetModel ()
 Create and set the specified moisture model. More...
 

Private Attributes

amrex::Vector< std::unique_ptr< NullMoist > > m_moist_model
 

Additional Inherited Members

- Static Public Member Functions inherited from Microphysics
static MoistureModelType modelType (const MoistureType a_moisture_type)
 query if a specified moisture model is Eulerian or Lagrangian More...
 

Detailed Description

Eulerian microphysics interface.

Constructor & Destructor Documentation

◆ EulerianMicrophysics() [1/2]

EulerianMicrophysics::EulerianMicrophysics ( )
inline

Null constructor.

19 { }

◆ ~EulerianMicrophysics()

EulerianMicrophysics::~EulerianMicrophysics ( )
default

default destructor

◆ EulerianMicrophysics() [2/2]

EulerianMicrophysics::EulerianMicrophysics ( const int &  nlev,
const MoistureType &  a_model_type 
)
inline

Constructor: create the moisture model.

Parameters
nlevNumber of AMR levels
a_model_typemoisture 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  }
amrex::Vector< std::unique_ptr< NullMoist > > m_moist_model
Definition: ERF_EulerianMicrophysics.H:131
static MoistureModelType modelType(const MoistureType a_moisture_type)
query if a specified moisture model is Eulerian or Lagrangian
Definition: ERF_Microphysics.H:64
Here is the call graph for this function:

Member Function Documentation

◆ Advance()

void EulerianMicrophysics::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 >> &   
)
inlineoverridevirtual

Advance the moisture model for one time step.

Parameters
levAMR level
dt_advanceTime step
solverChoiceSolver choice object

Implements Microphysics.

72  {
73  m_moist_model[lev]->Advance(dt_advance, solverChoice);
74  }

◆ Define()

void EulerianMicrophysics::Define ( const int &  lev,
SolverChoice sc 
)
inlineoverridevirtual

Define the moisture model.

Parameters
levAMR level
scSolver choice object

Implements Microphysics.

47  {
48  m_moist_model[lev]->Define(sc);
49  }

◆ Get_Qmoist_Ptr()

amrex::MultiFab* EulerianMicrophysics::Get_Qmoist_Ptr ( const int &  lev,
const int &  varIdx 
)
inlineoverridevirtual

get pointer to a moisture variable

Parameters
levAMR level
varIdxmoisture variable index

Implements Microphysics.

93  {
94  return m_moist_model[lev]->Qmoist_Ptr(varIdx);
95  }

◆ Get_Qmoist_Restart_Vars()

void EulerianMicrophysics::Get_Qmoist_Restart_Vars ( const int  a_lev,
const SolverChoice a_sc,
std::vector< int > &  a_idx,
std::vector< std::string > &  a_names 
) const
inlineoverridevirtual

get the indices and names of moisture model variables for restart at a given level

Parameters
a_levlevel
a_scSolver choice object
a_idxindices
a_namesnames

Implements Microphysics.

115  {
116  m_moist_model[a_lev]->Qmoist_Restart_Vars( a_sc, a_idx, a_names );
117  }

◆ Get_Qmoist_Size()

int EulerianMicrophysics::Get_Qmoist_Size ( const int &  )
inlineoverridevirtual

get the number of moisture model variables

Implements Microphysics.

99  {
100  return m_moist_model[0]->Qmoist_Size();
101  }

◆ Get_Qstate_Size()

int EulerianMicrophysics::Get_Qstate_Size ( )
inlineoverridevirtual

get the number of moisture-model-related conserved state variables

Implements Microphysics.

105  {
106  return m_moist_model[0]->Qstate_Size();
107  }

◆ Init()

void EulerianMicrophysics::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 
)
inlineoverridevirtual

Initialize the moisture model.

Parameters
levAMR level
cons_inConserved state variables
gridsGrids
geomGeometry
dt_advanceTime step

Implements Microphysics.

59  {
60  m_moist_model[lev]->Init(cons_in, grids, geom, dt_advance,
61  z_phys_nd, detJ_cc);
62  }

◆ SetModel()

template<class NewMoistModel >
void EulerianMicrophysics::SetModel ( )
inlineprotected

Create and set the specified moisture model.

124  {
125  for (int lev(0); lev<m_moist_model.size(); ++lev) {
126  m_moist_model[lev] = std::make_unique<NewMoistModel>();
127  }
128  }

◆ Update_Micro_Vars_Lev()

void EulerianMicrophysics::Update_Micro_Vars_Lev ( const int &  lev,
amrex::MultiFab &  cons_in 
)
inlineoverridevirtual

update microphysics variables from ERF state variables

Parameters
cons_inAMR level Conserved state variables

Implements Microphysics.

79  {
80  m_moist_model[lev]->Update_Micro_Vars(cons_in);
81  }

◆ Update_State_Vars_Lev()

void EulerianMicrophysics::Update_State_Vars_Lev ( const int &  lev,
amrex::MultiFab &  cons_in 
)
inlineoverridevirtual

update ERF state variables from microphysics variables

Parameters
levAMR level
cons_inConserved state variables

Implements Microphysics.

86  {
87  m_moist_model[lev]->Update_State_Vars(cons_in);
88  }

Member Data Documentation

◆ m_moist_model

amrex::Vector<std::unique_ptr<NullMoist> > EulerianMicrophysics::m_moist_model
private

The documentation for this class was generated from the following file: