ABLATE Source Documentation  0.12.33
chemistryModel.hpp
1 #ifndef ABLATELIBRARY_CHEMISTRYMODEL_HPP
2 #define ABLATELIBRARY_CHEMISTRYMODEL_HPP
3 
4 #include <petsc.h>
5 #include <memory>
6 #include <string>
7 #include <utility>
8 #include <vector>
9 #include "eos/eos.hpp"
10 #include "solver/cellSolver.hpp"
11 #include "solver/solver.hpp"
12 
13 namespace ablate::eos {
14 
18 class ChemistryModel : public eos::EOS {
19  public:
24  explicit ChemistryModel(std::string name) : eos::EOS(std::move(name)){};
25 
31  public:
32  virtual ~SourceCalculator() = default;
36  virtual void ComputeSource(const ablate::domain::Range& cellRange, PetscReal time, PetscReal dt, Vec solution) = 0;
37 
41  virtual void AddSource(const ablate::domain::Range& cellRange, Vec solution, Vec source) = 0;
42  };
43 
49  virtual std::shared_ptr<SourceCalculator> CreateSourceCalculator(const std::vector<domain::Field>& fields, const ablate::domain::Range& cellRange) = 0;
50 
54  virtual std::vector<std::tuple<ablate::solver::CellSolver::SolutionFieldUpdateFunction, void*, std::vector<std::string>>> GetSolutionFieldUpdates() { return {}; }
55 
56  virtual inline double GetEnthalpyOfFormation(std::string_view speciesName) const { return {}; };
57 
58  [[nodiscard]] virtual std::map<std::string, double> GetSpeciesMolecularMass() const { return {}; };
59 
60  [[nodiscard]] virtual std::map<std::string, double> GetElementInformation() const { return {}; };
61 
62  [[nodiscard]] virtual std::map<std::string, std::map<std::string, int>> GetSpeciesElementalInformation() const { return {}; };
68  PetscErrorCode (*function)(const PetscReal conserved[], const PetscReal yi[], PetscReal T, PetscReal* property, void* ctx) = nullptr;
70  std::shared_ptr<void> context = nullptr;
72  PetscInt propertySize = 1;
73  };
74 
75  [[nodiscard]] virtual ThermodynamicTemperatureMassFractionFunction GetThermodynamicTemperatureMassFractionFunction(ThermodynamicProperty property, const std::vector<domain::Field>& fields) const {
76  return {};
77  };
78 };
79 } // namespace ablate::eos
80 
81 #endif // ABLATELIBRARY_CHEMISTRYMODEL_HPP
Definition: chemistryModel.hpp:30
virtual void AddSource(const ablate::domain::Range &cellRange, Vec solution, Vec source)=0
virtual void ComputeSource(const ablate::domain::Range &cellRange, PetscReal time, PetscReal dt, Vec solution)=0
Definition: chemistryModel.hpp:18
virtual std::vector< std::tuple< ablate::solver::CellSolver::SolutionFieldUpdateFunction, void *, std::vector< std::string > > > GetSolutionFieldUpdates()
Definition: chemistryModel.hpp:54
ChemistryModel(std::string name)
Definition: chemistryModel.hpp:24
virtual std::shared_ptr< SourceCalculator > CreateSourceCalculator(const std::vector< domain::Field > &fields, const ablate::domain::Range &cellRange)=0
Definition: eos.hpp:60
Definition: range.hpp:11
std::shared_ptr< void > context
optional context to pass into the function
Definition: chemistryModel.hpp:70
PetscInt propertySize
the property size being set
Definition: chemistryModel.hpp:72