ABLATE Source Documentation  0.12.33
tChemBase.hpp
1 #ifndef ABLATELIBRARY_TCHEMBASE_HPP
2 #define ABLATELIBRARY_TCHEMBASE_HPP
3 
4 #include <filesystem>
5 #include <map>
6 #include <memory>
7 #include "TChem_KineticModelData.hpp"
8 #include "chemistryModel.hpp"
9 #include "eos.hpp"
10 #include "eos/tChem/pressure.hpp"
11 #include "eos/tChem/sensibleEnthalpy.hpp"
12 #include "eos/tChem/sensibleInternalEnergy.hpp"
13 #include "eos/tChem/sourceCalculator.hpp"
14 #include "eos/tChem/speedOfSound.hpp"
15 #include "eos/tChem/temperature.hpp"
16 #include "monitors/logs/log.hpp"
17 #include "parameters/parameters.hpp"
18 #include "utilities/intErrorChecker.hpp"
19 #include "utilities/vectorUtilities.hpp"
20 
21 namespace ablate::eos {
22 
23 namespace tChemLib = TChem;
24 
25 class TChemBase : public ChemistryModel {
26  protected:
29 
31  const std::filesystem::path mechanismFile;
32 
34  std::shared_ptr<ablate::monitors::logs::Log> log;
35 
39  tChemLib::KineticModelData kineticsModel;
40 
44  std::shared_ptr<tChemLib::KineticModelGasConstData<typename Tines::UseThisDevice<host_exec_space>::type>> kineticsModelDataHost;
45 
49  std::shared_ptr<tChemLib::KineticModelGasConstData<typename Tines::UseThisDevice<exec_space>::type>> kineticsModelDataDevice;
50 
54  std::vector<std::string> species;
55 
59  real_type_1d_view enthalpyReferenceDevice;
60 
64  real_type_1d_view_host enthalpyReferenceHost;
65 
66  public:
72  explicit TChemBase(const std::string& eosName, std::filesystem::path mechanismFile, const std::shared_ptr<ablate::monitors::logs::Log>& = {},
73  const std::shared_ptr<ablate::parameters::Parameters>& options = {});
74 
79  void View(std::ostream& stream) const override;
80 
84  tChemLib::KineticModelData& GetKineticModelData() { return kineticsModel; }
85 
89  real_type_1d_view GetEnthalpyOfFormation() { return enthalpyReferenceDevice; };
90 
96  [[nodiscard]] const std::vector<std::string>& GetSpeciesVariables() const override { return species; }
97 
102  [[nodiscard]] virtual const std::vector<std::string>& GetProgressVariables() const override { return ablate::utilities::VectorUtilities::Empty<std::string>; }
103 
109  inline double GetEnthalpyOfFormation(std::string_view speciesName) const override {
110  auto it = std::find_if(species.begin(), species.end(), [speciesName](const auto& component) { return component == speciesName; });
111  // If element was found
112  if (it != species.end()) {
113  return enthalpyReferenceHost[std::distance(species.begin(), it)];
114  } else {
115  throw std::invalid_argument(std::string("Cannot locate species ") + std::string(speciesName) + " in EOS " + type);
116  }
117  }
118 
119  protected:
123  static const inline std::array<std::string, 2> validFileExtensions = {".yaml", ".yml"};
124 
126  // memory access locations for fields
127  PetscInt dim;
128  PetscInt eulerOffset;
129  PetscInt densityYiOffset;
130 
132  real_type_2d_view_host stateHost;
134  real_type_2d_view_host perSpeciesHost;
136  real_type_1d_view_host mixtureHost;
137 
139  real_type_1d_view_host enthalpyReferenceHost;
140 
142  tChemLib::UseThisTeamPolicy<tChemLib::host_exec_space>::type policy;
143 
145  std::shared_ptr<tChemLib::KineticModelGasConstData<typename Tines::UseThisDevice<host_exec_space>::type>> kineticsModelDataHost;
146  };
147 
148  public:
149  // Private static helper functions
150  inline const static double TREF = 298.15;
151 };
152 
153 } // namespace ablate::eos
154 #endif // ABLATELIBRARY_TCHEMBASE_HPP
Definition: chemistryModel.hpp:18
Definition: tChemBase.hpp:25
real_type_1d_view GetEnthalpyOfFormation()
Definition: tChemBase.hpp:89
tChemLib::KineticModelData & GetKineticModelData()
Definition: tChemBase.hpp:84
tChemLib::KineticModelData kineticsModel
Definition: tChemBase.hpp:39
std::shared_ptr< tChemLib::KineticModelGasConstData< typename Tines::UseThisDevice< host_exec_space >::type > > kineticsModelDataHost
Definition: tChemBase.hpp:44
virtual const std::vector< std::string > & GetProgressVariables() const override
Definition: tChemBase.hpp:102
tChem::SourceCalculator::ChemistryConstraints constraints
hold a copy of the constrains that can be used for single or batch source calculation
Definition: tChemBase.hpp:28
std::shared_ptr< tChemLib::KineticModelGasConstData< typename Tines::UseThisDevice< exec_space >::type > > kineticsModelDataDevice
Definition: tChemBase.hpp:49
void View(std::ostream &stream) const override
Definition: tChemBase.cpp:79
real_type_1d_view_host enthalpyReferenceHost
Definition: tChemBase.hpp:64
std::shared_ptr< ablate::monitors::logs::Log > log
an optional log file for tchem echo redirection
Definition: tChemBase.hpp:34
std::vector< std::string > species
Definition: tChemBase.hpp:54
double GetEnthalpyOfFormation(std::string_view speciesName) const override
Definition: tChemBase.hpp:109
static const std::array< std::string, 2 > validFileExtensions
Definition: tChemBase.hpp:123
const std::filesystem::path mechanismFile
the mechanismFile may be chemkin or yaml based
Definition: tChemBase.hpp:31
const std::vector< std::string > & GetSpeciesVariables() const override
Definition: tChemBase.hpp:96
TChemBase(const std::string &eosName, std::filesystem::path mechanismFile, const std::shared_ptr< ablate::monitors::logs::Log > &={}, const std::shared_ptr< ablate::parameters::Parameters > &options={})
Definition: tChemBase.cpp:10
real_type_1d_view enthalpyReferenceDevice
Definition: tChemBase.hpp:59
Definition: tChemBase.hpp:125
tChemLib::UseThisTeamPolicy< tChemLib::host_exec_space >::type policy
the kokkos team policy for this function
Definition: tChemBase.hpp:142
real_type_1d_view_host enthalpyReferenceHost
store the enthalpyReferencePerSpecies
Definition: tChemBase.hpp:139
real_type_2d_view_host stateHost
per species state
Definition: tChemBase.hpp:132
real_type_2d_view_host perSpeciesHost
per species array
Definition: tChemBase.hpp:134
real_type_1d_view_host mixtureHost
mass weighted mixture
Definition: tChemBase.hpp:136
std::shared_ptr< tChemLib::KineticModelGasConstData< typename Tines::UseThisDevice< host_exec_space >::type > > kineticsModelDataHost
the kinetics data
Definition: tChemBase.hpp:145
hold a struct that can be used for chemistry constraints
Definition: sourceCalculator.hpp:26