ABLATE Source Documentation  0.12.34
chemTab.hpp
1 #ifndef ABLATELIBRARY_CHEMTAB_HPP
2 #define ABLATELIBRARY_CHEMTAB_HPP
3 
4 #include <petscmat.h>
5 #include <filesystem>
6 #include <istream>
7 #include "chemistryModel.hpp"
8 #include "eos/tChem.hpp"
9 #ifdef WITH_TENSORFLOW
10 #include <tensorflow/c/c_api.h>
11 #include "finiteVolume/compressibleFlowFields.hpp"
12 #include "utilities/vectorUtilities.hpp"
13 #endif
14 
15 namespace ablate::eos {
16 
17 #ifdef WITH_TENSORFLOW
18 class ChemTab : public ChemistryModel, public std::enable_shared_from_this<ChemTab>, public utilities::Loggable<ChemTab> {
19  public:
20  inline const static std::string DENSITY_YI_DECODE_FIELD = "DENSITY_YI_DECODE";
21 
22  private:
24  std::shared_ptr<ablate::eos::TChem> referenceEOS;
25 
26  // hold the required tensorflow information
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);
34 
35  PetscReal** Wmat = nullptr;
36  PetscReal* sourceEnergyScaler = nullptr;
37 
38  // Store any initializers specified by the metadata
39  std::map<std::string, std::map<std::string, double>> initializers;
40 
44  void ExtractMetaData(std::istream& inputStream);
45  static void LoadBasisVectors(std::istream& inputStream, std::size_t columns, PetscReal** W);
46 
55  void ChemTabModelComputeFunction(PetscReal density, const PetscReal densityProgressVariable[], PetscReal* predictedSourceEnergy, PetscReal* progressVariableSource,
56  PetscReal* densityMassFractions) const;
57 
59  [[nodiscard]] std::vector<std::string> GetFieldTags() const override { return std::vector<std::string>{ablate::finiteVolume::CompressibleFlowFields::MinusOneToOneRange}; }
60 
65  class ChemTabSourceCalculator : public ChemistryModel::SourceCalculator {
66  private:
68  const PetscInt densityOffset;
69  const PetscInt densityEnergyOffset;
70  const PetscInt densityProgressVariableOffset;
72  const std::shared_ptr<ChemTab> chemTabModel;
73 
74  public:
75  ChemTabSourceCalculator(PetscInt densityOffset, PetscInt densityEnergyOffset, PetscInt densityProgressVariableOffset, std::shared_ptr<ChemTab> chemTabModel);
76 
80  void ComputeSource(const ablate::domain::Range& cellRange, PetscReal time, PetscReal dt, Vec solution) override{};
81 
85  void AddSource(const ablate::domain::Range& cellRange, Vec locSolution, Vec locSource) override;
86  };
87 
95  void ComputeProgressVariables(const PetscReal* massFractions, PetscReal* progressVariables) const;
96 
107  static PetscErrorCode ComputeMassFractions(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], PetscScalar* u, void* ctx);
108 
109  public:
110  explicit ChemTab(const std::filesystem::path& path);
111  ~ChemTab() override;
112 
117  [[nodiscard]] const std::vector<std::string>& GetSpeciesVariables() const override { return ablate::utilities::VectorUtilities::Empty<std::string>; }
118 
123  [[nodiscard]] const std::vector<std::string>& GetSpeciesNames() const { return speciesNames; }
124 
129  [[nodiscard]] const std::vector<std::string>& GetFieldFunctionProperties() const override { return progressVariablesNames; }
130 
135  [[nodiscard]] const std::vector<std::string>& GetProgressVariables() const override { return progressVariablesNames; }
136 
141  std::vector<std::shared_ptr<domain::FieldDescriptor>> GetAdditionalFields() const override;
142 
149  void ChemistrySource(PetscReal density, const PetscReal densityProgressVariable[], PetscReal* densityEnergySource, PetscReal* progressVariableSource) const;
150 
158  std::shared_ptr<SourceCalculator> CreateSourceCalculator(const std::vector<domain::Field>& fields, const ablate::domain::Range& cellRange) override;
159 
164  void View(std::ostream& stream) const override;
165 
172  [[nodiscard]] eos::ThermodynamicFunction GetThermodynamicFunction(eos::ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override;
173 
180  [[nodiscard]] eos::ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(eos::ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override;
181 
188  [[nodiscard]] eos::EOSFunction GetFieldFunctionFunction(const std::string& field, eos::ThermodynamicProperty property1, eos::ThermodynamicProperty property2,
189  std::vector<std::string> otherProperties) const override;
190 
198  void ComputeProgressVariables(const std::vector<PetscReal>& massFractions, std::vector<PetscReal>& progressVariables) const;
199 
209  void ComputeMassFractions(const std::vector<PetscReal>& progressVariables, std::vector<PetscReal>& massFractions, PetscReal density = 1.0) const;
210 
217  void ComputeMassFractions(const PetscReal* progressVariables, PetscReal* massFractions, PetscReal density = 1.0) const;
218 
224  void GetInitializerProgressVariables(const std::string& name, std::vector<PetscReal>& progressVariables) const;
225 
230  [[nodiscard]] std::vector<std::tuple<ablate::solver::CellSolver::SolutionFieldUpdateFunction, void*, std::vector<std::string>>> GetSolutionFieldUpdates() override;
231 };
232 
233 #else
234 class ChemTab : public ChemistryModel {
235  public:
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); }
239 
240  [[nodiscard]] const std::vector<std::string>& GetSpeciesVariables() const override { throw std::runtime_error(errorMessage); }
241 
242  [[nodiscard]] const std::vector<std::string>& GetSpeciesNames() const { throw std::runtime_error(errorMessage); }
243 
244  [[nodiscard]] const std::vector<std::string>& GetFieldFunctionProperties() const override { throw std::runtime_error(errorMessage); }
245 
246  [[nodiscard]] const std::vector<std::string>& GetProgressVariables() const override { throw std::runtime_error(errorMessage); }
247 
248  [[nodiscard]] std::vector<std::string> GetFieldTags() const override { throw std::runtime_error(errorMessage); }
249 
250  std::shared_ptr<SourceCalculator> CreateSourceCalculator(const std::vector<domain::Field>& fields, const ablate::domain::Range& cellRange) override { throw std::runtime_error(errorMessage); }
251 
252  void View(std::ostream& stream) const override { throw std::runtime_error(errorMessage); }
253 
254  [[nodiscard]] eos::ThermodynamicFunction GetThermodynamicFunction(eos::ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override {
255  throw std::runtime_error(errorMessage);
256  }
257 
258  [[nodiscard]] eos::ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(eos::ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override {
259  throw std::runtime_error(errorMessage);
260  }
261 
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);
265  }
266 
267  [[nodiscard]] ThermodynamicTemperatureMassFractionFunction GetThermodynamicTemperatureMassFractionFunction(ThermodynamicProperty property,
268  const std::vector<domain::Field>& fields) const override {
269  throw std::runtime_error(errorMessage);
270  }
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); }
274 
275  void ChemistrySource(PetscReal density, const PetscReal densityProgressVariable[], PetscReal* densityEnergySource, PetscReal* progressVariableSource) const {
276  throw std::runtime_error(errorMessage);
277  }
278 
279  void ComputeProgressVariables(const std::vector<PetscReal>& massFractions, std::vector<PetscReal>& progressVariables) const { throw std::runtime_error(errorMessage); }
280 
281  void ComputeMassFractions(const std::vector<PetscReal>& progressVariables, std::vector<PetscReal>& massFractions, PetscReal density = 1.0) const { throw std::runtime_error(errorMessage); }
282 
283  void ComputeMassFractions(const PetscReal* progressVariables, PetscReal* massFractions, PetscReal density = 1.0) const { throw std::runtime_error(errorMessage); }
284 
285  void GetInitializerProgressVariables(const std::string& name, std::vector<PetscReal>& progressVariables) const { throw std::runtime_error(errorMessage); }
286 
287  [[nodiscard]] std::vector<std::tuple<ablate::solver::CellSolver::SolutionFieldUpdateFunction, void*, std::vector<std::string>>> GetSolutionFieldUpdates() override {
288  throw std::runtime_error(errorMessage);
289  }
290 };
291 #endif
292 } // namespace ablate::eos
293 #endif // ABLATELIBRARY_CHEMTAB_HPP
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
Definition: range.hpp:11