ABLATE Source Documentation  0.12.33
perfectGas.hpp
1 #ifndef ABLATELIBRARY_PERFECTGAS_HPP
2 #define ABLATELIBRARY_PERFECTGAS_HPP
3 #include <map>
4 #include <memory>
5 #include "eos.hpp"
6 #include "parameters/parameters.hpp"
7 #include "utilities/vectorUtilities.hpp"
8 
9 namespace ablate::eos {
10 
11 class PerfectGas : public EOS {
12  private:
13  // constant to be replaced with the number of species for thermodynamic properties
14  constexpr static inline int SPECIES_SIZE = -1;
15 
16  // the perfect gas does not allow species
17  const std::vector<std::string> species;
18  struct Parameters {
19  PetscReal gamma;
20  PetscReal rGas;
21  PetscInt numberSpecies;
22  };
23  Parameters parameters{};
24 
25  struct FunctionContext {
26  PetscInt dim;
27  PetscInt eulerOffset;
28  Parameters parameters;
29  };
30 
39  static PetscErrorCode DensityFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
40  static PetscErrorCode PressureFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
41  static PetscErrorCode TemperatureFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
42  static PetscErrorCode InternalSensibleEnergyFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
43  static PetscErrorCode SensibleEnthalpyFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
44  static PetscErrorCode SpecificHeatConstantVolumeFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
45  static PetscErrorCode SpecificHeatConstantPressureFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
46  static PetscErrorCode SpeedOfSoundFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
47  static PetscErrorCode SpeciesSensibleEnthalpyFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
59  static PetscErrorCode DensityTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
60  static PetscErrorCode PressureTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
61  static PetscErrorCode TemperatureTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
62  static PetscErrorCode InternalSensibleEnergyTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
63  static PetscErrorCode SensibleEnthalpyTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
64  static PetscErrorCode SpecificHeatConstantVolumeTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
65  static PetscErrorCode SpecificHeatConstantPressureTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
66  static PetscErrorCode SpeedOfSoundTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
67  static PetscErrorCode SpeciesSensibleEnthalpyTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
73  using ThermodynamicStaticFunction = PetscErrorCode (*)(const PetscReal conserved[], PetscReal* property, void* ctx);
74  using ThermodynamicTemperatureStaticFunction = PetscErrorCode (*)(const PetscReal conserved[], PetscReal temperature, PetscReal* property, void* ctx);
75  std::map<ThermodynamicProperty, std::tuple<ThermodynamicStaticFunction, ThermodynamicTemperatureStaticFunction, int>> thermodynamicFunctions = {
76  {ThermodynamicProperty::Density, {DensityFunction, DensityTemperatureFunction, 1}},
77  {ThermodynamicProperty::Pressure, {PressureFunction, PressureTemperatureFunction, 1}},
78  {ThermodynamicProperty::Temperature, {TemperatureFunction, TemperatureTemperatureFunction, 1}},
79  {ThermodynamicProperty::InternalSensibleEnergy, {InternalSensibleEnergyFunction, InternalSensibleEnergyTemperatureFunction, 1}},
80  {ThermodynamicProperty::SensibleEnthalpy, {SensibleEnthalpyFunction, SensibleEnthalpyTemperatureFunction, 1}},
81  {ThermodynamicProperty::SpecificHeatConstantVolume, {SpecificHeatConstantVolumeFunction, SpecificHeatConstantVolumeTemperatureFunction, 1}},
82  {ThermodynamicProperty::SpecificHeatConstantPressure, {SpecificHeatConstantPressureFunction, SpecificHeatConstantPressureTemperatureFunction, 1}},
83  {ThermodynamicProperty::SpeedOfSound, {SpeedOfSoundFunction, SpeedOfSoundTemperatureFunction, 1}},
84  {ThermodynamicProperty::SpeciesSensibleEnthalpy, {SpeciesSensibleEnthalpyFunction, SpeciesSensibleEnthalpyTemperatureFunction, SPECIES_SIZE}}};
85 
86  public:
87  explicit PerfectGas(const std::shared_ptr<ablate::parameters::Parameters>&, std::vector<std::string> species = {});
88  void View(std::ostream& stream) const override;
93  [[nodiscard]] PetscReal GetSpecificHeatRatio() const { return parameters.gamma; }
94 
99  [[nodiscard]] PetscReal GetGasConstant() const { return parameters.rGas; }
100 
107  [[nodiscard]] ThermodynamicFunction GetThermodynamicFunction(ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override;
108 
115  [[nodiscard]] ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override;
116 
123  [[nodiscard]] EOSFunction GetFieldFunctionFunction(const std::string& field, ThermodynamicProperty property1, ThermodynamicProperty property2,
124  std::vector<std::string> otherProperties) const override;
125 
130  [[nodiscard]] const std::vector<std::string>& GetSpeciesVariables() const override { return species; }
131 
136  [[nodiscard]] virtual const std::vector<std::string>& GetProgressVariables() const override { return ablate::utilities::VectorUtilities::Empty<std::string>; }
137 };
138 
139 } // namespace ablate::eos
140 #endif // ABLATELIBRARY_PERFECTGAS_HPP
Definition: eos.hpp:60
Definition: perfectGas.hpp:11
ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: perfectGas.cpp:35
virtual const std::vector< std::string > & GetProgressVariables() const override
Definition: perfectGas.hpp:136
PetscReal GetSpecificHeatRatio() const
Definition: perfectGas.hpp:93
const std::vector< std::string > & GetSpeciesVariables() const override
Definition: perfectGas.hpp:130
ThermodynamicFunction GetThermodynamicFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: perfectGas.cpp:23
void View(std::ostream &stream) const override
Definition: perfectGas.cpp:9
PetscReal GetGasConstant() const
Definition: perfectGas.hpp:99
EOSFunction GetFieldFunctionFunction(const std::string &field, ThermodynamicProperty property1, ThermodynamicProperty property2, std::vector< std::string > otherProperties) const override
Definition: perfectGas.cpp:199