ABLATE Source Documentation  0.12.34
sourceCalculatorZeroRK.hpp
1 #ifndef ABLATELIBRARY_ZERORK_SOURCECALCULATOR_HPP
2 #define ABLATELIBRARY_ZERORK_SOURCECALCULATOR_HPP
3 
4 #include "eos/chemistryModel.hpp"
5 #include "zerork/mechanism.h"
6 #include "zerork/utilities.h"
7 #include "zerork_cfd_plugin.h"
8 
9 namespace ablate::eos {
10 class zerorkEOS;
11 }
12 
13 namespace ablate::eos::zerorkeos {
14 
18 class SourceCalculator : public ChemistryModel::SourceCalculator, private utilities::Loggable<SourceCalculator> {
19  public:
23  enum class ReactorType { ConstantVolume, ConstantPressure };
24 
27  double relTolerance = 1.0E-6;
28  double absTolerance = 1.0E-10;
29 
30  // set the limiter on chemical rates of progress,
31  // for 1D flames set it to a value btw 1e18 and 1e20,
32  // 1e16 starts changing the flameshape
33  double stepLimiter = 1.0E200;
34 
35  // use cvode or SEULEX
36  int useSEULEX = 0;
37 
38  // zerork output state
39  int verbose = 0;
40 
41  // load balancing
42  int loadBalance = 1;
43 
44  // load balancing
45  int gpu = 0;
46 
47  // max iterations for cvode
48  int dumpfailed = 0;
49 
50  // max iterations for cvode
51  int maxiteration = 10000;
52 
53  // max iterations for cvode
54  int cvode_num_retries = 5;
55 
56  // max iterations for cvode
57  double cvode_retry_absolute_tolerance_adjustment = 0.1;
58 
59  // max iterations for cvode
60  double cvode_retry_relative_tolerance_adjustment = 0.1;
61 
62  // iterative solution for approximate jacobian recommended to be set to 0 whenever dense is selected
63  int iterative = 0;
64 
65  // For large mechanisms(~Nspec>100) it is recommended to use a sparse math jacobian option
66  bool sparseJacobian = false;
67 
68  // timing log information
69  bool timinglog = false;
70 
71  // store the reactor type in the chemistry constrains
72  ReactorType reactorType = ReactorType::ConstantVolume;
73 
74  // store an optional threshold temperature. Only compute the reactions if the temperature is above thresholdTemperature
75  double thresholdTemperature = 0.0;
76 
77  int errorhandle = 0;
78 
79  bool dumpreactor = false;
80 
81  void Set(const std::shared_ptr<ablate::parameters::Parameters>&);
82  };
89  SourceCalculator(const std::vector<domain::Field>& fields, std::shared_ptr<zerorkEOS> zerorkEos, ablate::eos::zerorkeos::SourceCalculator::ChemistryConstraints constraints,
90  const ablate::domain::Range& cellRange);
91 
95  void ComputeSource(const ablate::domain::Range& cellRange, PetscReal time, PetscReal dt, Vec globalSolution) override;
96 
100  void AddSource(const ablate::domain::Range& cellRange, Vec localXVec, Vec localFVec) override;
101 
102  private:
103  std::vector<double> sourceZeroRKAtI;
104  zerork_handle zrm_handle;
110  std::shared_ptr<eos::zerorkEOS> eos;
111 
112  const size_t numberSpecies;
113 
115  PetscInt eulerId;
116 
118  PetscInt densityYiId;
119 };
120 
127 std::ostream& operator<<(std::ostream& os, const SourceCalculator::ReactorType& v);
128 
135 std::istream& operator>>(std::istream& is, SourceCalculator::ReactorType& v);
136 
137 } // namespace ablate::eos::zerorkeos
138 
139 #endif // ABLATELIBRARY_BATCHSOURCE_HPP
Definition: chemistryModel.hpp:30
Definition: sourceCalculatorZeroRK.hpp:18
ReactorType
Definition: sourceCalculatorZeroRK.hpp:23
void AddSource(const ablate::domain::Range &cellRange, Vec localXVec, Vec localFVec) override
Definition: sourceCalculatorZeroRK.cpp:332
void ComputeSource(const ablate::domain::Range &cellRange, PetscReal time, PetscReal dt, Vec globalSolution) override
Definition: sourceCalculatorZeroRK.cpp:141
SourceCalculator(const std::vector< domain::Field > &fields, std::shared_ptr< zerorkEOS > zerorkEos, ablate::eos::zerorkeos::SourceCalculator::ChemistryConstraints constraints, const ablate::domain::Range &cellRange)
Definition: sourceCalculatorZeroRK.cpp:33
Definition: loggable.hpp:9
Definition: range.hpp:11
hold a struct that can be used for chemistry constraints
Definition: sourceCalculatorZeroRK.hpp:26