ABLATE Source Documentation  0.12.33
eos.hpp
1 #ifndef ABLATELIBRARY_EOS_HPP
2 #define ABLATELIBRARY_EOS_HPP
3 #include <petsc.h>
4 #include <functional>
5 #include <iostream>
6 #include <memory>
7 #include <string>
8 #include <utility>
9 #include <vector>
10 #include "domain/field.hpp"
11 #include "domain/fieldDescriptor.hpp"
12 
13 namespace ablate::eos {
14 
15 enum class ThermodynamicProperty {
16  Pressure,
17  Temperature,
18  InternalSensibleEnergy,
19  SensibleEnthalpy,
20  SpecificHeatConstantVolume,
21  SpecificHeatConstantPressure,
22  SpeedOfSound,
23  SpeciesSensibleEnthalpy,
24  Density
25 };
26 
32  PetscErrorCode (*function)(const PetscReal conserved[], PetscReal* property, void* ctx) = nullptr;
34  std::shared_ptr<void> context = nullptr;
36  PetscInt propertySize = 1;
37 };
38 
44  PetscErrorCode (*function)(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx) = nullptr;
46  std::shared_ptr<void> context = nullptr;
48  PetscInt propertySize = 1;
49 };
50 
54 using EOSFunction = std::function<void(PetscReal property1, PetscReal property2, PetscInt dim, const PetscReal velocity[], const PetscReal other[], PetscReal conserved[])>;
55 
60 class EOS {
61  protected:
62  const std::string type;
63 
64  public:
65  explicit EOS(std::string typeIn) : type(std::move(typeIn)){};
66  virtual ~EOS() = default;
67 
68  // Some known extra properties
69  // the conserved (density*yi) solution field for species mass fractions
70  inline const static std::string YI = "Yi";
71  // progress fields are used by the eos/chemistry model to transport required non species
72  inline const static std::string PROGRESS = "Progress";
73 
78  virtual void View(std::ostream& stream) const = 0;
79 
86  [[nodiscard]] virtual ThermodynamicFunction GetThermodynamicFunction(ThermodynamicProperty property, const std::vector<domain::Field>& fields) const = 0;
87 
94  [[nodiscard]] virtual ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(ThermodynamicProperty property, const std::vector<domain::Field>& fields) const = 0;
95 
102  [[nodiscard]] virtual EOSFunction GetFieldFunctionFunction(const std::string& field, ThermodynamicProperty property1, ThermodynamicProperty property2,
103  std::vector<std::string> otherProperties) const = 0;
104 
110  [[nodiscard]] virtual const std::vector<std::string>& GetSpeciesVariables() const = 0;
111 
116  [[nodiscard]] virtual const std::vector<std::string>& GetProgressVariables() const = 0;
117 
122  [[nodiscard]] virtual std::vector<std::string> GetFieldTags() const { return std::vector<std::string>{}; }
123 
128  [[nodiscard]] virtual std::vector<std::shared_ptr<domain::FieldDescriptor>> GetAdditionalFields() const { return {}; }
129 
135  [[nodiscard]] virtual const std::vector<std::string>& GetFieldFunctionProperties() const { return GetSpeciesVariables(); }
136 
143  friend std::ostream& operator<<(std::ostream& out, const EOS& eos) {
144  eos.View(out);
145  return out;
146  }
147 };
148 
154 constexpr std::string_view to_string(const ThermodynamicProperty& prop) {
155  switch (prop) {
156  case ThermodynamicProperty::Pressure:
157  return "pressure";
158  case ThermodynamicProperty::Temperature:
159  return "temperature";
160  case ThermodynamicProperty::InternalSensibleEnergy:
161  return "internalSensibleEnergy";
162  case ThermodynamicProperty::SensibleEnthalpy:
163  return "sensibleEnthalpy";
164  case ThermodynamicProperty::SpecificHeatConstantVolume:
165  return "specificHeatConstantVolume";
166  case ThermodynamicProperty::SpecificHeatConstantPressure:
167  return "specificHeatConstantPressure";
168  case ThermodynamicProperty::SpeedOfSound:
169  return "speedOfSound";
170  case ThermodynamicProperty::SpeciesSensibleEnthalpy:
171  return "speciesSensibleEnthalpy";
172  case ThermodynamicProperty::Density:
173  return "density";
174  }
175  return "";
176 }
177 
183 constexpr ThermodynamicProperty from_string(const std::string_view& prop) {
184  if (prop == "pressure") return ThermodynamicProperty::Pressure;
185  if (prop == "temperature") return ThermodynamicProperty::Temperature;
186  if (prop == "internalSensibleEnergy") return ThermodynamicProperty::InternalSensibleEnergy;
187  if (prop == "sensibleEnthalpy") return ThermodynamicProperty::SensibleEnthalpy;
188  if (prop == "specificHeatConstantVolume") return ThermodynamicProperty::SpecificHeatConstantVolume;
189  if (prop == "specificHeatConstantPressure") return ThermodynamicProperty::SpecificHeatConstantPressure;
190  if (prop == "speedOfSound") return ThermodynamicProperty::SpeedOfSound;
191  if (prop == "speciesSensibleEnthalpy") return ThermodynamicProperty::SpeciesSensibleEnthalpy;
192  if (prop == "density") return ThermodynamicProperty::Density;
193  throw std::invalid_argument("No known ThermodynamicProperty for " + std::string(prop));
194 }
195 
202 inline std::ostream& operator<<(std::ostream& out, const ThermodynamicProperty& prop) {
203  auto string = to_string(prop);
204  out << string;
205  return out;
206 }
207 
214 inline std::istream& operator>>(std::istream& in, ThermodynamicProperty& prop) {
215  std::string propString;
216  in >> propString;
217  prop = from_string(propString);
218  return in;
219 }
220 
221 } // namespace ablate::eos
222 
223 #endif // ABLATELIBRARY_EOS_HPP
Definition: eos.hpp:60
virtual ThermodynamicFunction GetThermodynamicFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const =0
virtual EOSFunction GetFieldFunctionFunction(const std::string &field, ThermodynamicProperty property1, ThermodynamicProperty property2, std::vector< std::string > otherProperties) const =0
virtual std::vector< std::string > GetFieldTags() const
Definition: eos.hpp:122
virtual const std::vector< std::string > & GetFieldFunctionProperties() const
Definition: eos.hpp:135
virtual std::vector< std::shared_ptr< domain::FieldDescriptor > > GetAdditionalFields() const
Definition: eos.hpp:128
virtual const std::vector< std::string > & GetProgressVariables() const =0
friend std::ostream & operator<<(std::ostream &out, const EOS &eos)
Definition: eos.hpp:143
virtual void View(std::ostream &stream) const =0
virtual ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const =0
virtual const std::vector< std::string > & GetSpeciesVariables() const =0
PetscInt propertySize
the property size being set
Definition: eos.hpp:36
std::shared_ptr< void > context
optional context to pass into the function
Definition: eos.hpp:34
std::shared_ptr< void > context
optional context to pass into the function
Definition: eos.hpp:46
PetscInt propertySize
the property size being set
Definition: eos.hpp:48