ABLATE Source Documentation  0.12.33
twoPhase.hpp
1 #ifndef ABLATELIBRARY_TWOPHASE_HPP
2 #define ABLATELIBRARY_TWOPHASE_HPP
3 
4 #include <memory>
5 #include "eos.hpp"
6 #include "eos/perfectGas.hpp"
7 #include "eos/stiffenedGas.hpp"
8 #include "parameters/parameters.hpp"
9 namespace ablate::eos {
10 class TwoPhase : public EOS { // , public std::enabled_shared_from_this<TwoPhase>
11  public:
12  inline const static std::string VF = "volumeFraction";
13 
14  private:
15  const std::shared_ptr<eos::EOS> eos1;
16  const std::shared_ptr<eos::EOS> eos2;
17  // this mixed eos does not allow species, get species from eos1 eos2
18  std::vector<std::string> species;
19  std::vector<std::string> otherPropertiesList = {"VF"}; // otherProperties must include volumeFraction for Field Function initialization
20  struct Parameters {
21  PetscReal gamma1;
22  PetscReal rGas1;
23  PetscReal Cp1;
24  PetscReal p01;
25  PetscReal gamma2;
26  PetscReal rGas2;
27  PetscReal Cp2;
28  PetscReal p02;
29  PetscInt numberSpecies1;
30  std::vector<std::string> species1;
31  PetscInt numberSpecies2;
32  std::vector<std::string> species2;
33  };
34  Parameters parameters;
35  struct FunctionContext {
36  PetscInt dim;
37  PetscInt eulerOffset;
38  PetscInt densityVFOffset;
39  PetscInt volumeFractionOffset;
40  Parameters parameters;
41  };
42 
43  public:
44  struct DecodeIn {
45  PetscReal alpha;
46  PetscReal alphaRho1;
47  PetscReal rho;
48  PetscReal e;
49  Parameters parameters;
50  };
51  struct DecodeOut {
52  PetscReal rho1;
53  PetscReal rho2;
54  PetscReal e1;
55  PetscReal e2;
56  PetscReal p;
57  PetscReal T;
58  };
59 
60  private:
61  // functions for all cases
62  static PetscErrorCode DensityFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
63  static PetscErrorCode InternalSensibleEnergyFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
64  static PetscErrorCode SpeciesSensibleEnthalpyFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
65 
66  static PetscErrorCode DensityTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
67  static PetscErrorCode InternalSensibleEnergyTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
68  static PetscErrorCode SpeciesSensibleEnthalpyTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
69 
70  // functions for GasGas case
71  static PetscErrorCode PressureFunctionGasGas(const PetscReal conserved[], PetscReal* property, void* ctx);
72  static PetscErrorCode TemperatureFunctionGasGas(const PetscReal conserved[], PetscReal* property, void* ctx);
73  static PetscErrorCode SensibleEnthalpyFunctionGasGas(const PetscReal conserved[], PetscReal* property, void* ctx);
74  static PetscErrorCode SpecificHeatConstantVolumeFunctionGasGas(const PetscReal conserved[], PetscReal* property, void* ctx);
75  static PetscErrorCode SpecificHeatConstantPressureFunctionGasGas(const PetscReal conserved[], PetscReal* property, void* ctx);
76  static PetscErrorCode SpeedOfSoundFunctionGasGas(const PetscReal conserved[], PetscReal* property, void* ctx);
77 
78  static PetscErrorCode PressureTemperatureFunctionGasGas(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
79  static PetscErrorCode TemperatureTemperatureFunctionGasGas(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
80  static PetscErrorCode SensibleEnthalpyTemperatureFunctionGasGas(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
81  static PetscErrorCode SpecificHeatConstantVolumeTemperatureFunctionGasGas(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
82  static PetscErrorCode SpecificHeatConstantPressureTemperatureFunctionGasGas(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
83  static PetscErrorCode SpeedOfSoundTemperatureFunctionGasGas(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
84 
85  // functions for GasLiquid
86  static PetscErrorCode PressureFunctionGasLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
87  static PetscErrorCode TemperatureFunctionGasLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
88  static PetscErrorCode SensibleEnthalpyFunctionGasLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
89  static PetscErrorCode SpecificHeatConstantVolumeFunctionGasLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
90  static PetscErrorCode SpecificHeatConstantPressureFunctionGasLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
91  static PetscErrorCode SpeedOfSoundFunctionGasLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
92 
93  static PetscErrorCode PressureTemperatureFunctionGasLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
94  static PetscErrorCode TemperatureTemperatureFunctionGasLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
95  static PetscErrorCode SensibleEnthalpyTemperatureFunctionGasLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
96  static PetscErrorCode SpecificHeatConstantVolumeTemperatureFunctionGasLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
97  static PetscErrorCode SpecificHeatConstantPressureTemperatureFunctionGasLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
98  static PetscErrorCode SpeedOfSoundTemperatureFunctionGasLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
99 
100  // functions for LiquidLiquid
101  static PetscErrorCode PressureFunctionLiquidLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
102  static PetscErrorCode TemperatureFunctionLiquidLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
103  static PetscErrorCode SensibleEnthalpyFunctionLiquidLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
104  static PetscErrorCode SpecificHeatConstantVolumeFunctionLiquidLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
105  static PetscErrorCode SpecificHeatConstantPressureFunctionLiquidLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
106  static PetscErrorCode SpeedOfSoundFunctionLiquidLiquid(const PetscReal conserved[], PetscReal* property, void* ctx);
107 
108  static PetscErrorCode PressureTemperatureFunctionLiquidLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
109  static PetscErrorCode TemperatureTemperatureFunctionLiquidLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
110  static PetscErrorCode SensibleEnthalpyTemperatureFunctionLiquidLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
111  static PetscErrorCode SpecificHeatConstantVolumeTemperatureFunctionLiquidLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
112  static PetscErrorCode SpecificHeatConstantPressureTemperatureFunctionLiquidLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
113  static PetscErrorCode SpeedOfSoundTemperatureFunctionLiquidLiquid(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
114 
115  using ThermodynamicStaticFunction = PetscErrorCode (*)(const PetscReal conserved[], PetscReal* property, void* ctx);
116  using ThermodynamicTemperatureStaticFunction = PetscErrorCode (*)(const PetscReal conserved[], PetscReal temperature, PetscReal* property, void* ctx);
117  // map for GasGas case
118  std::map<ThermodynamicProperty, std::pair<ThermodynamicStaticFunction, ThermodynamicTemperatureStaticFunction>> thermodynamicFunctionsGasGas = {
119  {ThermodynamicProperty::Density, {DensityFunction, DensityTemperatureFunction}},
120  {ThermodynamicProperty::Pressure, {PressureFunctionGasGas, PressureTemperatureFunctionGasGas}},
121  {ThermodynamicProperty::Temperature, {TemperatureFunctionGasGas, TemperatureTemperatureFunctionGasGas}},
122  {ThermodynamicProperty::InternalSensibleEnergy, {InternalSensibleEnergyFunction, InternalSensibleEnergyTemperatureFunction}},
123  {ThermodynamicProperty::SensibleEnthalpy, {SensibleEnthalpyFunctionGasGas, SensibleEnthalpyTemperatureFunctionGasGas}},
124  {ThermodynamicProperty::SpecificHeatConstantVolume, {SpecificHeatConstantVolumeFunctionGasGas, SpecificHeatConstantVolumeTemperatureFunctionGasGas}},
125  {ThermodynamicProperty::SpecificHeatConstantPressure, {SpecificHeatConstantPressureFunctionGasGas, SpecificHeatConstantPressureTemperatureFunctionGasGas}},
126  {ThermodynamicProperty::SpeedOfSound, {SpeedOfSoundFunctionGasGas, SpeedOfSoundTemperatureFunctionGasGas}},
127  {ThermodynamicProperty::SpeciesSensibleEnthalpy, {SpeciesSensibleEnthalpyFunction, SpeciesSensibleEnthalpyTemperatureFunction}}};
128  // map for GasLiquid case
129  std::map<ThermodynamicProperty, std::pair<ThermodynamicStaticFunction, ThermodynamicTemperatureStaticFunction>> thermodynamicFunctionsGasLiquid = {
130  {ThermodynamicProperty::Density, {DensityFunction, DensityTemperatureFunction}},
131  {ThermodynamicProperty::Pressure, {PressureFunctionGasLiquid, PressureTemperatureFunctionGasLiquid}},
132  {ThermodynamicProperty::Temperature, {TemperatureFunctionGasLiquid, TemperatureTemperatureFunctionGasLiquid}},
133  {ThermodynamicProperty::InternalSensibleEnergy, {InternalSensibleEnergyFunction, InternalSensibleEnergyTemperatureFunction}},
134  {ThermodynamicProperty::SensibleEnthalpy, {SensibleEnthalpyFunctionGasLiquid, SensibleEnthalpyTemperatureFunctionGasLiquid}},
135  {ThermodynamicProperty::SpecificHeatConstantVolume, {SpecificHeatConstantVolumeFunctionGasLiquid, SpecificHeatConstantVolumeTemperatureFunctionGasLiquid}},
136  {ThermodynamicProperty::SpecificHeatConstantPressure, {SpecificHeatConstantPressureFunctionGasLiquid, SpecificHeatConstantPressureTemperatureFunctionGasLiquid}},
137  {ThermodynamicProperty::SpeedOfSound, {SpeedOfSoundFunctionGasLiquid, SpeedOfSoundTemperatureFunctionGasLiquid}},
138  {ThermodynamicProperty::SpeciesSensibleEnthalpy, {SpeciesSensibleEnthalpyFunction, SpeciesSensibleEnthalpyTemperatureFunction}}};
139  // map for LiquidLiquid case
140  std::map<ThermodynamicProperty, std::pair<ThermodynamicStaticFunction, ThermodynamicTemperatureStaticFunction>> thermodynamicFunctionsLiquidLiquid = {
141  {ThermodynamicProperty::Density, {DensityFunction, DensityTemperatureFunction}},
142  {ThermodynamicProperty::Pressure, {PressureFunctionLiquidLiquid, PressureTemperatureFunctionLiquidLiquid}},
143  {ThermodynamicProperty::Temperature, {TemperatureFunctionLiquidLiquid, TemperatureTemperatureFunctionLiquidLiquid}},
144  {ThermodynamicProperty::InternalSensibleEnergy, {InternalSensibleEnergyFunction, InternalSensibleEnergyTemperatureFunction}},
145  {ThermodynamicProperty::SensibleEnthalpy, {SensibleEnthalpyFunctionLiquidLiquid, SensibleEnthalpyTemperatureFunctionLiquidLiquid}},
146  {ThermodynamicProperty::SpecificHeatConstantVolume, {SpecificHeatConstantVolumeFunctionLiquidLiquid, SpecificHeatConstantVolumeTemperatureFunctionLiquidLiquid}},
147  {ThermodynamicProperty::SpecificHeatConstantPressure, {SpecificHeatConstantPressureFunctionLiquidLiquid, SpecificHeatConstantPressureTemperatureFunctionLiquidLiquid}},
148  {ThermodynamicProperty::SpeedOfSound, {SpeedOfSoundFunctionLiquidLiquid, SpeedOfSoundTemperatureFunctionLiquidLiquid}},
149  {ThermodynamicProperty::SpeciesSensibleEnthalpy, {SpeciesSensibleEnthalpyFunction, SpeciesSensibleEnthalpyTemperatureFunction}}};
150 
154  const std::set<ThermodynamicProperty> speciesSizedProperties = {ThermodynamicProperty::SpeciesSensibleEnthalpy};
155 
156  public:
157  explicit TwoPhase(std::shared_ptr<eos::EOS> eos1, std::shared_ptr<eos::EOS> eos2);
158  void View(std::ostream& stream) const override;
159 
160  const std::shared_ptr<ablate::eos::EOS> GetEOSGas() const { return eos1; }
161  const std::shared_ptr<ablate::eos::EOS> GetEOSLiquid() const { return eos2; }
162 
163  ThermodynamicFunction GetThermodynamicFunction(ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override;
164 
165  ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override;
166 
167  EOSFunction GetFieldFunctionFunction(const std::string& field, ThermodynamicProperty property1, ThermodynamicProperty property2, std::vector<std::string> otherProperties) const override;
168  const std::vector<std::string>& GetFieldFunctionProperties() const override { return otherPropertiesList; } // list of other properties i.e. VF;
169 
170  const std::vector<std::string>& GetSpeciesVariables() const override { return species; } // lists species of eos1 first, then eos2, no distinction for which fluid the species exists in
171  [[nodiscard]] virtual const std::vector<std::string>& GetProgressVariables() const override { return ablate::utilities::VectorUtilities::Empty<std::string>; }
172 };
173 } // namespace ablate::eos
174 
175 #endif // ABLATELIBRARY_TWOPHASE_HPP
Definition: eos.hpp:60
Definition: twoPhase.hpp:10
const std::vector< std::string > & GetFieldFunctionProperties() const override
Definition: twoPhase.hpp:168
void View(std::ostream &stream) const override
Definition: twoPhase.cpp:212
ThermodynamicFunction GetThermodynamicFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: twoPhase.cpp:220
const std::vector< std::string > & GetSpeciesVariables() const override
Definition: twoPhase.hpp:170
ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: twoPhase.cpp:268
EOSFunction GetFieldFunctionFunction(const std::string &field, ThermodynamicProperty property1, ThermodynamicProperty property2, std::vector< std::string > otherProperties) const override
Definition: twoPhase.cpp:308
virtual const std::vector< std::string > & GetProgressVariables() const override
Definition: twoPhase.hpp:171
Definition: twoPhase.hpp:44
Definition: twoPhase.hpp:51