ABLATE Source Documentation  0.12.34
tChem.hpp
1 #ifndef ABLATELIBRARY_TCHEM_HPP
2 #define ABLATELIBRARY_TCHEM_HPP
3 
4 #include <filesystem>
5 #include <map>
6 #include <memory>
7 #include "TChem_KineticModelData.hpp"
8 #include "chemistryModel.hpp"
9 #include "eos.hpp"
10 #include "eos/tChem/pressure.hpp"
11 #include "eos/tChem/sensibleEnthalpy.hpp"
12 #include "eos/tChem/sensibleInternalEnergy.hpp"
13 #include "eos/tChem/sourceCalculator.hpp"
14 #include "eos/tChem/speedOfSound.hpp"
15 #include "eos/tChem/temperature.hpp"
16 #include "monitors/logs/log.hpp"
17 #include "parameters/parameters.hpp"
18 #include "tChemBase.hpp"
19 #include "utilities/intErrorChecker.hpp"
20 #include "utilities/vectorUtilities.hpp"
21 
22 namespace ablate::eos {
23 
24 namespace tChemLib = TChem;
25 
26 class TChem : public TChemBase, public std::enable_shared_from_this<ablate::eos::TChem> {
27  public:
33  PetscErrorCode (*function)(const PetscReal conserved[], const PetscReal yi[], PetscReal* property, void* ctx) = nullptr;
35  std::shared_ptr<void> context = nullptr;
37  PetscInt propertySize = 1;
38  };
39 
40  public:
46  explicit TChem(std::filesystem::path mechanismFile, std::shared_ptr<ablate::monitors::logs::Log> = {}, const std::shared_ptr<ablate::parameters::Parameters>& options = {});
47 
54  [[nodiscard]] ThermodynamicFunction GetThermodynamicFunction(ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override;
55 
62  [[nodiscard]] ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override;
63 
70  [[nodiscard]] ThermodynamicMassFractionFunction GetThermodynamicMassFractionFunction(ThermodynamicProperty property, const std::vector<domain::Field>& fields) const;
71 
78  [[nodiscard]] ThermodynamicTemperatureMassFractionFunction GetThermodynamicTemperatureMassFractionFunction(ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override;
79 
86  [[nodiscard]] EOSFunction GetFieldFunctionFunction(const std::string& field, ThermodynamicProperty property1, ThermodynamicProperty property2,
87  std::vector<std::string> otherProperties) const override;
88 
93  [[nodiscard]] std::map<std::string, double> GetElementInformation() const override;
94 
99  [[nodiscard]] std::map<std::string, std::map<std::string, int>> GetSpeciesElementalInformation() const override;
100 
105  [[nodiscard]] std::map<std::string, double> GetSpeciesMolecularMass() const override;
106 
113  std::shared_ptr<SourceCalculator> CreateSourceCalculator(const std::vector<domain::Field>& fields, const ablate::domain::Range& cellRange) override;
114 
115  protected:
124  [[nodiscard]] std::shared_ptr<FunctionContext> BuildFunctionContext(ablate::eos::ThermodynamicProperty property, const std::vector<domain::Field>& fields, bool checkDensityYi = true) const;
125 
134  static PetscErrorCode DensityFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
135  static PetscErrorCode TemperatureFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
136  static PetscErrorCode PressureFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
137  static PetscErrorCode InternalSensibleEnergyFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
138  static PetscErrorCode SensibleEnthalpyFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
139  static PetscErrorCode SpecificHeatConstantVolumeFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
140  static PetscErrorCode SpecificHeatConstantPressureFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
141  static PetscErrorCode SpeedOfSoundFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
142  static PetscErrorCode SpeciesSensibleEnthalpyFunction(const PetscReal conserved[], PetscReal* property, void* ctx);
154  static PetscErrorCode DensityTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
155  static PetscErrorCode TemperatureTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
156  static PetscErrorCode PressureTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
157  static PetscErrorCode InternalSensibleEnergyTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
158  static PetscErrorCode SensibleEnthalpyTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
159  static PetscErrorCode SpecificHeatConstantVolumeTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
160  static PetscErrorCode SpecificHeatConstantPressureTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
161  static PetscErrorCode SpeedOfSoundTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
162  static PetscErrorCode SpeciesSensibleEnthalpyTemperatureFunction(const PetscReal conserved[], PetscReal T, PetscReal* property, void* ctx);
173  static PetscErrorCode DensityMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal* property, void* ctx);
174  static PetscErrorCode TemperatureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal* property, void* ctx);
175  static PetscErrorCode PressureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal* property, void* ctx);
176  static PetscErrorCode InternalSensibleEnergyMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal* property, void* ctx);
177  static PetscErrorCode SensibleEnthalpyMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal* property, void* ctx);
178  static PetscErrorCode SpecificHeatConstantVolumeMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal* property, void* ctx);
179  static PetscErrorCode SpecificHeatConstantPressureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal* property, void* ctx);
180  static PetscErrorCode SpeedOfSoundMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal* property, void* ctx);
181  static PetscErrorCode SpeciesSensibleEnthalpyMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal* property, void* ctx);
193  static PetscErrorCode DensityTemperatureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal T, PetscReal* property, void* ctx);
194  static PetscErrorCode TemperatureTemperatureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal T, PetscReal* property, void* ctx);
195  static PetscErrorCode PressureTemperatureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal T, PetscReal* property, void* ctx);
196  static PetscErrorCode InternalSensibleEnergyTemperatureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal T, PetscReal* property, void* ctx);
197  static PetscErrorCode SensibleEnthalpyTemperatureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal T, PetscReal* property, void* ctx);
198  static PetscErrorCode SpecificHeatConstantVolumeTemperatureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal T, PetscReal* property, void* ctx);
199  static PetscErrorCode SpecificHeatConstantPressureTemperatureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal T, PetscReal* property, void* ctx);
200  static PetscErrorCode SpeedOfSoundTemperatureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal T, PetscReal* property, void* ctx);
201  static PetscErrorCode SpeciesSensibleEnthalpyTemperatureMassFractionFunction(const PetscReal conserved[], const PetscReal yi[], PetscReal T, PetscReal* property, void* ctx);
211  using ThermodynamicStaticFunction = PetscErrorCode (*)(const PetscReal conserved[], PetscReal* property, void* ctx);
212  using ThermodynamicTemperatureStaticFunction = PetscErrorCode (*)(const PetscReal conserved[], PetscReal temperature, PetscReal* property, void* ctx);
213  std::map<ThermodynamicProperty, std::tuple<ThermodynamicStaticFunction, ThermodynamicTemperatureStaticFunction, std::function<ordinal_type(ordinal_type)>>> thermodynamicFunctions = {
214  {ThermodynamicProperty::Density, {DensityFunction, DensityTemperatureFunction, [](auto) { return 0; }}},
215  {ThermodynamicProperty::Pressure,
216  {PressureFunction, PressureTemperatureFunction, ablate::eos::tChem::Temperature::getWorkSpaceSize}} ,
217  {ThermodynamicProperty::Temperature, {TemperatureFunction, TemperatureTemperatureFunction, ablate::eos::tChem::Temperature::getWorkSpaceSize}},
218  {ThermodynamicProperty::InternalSensibleEnergy, {InternalSensibleEnergyFunction, InternalSensibleEnergyTemperatureFunction, ablate::eos::tChem::SensibleInternalEnergy::getWorkSpaceSize}},
219  {ThermodynamicProperty::SensibleEnthalpy, {SensibleEnthalpyFunction, SensibleEnthalpyTemperatureFunction, ablate::eos::tChem::SensibleEnthalpy::getWorkSpaceSize}},
220  {ThermodynamicProperty::SpecificHeatConstantVolume, {SpecificHeatConstantVolumeFunction, SpecificHeatConstantVolumeTemperatureFunction, [](auto nSpec) { return nSpec; }}},
221  {ThermodynamicProperty::SpecificHeatConstantPressure,
222  {SpecificHeatConstantPressureFunction,
223  SpecificHeatConstantPressureTemperatureFunction,
224  ablate::eos::tChem::Temperature::getWorkSpaceSize}} ,
225  {ThermodynamicProperty::SpeedOfSound, {SpeedOfSoundFunction, SpeedOfSoundTemperatureFunction, ablate::eos::tChem::SpeedOfSound::getWorkSpaceSize}},
226  {ThermodynamicProperty::SpeciesSensibleEnthalpy,
227  {SpeciesSensibleEnthalpyFunction,
228  SpeciesSensibleEnthalpyTemperatureFunction,
229  ablate::eos::tChem::Temperature::getWorkSpaceSize}}
230  };
231 
235  using ThermodynamicStaticMassFractionFunction = PetscErrorCode (*)(const PetscReal conserved[], const PetscReal yi[], PetscReal* property, void* ctx);
236  using ThermodynamicTemperatureStaticMassFractionFunction = PetscErrorCode (*)(const PetscReal conserved[], const PetscReal yi[], PetscReal temperature, PetscReal* property, void* ctx);
237  std::map<ThermodynamicProperty, std::tuple<ThermodynamicStaticMassFractionFunction, ThermodynamicTemperatureStaticMassFractionFunction, std::function<ordinal_type(ordinal_type)>>>
238  thermodynamicMassFractionFunctions = {
239  {ThermodynamicProperty::Density, {DensityMassFractionFunction, DensityTemperatureMassFractionFunction, [](auto) { return 0; }}},
240  {ThermodynamicProperty::Pressure,
241  {PressureMassFractionFunction,
242  PressureTemperatureMassFractionFunction,
243  ablate::eos::tChem::Temperature::getWorkSpaceSize}} ,
244  {ThermodynamicProperty::Temperature, {TemperatureMassFractionFunction, TemperatureTemperatureMassFractionFunction, ablate::eos::tChem::Temperature::getWorkSpaceSize}},
245  {ThermodynamicProperty::InternalSensibleEnergy,
246  {InternalSensibleEnergyMassFractionFunction, InternalSensibleEnergyTemperatureMassFractionFunction, ablate::eos::tChem::SensibleInternalEnergy::getWorkSpaceSize}},
247  {ThermodynamicProperty::SensibleEnthalpy, {SensibleEnthalpyMassFractionFunction, SensibleEnthalpyTemperatureMassFractionFunction, ablate::eos::tChem::SensibleEnthalpy::getWorkSpaceSize}},
248  {ThermodynamicProperty::SpecificHeatConstantVolume,
249  {SpecificHeatConstantVolumeMassFractionFunction, SpecificHeatConstantVolumeTemperatureMassFractionFunction, [](auto nSpec) { return nSpec; }}},
250  {ThermodynamicProperty::SpecificHeatConstantPressure,
251  {SpecificHeatConstantPressureMassFractionFunction,
252  SpecificHeatConstantPressureTemperatureMassFractionFunction,
253  ablate::eos::tChem::Temperature::getWorkSpaceSize}} ,
254  {ThermodynamicProperty::SpeedOfSound, {SpeedOfSoundMassFractionFunction, SpeedOfSoundTemperatureMassFractionFunction, ablate::eos::tChem::SpeedOfSound::getWorkSpaceSize}},
255  {ThermodynamicProperty::SpeciesSensibleEnthalpy,
256  {SpeciesSensibleEnthalpyMassFractionFunction,
257  SpeciesSensibleEnthalpyTemperatureMassFractionFunction,
258  ablate::eos::tChem::Temperature::getWorkSpaceSize}}
259  };
260 
264  const std::set<ThermodynamicProperty> speciesSizedProperties = {ThermodynamicProperty::SpeciesSensibleEnthalpy};
265 
271  static void FillWorkingVectorFromDensityMassFractions(double density, double temperature, const double* densityYi, const tChemLib::Impl::StateVector<real_type_1d_view_host>& stateVector);
272 
278  static void FillWorkingVectorFromMassFractions(double density, double temperature, const double* densityYi, const tChemLib::Impl::StateVector<real_type_1d_view_host>& stateVector);
279 };
280 
281 } // namespace ablate::eos
282 #endif // ABLATELIBRARY_TCHEM_HPP
Definition: tChemBase.hpp:25
const std::filesystem::path mechanismFile
the mechanismFile may be chemkin or yaml based
Definition: tChemBase.hpp:31
Definition: tChem.hpp:26
PetscErrorCode(*)(const PetscReal conserved[], PetscReal *property, void *ctx) ThermodynamicStaticFunction
Definition: tChem.hpp:211
ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: tChem.cpp:73
std::shared_ptr< SourceCalculator > CreateSourceCalculator(const std::vector< domain::Field > &fields, const ablate::domain::Range &cellRange) override
Definition: tChem.cpp:934
ThermodynamicMassFractionFunction GetThermodynamicMassFractionFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const
Definition: tChem.cpp:79
std::shared_ptr< FunctionContext > BuildFunctionContext(ablate::eos::ThermodynamicProperty property, const std::vector< domain::Field > &fields, bool checkDensityYi=true) const
Definition: tChem.cpp:16
PetscErrorCode(*)(const PetscReal conserved[], const PetscReal yi[], PetscReal *property, void *ctx) ThermodynamicStaticMassFractionFunction
Definition: tChem.hpp:235
static void FillWorkingVectorFromMassFractions(double density, double temperature, const double *densityYi, const tChemLib::Impl::StateVector< real_type_1d_view_host > &stateVector)
Definition: tChem.cpp:604
std::map< std::string, double > GetElementInformation() const override
Definition: tChem.cpp:882
ThermodynamicTemperatureMassFractionFunction GetThermodynamicTemperatureMassFractionFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: tChem.cpp:86
const std::set< ThermodynamicProperty > speciesSizedProperties
Definition: tChem.hpp:264
TChem(std::filesystem::path mechanismFile, std::shared_ptr< ablate::monitors::logs::Log >={}, const std::shared_ptr< ablate::parameters::Parameters > &options={})
Definition: tChem.cpp:13
static void FillWorkingVectorFromDensityMassFractions(double density, double temperature, const double *densityYi, const tChemLib::Impl::StateVector< real_type_1d_view_host > &stateVector)
Definition: tChem.cpp:581
std::map< std::string, std::map< std::string, int > > GetSpeciesElementalInformation() const override
Definition: tChem.cpp:897
ThermodynamicFunction GetThermodynamicFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: tChem.cpp:67
EOSFunction GetFieldFunctionFunction(const std::string &field, ThermodynamicProperty property1, ThermodynamicProperty property2, std::vector< std::string > otherProperties) const override
Definition: tChem.cpp:627
std::map< std::string, double > GetSpeciesMolecularMass() const override
Definition: tChem.cpp:922
Definition: range.hpp:11
std::shared_ptr< void > context
optional context to pass into the function
Definition: tChem.hpp:35
PetscInt propertySize
the property size being set
Definition: tChem.hpp:37