1 #ifndef ABLATELIBRARY_CHEMTAB_HPP
2 #define ABLATELIBRARY_CHEMTAB_HPP
7 #include "chemistryModel.hpp"
8 #include "eos/tChem.hpp"
10 #include <tensorflow/c/c_api.h>
11 #include "finiteVolume/compressibleFlowFields.hpp"
12 #include "utilities/vectorUtilities.hpp"
15 namespace ablate::eos {
17 #ifdef WITH_TENSORFLOW
18 class ChemTab :
public ChemistryModel,
public std::enable_shared_from_this<ChemTab>,
public utilities::Loggable<ChemTab> {
20 inline const static std::string DENSITY_YI_DECODE_FIELD =
"DENSITY_YI_DECODE";
24 std::shared_ptr<ablate::eos::TChem> referenceEOS;
27 TF_Graph* graph =
nullptr;
28 TF_Status* status =
nullptr;
29 TF_SessionOptions* sessionOpts =
nullptr;
30 TF_Buffer* runOpts =
nullptr;
31 TF_Session* session =
nullptr;
32 std::vector<std::string> speciesNames = std::vector<std::string>(0);
33 std::vector<std::string> progressVariablesNames = std::vector<std::string>(0);
35 PetscReal** Wmat =
nullptr;
36 PetscReal* sourceEnergyScaler =
nullptr;
39 std::map<std::string, std::map<std::string, double>> initializers;
44 void ExtractMetaData(std::istream& inputStream);
45 static void LoadBasisVectors(std::istream& inputStream, std::size_t columns, PetscReal** W);
55 void ChemTabModelComputeFunction(PetscReal density,
const PetscReal densityProgressVariable[], PetscReal* predictedSourceEnergy, PetscReal* progressVariableSource,
56 PetscReal* densityMassFractions)
const;
59 [[nodiscard]] std::vector<std::string>
GetFieldTags()
const override {
return std::vector<std::string>{ablate::finiteVolume::CompressibleFlowFields::MinusOneToOneRange}; }
65 class ChemTabSourceCalculator :
public ChemistryModel::SourceCalculator {
68 const PetscInt densityOffset;
69 const PetscInt densityEnergyOffset;
70 const PetscInt densityProgressVariableOffset;
72 const std::shared_ptr<ChemTab> chemTabModel;
75 ChemTabSourceCalculator(PetscInt densityOffset, PetscInt densityEnergyOffset, PetscInt densityProgressVariableOffset, std::shared_ptr<ChemTab> chemTabModel);
80 void ComputeSource(
const ablate::domain::Range& cellRange, PetscReal time, PetscReal dt, Vec solution)
override{};
95 void ComputeProgressVariables(
const PetscReal* massFractions, PetscReal* progressVariables)
const;
107 static PetscErrorCode ComputeMassFractions(PetscReal time, PetscInt dim,
const PetscFVCellGeom* cellGeom,
const PetscInt uOff[], PetscScalar* u,
void* ctx);
110 explicit ChemTab(
const std::filesystem::path& path);
117 [[nodiscard]]
const std::vector<std::string>&
GetSpeciesVariables()
const override {
return ablate::utilities::VectorUtilities::Empty<std::string>; }
123 [[nodiscard]]
const std::vector<std::string>& GetSpeciesNames()
const {
return speciesNames; }
135 [[nodiscard]]
const std::vector<std::string>&
GetProgressVariables()
const override {
return progressVariablesNames; }
141 std::vector<std::shared_ptr<domain::FieldDescriptor>>
GetAdditionalFields()
const override;
149 void ChemistrySource(PetscReal density,
const PetscReal densityProgressVariable[], PetscReal* densityEnergySource, PetscReal* progressVariableSource)
const;
164 void View(std::ostream& stream)
const override;
172 [[nodiscard]] eos::ThermodynamicFunction
GetThermodynamicFunction(eos::ThermodynamicProperty property,
const std::vector<domain::Field>& fields)
const override;
180 [[nodiscard]] eos::ThermodynamicTemperatureFunction
GetThermodynamicTemperatureFunction(eos::ThermodynamicProperty property,
const std::vector<domain::Field>& fields)
const override;
188 [[nodiscard]] eos::EOSFunction
GetFieldFunctionFunction(
const std::string& field, eos::ThermodynamicProperty property1, eos::ThermodynamicProperty property2,
189 std::vector<std::string> otherProperties)
const override;
198 void ComputeProgressVariables(
const std::vector<PetscReal>& massFractions, std::vector<PetscReal>& progressVariables)
const;
209 void ComputeMassFractions(
const std::vector<PetscReal>& progressVariables, std::vector<PetscReal>& massFractions, PetscReal density = 1.0)
const;
217 void ComputeMassFractions(
const PetscReal* progressVariables, PetscReal* massFractions, PetscReal density = 1.0)
const;
224 void GetInitializerProgressVariables(
const std::string& name, std::vector<PetscReal>& progressVariables)
const;
230 [[nodiscard]] std::vector<std::tuple<ablate::solver::CellSolver::SolutionFieldUpdateFunction, void*, std::vector<std::string>>>
GetSolutionFieldUpdates()
override;
236 inline const static std::string DENSITY_YI_DECODE_FIELD =
"DENSITY_YI_DECODE";
237 static inline const std::string errorMessage =
"Using the ChemTab requires Tensorflow to be compile with ABLATE.";
238 ChemTab(std::filesystem::path path) :
ChemistryModel(
"ablate::chemistry::ChemTabModel") {
throw std::runtime_error(errorMessage); }
240 [[nodiscard]]
const std::vector<std::string>&
GetSpeciesVariables()
const override {
throw std::runtime_error(errorMessage); }
242 [[nodiscard]]
const std::vector<std::string>& GetSpeciesNames()
const {
throw std::runtime_error(errorMessage); }
246 [[nodiscard]]
const std::vector<std::string>&
GetProgressVariables()
const override {
throw std::runtime_error(errorMessage); }
248 [[nodiscard]] std::vector<std::string>
GetFieldTags()
const override {
throw std::runtime_error(errorMessage); }
252 void View(std::ostream& stream)
const override {
throw std::runtime_error(errorMessage); }
255 throw std::runtime_error(errorMessage);
259 throw std::runtime_error(errorMessage);
262 [[nodiscard]] eos::EOSFunction
GetFieldFunctionFunction(
const std::string& field, eos::ThermodynamicProperty property1, eos::ThermodynamicProperty property2,
263 std::vector<std::string> otherProperties)
const override {
264 throw std::runtime_error(errorMessage);
267 [[nodiscard]] ThermodynamicTemperatureMassFractionFunction GetThermodynamicTemperatureMassFractionFunction(ThermodynamicProperty property,
268 const std::vector<domain::Field>& fields)
const override {
269 throw std::runtime_error(errorMessage);
271 std::map<std::string, double> GetSpeciesMolecularMass()
const override {
throw std::runtime_error(errorMessage); }
272 std::map<std::string, std::map<std::string, int>> GetSpeciesElementalInformation()
const override {
throw std::runtime_error(errorMessage); }
273 std::map<std::string, double> GetElementInformation()
const override {
throw std::runtime_error(errorMessage); }
275 void ChemistrySource(PetscReal density,
const PetscReal densityProgressVariable[], PetscReal* densityEnergySource, PetscReal* progressVariableSource)
const {
276 throw std::runtime_error(errorMessage);
279 void ComputeProgressVariables(
const std::vector<PetscReal>& massFractions, std::vector<PetscReal>& progressVariables)
const {
throw std::runtime_error(errorMessage); }
281 void ComputeMassFractions(
const std::vector<PetscReal>& progressVariables, std::vector<PetscReal>& massFractions, PetscReal density = 1.0)
const {
throw std::runtime_error(errorMessage); }
283 void ComputeMassFractions(
const PetscReal* progressVariables, PetscReal* massFractions, PetscReal density = 1.0)
const {
throw std::runtime_error(errorMessage); }
285 void GetInitializerProgressVariables(
const std::string& name, std::vector<PetscReal>& progressVariables)
const {
throw std::runtime_error(errorMessage); }
287 [[nodiscard]] std::vector<std::tuple<ablate::solver::CellSolver::SolutionFieldUpdateFunction, void*, std::vector<std::string>>>
GetSolutionFieldUpdates()
override {
288 throw std::runtime_error(errorMessage);
Definition: chemTab.hpp:234
const std::vector< std::string > & GetSpeciesVariables() const override
Definition: chemTab.hpp:240
std::vector< std::tuple< ablate::solver::CellSolver::SolutionFieldUpdateFunction, void *, std::vector< std::string > > > GetSolutionFieldUpdates() override
Definition: chemTab.hpp:287
void View(std::ostream &stream) const override
Definition: chemTab.hpp:252
eos::ThermodynamicFunction GetThermodynamicFunction(eos::ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: chemTab.hpp:254
std::shared_ptr< SourceCalculator > CreateSourceCalculator(const std::vector< domain::Field > &fields, const ablate::domain::Range &cellRange) override
Definition: chemTab.hpp:250
std::vector< std::string > GetFieldTags() const override
Definition: chemTab.hpp:248
const std::vector< std::string > & GetProgressVariables() const override
Definition: chemTab.hpp:246
const std::vector< std::string > & GetFieldFunctionProperties() const override
Definition: chemTab.hpp:244
eos::EOSFunction GetFieldFunctionFunction(const std::string &field, eos::ThermodynamicProperty property1, eos::ThermodynamicProperty property2, std::vector< std::string > otherProperties) const override
Definition: chemTab.hpp:262
eos::ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(eos::ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: chemTab.hpp:258
Definition: chemistryModel.hpp:18
ChemistryModel(std::string name)
Definition: chemistryModel.hpp:24
virtual std::vector< std::shared_ptr< domain::FieldDescriptor > > GetAdditionalFields() const
Definition: eos.hpp:128