ABLATE Source Documentation  0.12.36
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 
37  // Store any initializers specified by the metadata
38  std::map<std::string, std::map<std::string, double>> initializers;
39 
43  void ExtractMetaData(std::istream& inputStream);
44  static void LoadBasisVectors(std::istream& inputStream, std::size_t columns, PetscReal** W);
45 
54  void ChemTabModelComputeFunction(PetscReal density, const PetscReal densityProgressVariables[], PetscReal* densityEnergySource, PetscReal* densityProgressVariableSource,
55  PetscReal* densityMassFractions) const;
56 
67  void ChemTabModelComputeFunction(const PetscReal density[], const PetscReal* const* const densityProgressVariables, PetscReal** densityEnergySource, PetscReal** densityProgressVariableSource,
68  PetscReal** densityMassFractions, size_t batch_size) const;
69 
71  [[nodiscard]] std::vector<std::string> GetFieldTags() const override { return std::vector<std::string>{ablate::finiteVolume::CompressibleFlowFields::MinusOneToOneRange}; }
72 
77  class ChemTabSourceCalculator : public ChemistryModel::SourceCalculator {
78  private:
80  const PetscInt densityOffset;
81  const PetscInt densityEnergyOffset;
82  const PetscInt densityProgressVariableOffset;
84  const std::shared_ptr<ChemTab> chemTabModel;
85 
86  public:
87  ChemTabSourceCalculator(PetscInt densityOffset, PetscInt densityEnergyOffset, PetscInt densityProgressVariableOffset, std::shared_ptr<ChemTab> chemTabModel);
88 
92  void ComputeSource(const ablate::domain::Range& cellRange, PetscReal time, PetscReal dt, Vec solution) override{};
93 
97  void AddSource(const ablate::domain::Range& cellRange, Vec locSolution, Vec locSource) override;
98  };
99 
105  void ComputeProgressVariables(const PetscReal* massFractions, PetscReal* progressVariables) const;
106 
114  void ComputeProgressVariables(const PetscReal* const* massFractions, PetscReal* const* progressVariables, size_t n) const;
115 
126  static PetscErrorCode ComputeMassFractions(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], PetscScalar* u, void* ctx);
127 
128  public:
129  explicit ChemTab(const std::filesystem::path& path);
130  ~ChemTab() override;
131 
136  [[nodiscard]] const std::vector<std::string>& GetSpeciesVariables() const override { return ablate::utilities::VectorUtilities::Empty<std::string>; }
137 
142  [[nodiscard]] const std::vector<std::string>& GetSpeciesNames() const { return speciesNames; }
143 
148  [[nodiscard]] const std::vector<std::string>& GetFieldFunctionProperties() const override { return progressVariablesNames; }
149 
154  [[nodiscard]] const std::vector<std::string>& GetProgressVariables() const override { return progressVariablesNames; }
155 
160  std::vector<std::shared_ptr<domain::FieldDescriptor>> GetAdditionalFields() const override;
161 
169  void ChemistrySource(const PetscReal density, const PetscReal densityProgressVariable[], PetscReal* densityEnergySource, PetscReal* densityProgressVariableSource) const;
170 
180  void ChemistrySource(const PetscReal* const density, const PetscReal* const* const densityProgressVariable, PetscReal** densityEnergySource, PetscReal** progressVariableSource, size_t n) const;
181 
189  std::shared_ptr<SourceCalculator> CreateSourceCalculator(const std::vector<domain::Field>& fields, const ablate::domain::Range& cellRange) override;
190 
195  void View(std::ostream& stream) const override;
196 
203  [[nodiscard]] eos::ThermodynamicFunction GetThermodynamicFunction(eos::ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override;
204 
211  [[nodiscard]] eos::ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(eos::ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override;
212 
219  [[nodiscard]] eos::EOSFunction GetFieldFunctionFunction(const std::string& field, eos::ThermodynamicProperty property1, eos::ThermodynamicProperty property2,
220  std::vector<std::string> otherProperties) const override;
221 
229  void ComputeProgressVariables(const std::vector<PetscReal>& massFractions, std::vector<PetscReal>& progressVariables) const;
230 
239  void ComputeMassFractions(std::vector<PetscReal>& progressVariables, std::vector<PetscReal>& massFractions, PetscReal density = 1.0) const;
240 
247  void ComputeMassFractions(const PetscReal* densityProgressVariables, PetscReal* densityMassFractions, const PetscReal density = 1.0) const;
248 
256  void ComputeMassFractions(const PetscReal* const* densityProgressVariables, PetscReal** densityMassFractions, const PetscReal density[], size_t n) const;
257 
263  void GetInitializerProgressVariables(const std::string& name, std::vector<PetscReal>& progressVariables) const;
264 
269  [[nodiscard]] std::vector<std::tuple<ablate::solver::CellSolver::SolutionFieldUpdateFunction, void*, std::vector<std::string>>> GetSolutionFieldUpdates() override;
270 
271  void ExtractModelOutputsAtPoint(const PetscReal density, PetscReal* densityEnergySource, PetscReal* densityProgressVariableSource, PetscReal* densityMassFractions,
272  const std::array<TF_Tensor*, 2>& outputValues, size_t id = 0) const;
273 };
274 
275 #else
276 class ChemTab : public ChemistryModel {
277  public:
278  inline const static std::string DENSITY_YI_DECODE_FIELD = "DENSITY_YI_DECODE";
279  static inline const std::string errorMessage = "Using the ChemTab requires Tensorflow to be compile with ABLATE.";
280  ChemTab(std::filesystem::path path) : ChemistryModel("ablate::chemistry::ChemTabModel") { throw std::runtime_error(errorMessage); }
281 
282  [[nodiscard]] const std::vector<std::string>& GetSpeciesVariables() const override { throw std::runtime_error(errorMessage); }
283 
284  [[nodiscard]] const std::vector<std::string>& GetSpeciesNames() const { throw std::runtime_error(errorMessage); }
285 
286  [[nodiscard]] const std::vector<std::string>& GetFieldFunctionProperties() const override { throw std::runtime_error(errorMessage); }
287 
288  [[nodiscard]] const std::vector<std::string>& GetProgressVariables() const override { throw std::runtime_error(errorMessage); }
289 
290  [[nodiscard]] std::vector<std::string> GetFieldTags() const override { throw std::runtime_error(errorMessage); }
291 
292  std::shared_ptr<SourceCalculator> CreateSourceCalculator(const std::vector<domain::Field>& fields, const ablate::domain::Range& cellRange) override { throw std::runtime_error(errorMessage); }
293 
294  void View(std::ostream& stream) const override { throw std::runtime_error(errorMessage); }
295 
296  [[nodiscard]] eos::ThermodynamicFunction GetThermodynamicFunction(eos::ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override {
297  throw std::runtime_error(errorMessage);
298  }
299 
300  [[nodiscard]] eos::ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(eos::ThermodynamicProperty property, const std::vector<domain::Field>& fields) const override {
301  throw std::runtime_error(errorMessage);
302  }
303 
304  [[nodiscard]] eos::EOSFunction GetFieldFunctionFunction(const std::string& field, eos::ThermodynamicProperty property1, eos::ThermodynamicProperty property2,
305  std::vector<std::string> otherProperties) const override {
306  throw std::runtime_error(errorMessage);
307  }
308 
309  [[nodiscard]] ThermodynamicTemperatureMassFractionFunction GetThermodynamicTemperatureMassFractionFunction(ThermodynamicProperty property,
310  const std::vector<domain::Field>& fields) const override {
311  throw std::runtime_error(errorMessage);
312  }
313  std::map<std::string, double> GetSpeciesMolecularMass() const override { throw std::runtime_error(errorMessage); }
314  std::map<std::string, std::map<std::string, int>> GetSpeciesElementalInformation() const override { throw std::runtime_error(errorMessage); }
315  std::map<std::string, double> GetElementInformation() const override { throw std::runtime_error(errorMessage); }
316 
317  void ChemistrySource(PetscReal density, const PetscReal densityProgressVariable[], PetscReal* densityEnergySource, PetscReal* progressVariableSource) const {
318  throw std::runtime_error(errorMessage);
319  }
320 
321  void ComputeProgressVariables(const std::vector<PetscReal>& massFractions, std::vector<PetscReal>& progressVariables) const { throw std::runtime_error(errorMessage); }
322 
323  void ComputeMassFractions(const std::vector<PetscReal>& progressVariables, std::vector<PetscReal>& massFractions, PetscReal density = 1.0) const { throw std::runtime_error(errorMessage); }
324 
325  void ComputeMassFractions(const PetscReal* progressVariables, PetscReal* massFractions, PetscReal density = 1.0) const { throw std::runtime_error(errorMessage); }
326 
327  void GetInitializerProgressVariables(const std::string& name, std::vector<PetscReal>& progressVariables) const { throw std::runtime_error(errorMessage); }
328 
329  [[nodiscard]] std::vector<std::tuple<ablate::solver::CellSolver::SolutionFieldUpdateFunction, void*, std::vector<std::string>>> GetSolutionFieldUpdates() override {
330  throw std::runtime_error(errorMessage);
331  }
332 };
333 #endif
334 } // namespace ablate::eos
335 #endif // ABLATELIBRARY_CHEMTAB_HPP
Definition: chemTab.hpp:276
const std::vector< std::string > & GetSpeciesVariables() const override
Definition: chemTab.hpp:282
std::vector< std::tuple< ablate::solver::CellSolver::SolutionFieldUpdateFunction, void *, std::vector< std::string > > > GetSolutionFieldUpdates() override
Definition: chemTab.hpp:329
void View(std::ostream &stream) const override
Definition: chemTab.hpp:294
eos::ThermodynamicFunction GetThermodynamicFunction(eos::ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: chemTab.hpp:296
std::shared_ptr< SourceCalculator > CreateSourceCalculator(const std::vector< domain::Field > &fields, const ablate::domain::Range &cellRange) override
Definition: chemTab.hpp:292
std::vector< std::string > GetFieldTags() const override
Definition: chemTab.hpp:290
const std::vector< std::string > & GetProgressVariables() const override
Definition: chemTab.hpp:288
const std::vector< std::string > & GetFieldFunctionProperties() const override
Definition: chemTab.hpp:286
eos::EOSFunction GetFieldFunctionFunction(const std::string &field, eos::ThermodynamicProperty property1, eos::ThermodynamicProperty property2, std::vector< std::string > otherProperties) const override
Definition: chemTab.hpp:304
eos::ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(eos::ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: chemTab.hpp:300
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