ABLATE Source Documentation  0.12.33
mixtureFractionCalculator.hpp
1 #ifndef ABLATELIBRARY_MIXTUREFRACTIONCALCULATOR_HPP
2 #define ABLATELIBRARY_MIXTUREFRACTIONCALCULATOR_HPP
3 
4 #include <map>
5 #include <memory>
6 #include <string>
7 #include <vector>
8 #include "eos/tChem.hpp"
9 #include "mathFunctions/fieldFunction.hpp"
10 
11 namespace ablate::monitors {
12 
19  private:
21  const std::shared_ptr<ablate::eos::ChemistryModel> eos;
22 
24  const std::vector<std::string> trackingElements;
25 
27  double zMixFuel, zMixOxidizer;
28 
30  std::vector<double> zMixCoefficients;
31 
38  static std::map<std::string, double> ToMassFractionMap(const std::shared_ptr<ablate::eos::EOS>& eos, const std::shared_ptr<ablate::mathFunctions::FieldFunction>& massFractions);
39 
40  public:
48  MixtureFractionCalculator(const std::shared_ptr<ablate::eos::EOS>& eos, std::map<std::string, double> massFractionsFuel, std::map<std::string, double> massFractionsOxidizer,
49  const std::vector<std::string>& trackingElements = {});
50 
59  MixtureFractionCalculator(const std::shared_ptr<ablate::eos::EOS>& eos, const std::shared_ptr<ablate::mathFunctions::FieldFunction>& massFractionsFuel,
60  const std::shared_ptr<ablate::mathFunctions::FieldFunction>& massFractionsOxidizer, const std::vector<std::string>& trackingElements = {});
67  template <class T>
68  inline T Calculate(T* yi) const {
69  double zMix = 0.;
70  for (std::size_t s = 0; s < zMixCoefficients.size(); s++) {
71  zMix += zMixCoefficients[s] * yi[s];
72  }
73  return (zMix - zMixOxidizer) / (zMixFuel - zMixOxidizer);
74  }
75 
80  std::shared_ptr<ablate::eos::ChemistryModel> GetEos() { return eos; }
81 };
82 
83 } // namespace ablate::monitors
84 #endif // ABLATELIBRARY_MIXTUREFRACTIONCALCULATOR_HPP
Definition: mixtureFractionCalculator.hpp:18
MixtureFractionCalculator(const std::shared_ptr< ablate::eos::EOS > &eos, std::map< std::string, double > massFractionsFuel, std::map< std::string, double > massFractionsOxidizer, const std::vector< std::string > &trackingElements={})
Definition: mixtureFractionCalculator.cpp:4
std::shared_ptr< ablate::eos::ChemistryModel > GetEos()
Definition: mixtureFractionCalculator.hpp:80
T Calculate(T *yi) const
Definition: mixtureFractionCalculator.hpp:68
Definition: boundarySolverMonitor.hpp:16