Canoe
Comprehensive Atmosphere N' Ocean Engine
microphysics.hpp
Go to the documentation of this file.
1 #ifndef SRC_MICROPHYSICS_MICROPHYSICS_HPP_
2 #define SRC_MICROPHYSICS_MICROPHYSICS_HPP_
3 
4 // C/C++
5 #include <map>
6 #include <memory>
7 #include <vector>
8 
9 // athena
10 #include <athena/athena.hpp>
11 
12 // canoe
13 #include <air_parcel.hpp>
14 #include <schedulers.hpp>
15 #include <virtual_groups.hpp>
16 
17 class MeshBlock;
18 class ParameterInput;
20 class Hydro;
21 
23 class Microphysics {
24  public:
28 
30  static const std::string input_key;
31 
34 
37 
40 
43 
44  public:
45  Microphysics(MeshBlock *pmb, ParameterInput *pin);
46  ~Microphysics();
47 
48  public:
49  size_t GetNumSystems() const { return systems_.size(); }
50  std::shared_ptr<MicrophysicalSchemeBase> GetSystem(int i) const {
51  return systems_[i];
52  }
53  // void AddFrictionalHeating(std::vector<AirParcel> &air_column) const;
54 
60  void EvolveSystems(AirColumn &ac, Real time, Real dt);
61 
62  public:
63  void SetVsedFromConserved(Hydro const *phydro);
64 
65  protected:
68 
70  std::vector<std::shared_ptr<MicrophysicalSchemeBase>> systems_;
71 
72  private:
74  MeshBlock const *pmy_block_;
75 };
76 
77 using MicrophysicsPtr = std::shared_ptr<Microphysics>;
78 
79 namespace AllTasks {
80 
81 // hydro tasks should be move into hydro in the future
82 bool hydro_implicit_correction(MeshBlock *pmb, IntegrationStage stage);
83 bool hydro_calculate_flux(MeshBlock *pmb, IntegrationStage stage);
84 
85 bool microphysics_set_sedimentaton_velocity(MeshBlock *pmb,
86  IntegrationStage stage);
87 bool microphysics_set_mass_flux(MeshBlock *pmb, IntegrationStage stage);
88 bool microphysics_evolve_system(MeshBlock *pmb, IntegrationStage stage);
89 
90 // scalar tasks should be move into scalars in the future
91 bool scalar_calculate_flux(MeshBlock *pmb, IntegrationStage stage);
92 
93 } // namespace AllTasks
94 
95 #endif // SRC_MICROPHYSICS_MICROPHYSICS_HPP_
std::vector< AirParcel > AirColumn
Definition: air_parcel.hpp:121
virtual base class for all microphysical schemes
root-level management class for microphysics
AthenaArray< Real > w
primitive variables: mass fraction [kg/kg]
void SetVsedFromConserved(Hydro const *phydro)
inbound functions
Microphysics(MeshBlock *pmb, ParameterInput *pin)
constructor and destructor
std::vector< std::shared_ptr< MicrophysicalSchemeBase > > systems_
pointers of microphysical systems
AthenaArray< Real > u
conserved variables: mass concentration [kg/m^3]
static const std::string input_key
microphysics input key in the input file [microphysics_config]
MeshBlock const * pmy_block_
meshblock pointer
AthenaArray< Real > vsed_[3]
sedimentation velocity at cell center [m/s]
size_t GetNumSystems() const
functions
void EvolveSystems(AirColumn &ac, Real time, Real dt)
Evolve all microphysical systems.
AthenaArray< Real > mass_flux[3]
mass flux of the dry fluid [kg/m^2/s]
AthenaArray< Real > vsedf[3]
sedimentation velocity at cell interface [m/s]
std::shared_ptr< MicrophysicalSchemeBase > GetSystem(int i) const
std::shared_ptr< Microphysics > MicrophysicsPtr
bool microphysics_set_mass_flux(MeshBlock *pmb, IntegrationStage stage)
bool hydro_calculate_flux(MeshBlock *pmb, IntegrationStage stage)
bool microphysics_evolve_system(MeshBlock *pmb, IntegrationStage stage)
bool scalar_calculate_flux(MeshBlock *pmb, IntegrationStage stage)
bool hydro_implicit_correction(MeshBlock *pmb, IntegrationStage stage)
bool microphysics_set_sedimentaton_velocity(MeshBlock *pmb, IntegrationStage stage)
std::pair< int, int > IntegrationStage
Definition: schedulers.hpp:12