ABLATE Source Documentation  0.12.34
zimmer.hpp
1 #ifndef ABLATELIBRARY_RADIATIONPROPERTIESZIMMER_H
2 #define ABLATELIBRARY_RADIATIONPROPERTIESZIMMER_H
3 
4 #include <array>
5 #include "eos/radiationProperties/radiationProperties.hpp"
6 #include "finiteVolume/compressibleFlowFields.hpp"
7 #include "radiation/radiation.hpp"
8 #include "radiationProperties.hpp"
9 #include "solver/cellSolver.hpp"
10 #include "utilities/mathUtilities.hpp"
11 
12 namespace ablate::eos::radiationProperties {
14 class Zimmer : public RadiationModel {
15  private:
16  struct FunctionContext {
17  PetscInt densityYiH2OOffset;
18  PetscInt densityYiCO2Offset;
19  PetscInt densityYiCOOffset;
20  PetscInt densityYiCH4Offset;
21 
22  PetscReal upperLimit;
23  PetscReal lowerLimit;
24 
25  const ThermodynamicFunction temperatureFunction;
26  const ThermodynamicTemperatureFunction densityFunction;
27  };
28 
32  const std::shared_ptr<eos::EOS> eos;
33 
35  constexpr static double kapparef = 1;
36  constexpr static double Tsurf = 300.;
37 
39  constexpr static std::array<double, 7> H2O_coeff = {0.22317e1, -.15829e1, .1329601e1, -.50707, .93334e-1, -0.83108e-2, 0.28834e-3};
40  constexpr static std::array<double, 7> CO2_coeff = {0.38041E1, -0.27808E1, 0.11672E1, -0.284910E0, 0.38163E-1, -0.26292E-2, 0.73662E-4};
41  constexpr static std::array<double, 5> CH4_coeff = {6.6334E0, -3.5686E-3, 1.6682E-8, 2.5611E-10, -2.6558E-14};
42  constexpr static std::array<double, 5> CO_1_coeff = {4.7869E0, -6.953E-2, 2.95775E-4, -4.25732E-7, 2.02894E-10};
43  constexpr static std::array<double, 5> CO_2_coeff = {1.0091E1, -1.183E-2, 4.7753E-6, -5.87209E-10, -2.5334E-14};
44  constexpr static double MWC = 1.2010700e+01;
45  constexpr static double MWO = 1.599940e+01;
46  constexpr static double MWH = 1.007940e+00;
47  constexpr static double UGC = 8314.4;
48  constexpr static double MWCO = MWC + MWO;
49  constexpr static double MWCO2 = MWC + 2. * MWO;
50  constexpr static double MWCH4 = MWC + 4. * MWH;
51  constexpr static double MWH2O = 2. * MWH + MWO;
52 
61  static PetscErrorCode ZimmerEmissionTemperatureFunction(const PetscReal conserved[], PetscReal temperature, PetscReal* epsilon, void* ctx);
62 
69  static PetscErrorCode ZimmerAbsorptionTemperatureFunction(const PetscReal conserved[], PetscReal temperature, PetscReal* kappa, void* ctx);
70 
71  public:
72  explicit Zimmer(std::shared_ptr<eos::EOS> eosIn, PetscReal upperLimitIn = 0, PetscReal lowerLimitIn = 0);
73  explicit Zimmer(const Zimmer&) = delete;
74  void operator=(const Zimmer&) = delete;
75 
82  [[nodiscard]] ThermodynamicTemperatureFunction GetRadiationPropertiesTemperatureFunction(RadiationProperty property, const std::vector<domain::Field>& fields) const override;
83 
84  PetscInt GetFieldComponentOffset(const std::string& str, const domain::Field& field) const;
85 
86  PetscReal upperLimitStored;
87  PetscReal lowerLimitStored;
88 };
89 
90 } // namespace ablate::eos::radiationProperties
91 
92 #endif // ABLATELIBRARY_RADIATIONPROPERTIESZIMMER_H
Definition: radiationProperties.hpp:11
ThermodynamicTemperatureFunction GetRadiationPropertiesTemperatureFunction(RadiationProperty property, const std::vector< domain::Field > &fields) const override
Definition: zimmer.cpp:83
PetscInt GetFieldComponentOffset(const std::string &str, const domain::Field &field) const
Definition: zimmer.cpp:136
Definition: field.hpp:17