ABLATE Source Documentation  0.12.34
sourceCalculator.hpp
1 #ifndef ABLATELIBRARY_TCHEM_SOURCECALCULATOR_HPP
2 #define ABLATELIBRARY_TCHEM_SOURCECALCULATOR_HPP
3 
4 #include <TChem_KineticModelGasConstData.hpp>
5 #include "eos/chemistryModel.hpp"
6 
7 namespace tChemLib = TChem;
8 
9 namespace ablate::eos {
10 class TChem;
11 }
12 
13 namespace ablate::eos::tChem {
14 
18 class SourceCalculator : public ChemistryModel::SourceCalculator, private utilities::Loggable<SourceCalculator> {
19  public:
23  enum class ReactorType { ConstantPressure, ConstantVolume };
24 
27  double dtMin = 1.0E-12;
28  double dtMax = 1.0E-1;
29  double dtDefault = 1E-4;
30  double dtEstimateFactor = 1.5;
31  double relToleranceTime = 1.0E-4;
32  double absToleranceTime = 1.0E-8;
33  double relToleranceNewton = 1.0E-6;
34  double absToleranceNewton = 1.0E-10;
35 
36  int maxNumNewtonIterations = 100;
37  int numTimeIterationsPerInterval = 100000;
38  int jacobianInterval = 1;
39  int maxAttempts = 4;
40 
41  // store the reactor type in the chemistry constrains
42  ReactorType reactorType = ReactorType::ConstantPressure;
43 
44  // store an optional threshold temperature. Only compute the reactions if the temperature is above thresholdTemperature
45  double thresholdTemperature = 0.0;
46 
47  void Set(const std::shared_ptr<ablate::parameters::Parameters>&);
48  };
49 
56  SourceCalculator(const std::vector<domain::Field>& fields, std::shared_ptr<TChem> tChemEos, ChemistryConstraints constraints, const ablate::domain::Range& cellRange);
57 
61  void ComputeSource(const ablate::domain::Range& cellRange, PetscReal time, PetscReal dt, Vec globalSolution) override;
62 
66  // static void ComputeSource(SourceCalculator& sourceCalculator, const ablate::domain::Range& cellRange, PetscReal time, PetscReal dt, Vec globalSolution);
67 
71  void AddSource(const ablate::domain::Range& cellRange, Vec localXVec, Vec localFVec) override;
72 
73  private:
75  ChemistryConstraints chemistryConstraints;
76 
80  std::shared_ptr<eos::TChem> eos;
81  const size_t numberSpecies;
82 
84  PetscInt eulerId;
85 
87  PetscInt densityYiId;
88 
89  // tchem memory storage on host/device. These will be sized for the number of active nodes in the domain
90  real_type_2d_view stateDevice;
91  real_type_2d_view_host stateHost;
92 
93  // store the end state for the device/host
94  real_type_2d_view endStateDevice;
95 
96  // the time advance information
97  time_advance_type_1d_view timeAdvanceDevice;
98  time_advance_type timeAdvanceDefault{};
99 
100  // store host/device memory for computing state
101  real_type_1d_view internalEnergyRefDevice;
102  real_type_1d_view_host internalEnergyRefHost;
103  real_type_2d_view perSpeciesScratchDevice;
104 
105  // store the source terms (density* energy + density*species)
106  real_type_2d_view_host sourceTermsHost;
107  real_type_2d_view sourceTermsDevice;
108 
109  // tolerance constraints
110  real_type_2d_view tolTimeDevice;
111  real_type_1d_view tolNewtonDevice;
112  real_type_2d_view facDevice;
113 
114  // store the time and delta for the ode solver
115  real_type_1d_view timeViewDevice;
116  real_type_1d_view dtViewDevice;
117 
118  // Hard code some values needed for the constant volume reactor
119  static inline constexpr bool solveTla = false; // do not calculate tangent linear approximation (TLA) for the const volume reactions
120  static inline constexpr real_type thetaTla = 0; // this is not used when solveTla is false
121 
122  // store device specific kineticModelGasConstants
123  tChemLib::KineticModelConstData<typename Tines::UseThisDevice<exec_space>::type> kineticModelGasConstDataDevice;
124  kmd_type_1d_view_host kineticModelDataClone;
125  Kokkos::View<KineticModelGasConstData<typename Tines::UseThisDevice<exec_space>::type>*, typename Tines::UseThisDevice<exec_space>::type> kineticModelGasConstDataDevices;
126 };
127 
134 std::ostream& operator<<(std::ostream& os, const SourceCalculator::ReactorType& v);
135 
142 std::istream& operator>>(std::istream& is, SourceCalculator::ReactorType& v);
143 
144 } // namespace ablate::eos::tChem
145 
146 #endif // ABLATELIBRARY_BATCHSOURCE_HPP
Definition: chemistryModel.hpp:30
Definition: sourceCalculator.hpp:18
void ComputeSource(const ablate::domain::Range &cellRange, PetscReal time, PetscReal dt, Vec globalSolution) override
Definition: sourceCalculator.cpp:120
void AddSource(const ablate::domain::Range &cellRange, Vec localXVec, Vec localFVec) override
Definition: sourceCalculator.cpp:404
SourceCalculator(const std::vector< domain::Field > &fields, std::shared_ptr< TChem > tChemEos, ChemistryConstraints constraints, const ablate::domain::Range &cellRange)
Definition: sourceCalculator.cpp:31
ReactorType
Definition: sourceCalculator.hpp:23
Definition: loggable.hpp:9
Definition: range.hpp:11
hold a struct that can be used for chemistry constraints
Definition: sourceCalculator.hpp:26