ABLATE Source Documentation  0.12.34
sootSpectrumProperties.hpp
1 #ifndef ABLATELIBRARY_SOOTSPECTRUMPROPERTIES_HPP
2 #define ABLATELIBRARY_SOOTSPECTRUMPROPERTIES_HPP
3 
4 #include "eos/tChemSoot.hpp"
5 #include "finiteVolume/compressibleFlowFields.hpp"
6 #include "radiation/radiation.hpp"
7 #include "radiationProperties.hpp"
8 #include "utilities/constants.hpp"
9 
10 namespace ablate::eos::radiationProperties {
12  private:
13  struct FunctionContext {
14  PetscInt densityYiCSolidCOffset;
15  const ThermodynamicFunction temperatureFunction;
16  const ThermodynamicTemperatureFunction densityFunction;
17  const std::vector<PetscReal> wavelengths;
18  const std::vector<PetscReal> bandwidths;
19  };
20  const std::shared_ptr<eos::EOS> eos;
21  constexpr static PetscReal rhoC = 2000;
22 
23  std::vector<PetscReal> wavelengthsIn;
24  std::vector<PetscReal> bandwidthsIn;
25 
26  public:
27  SootSpectrumProperties(std::shared_ptr<eos::EOS> eosIn, int num = 0, double min = 0.4E-6, double max = 30E-6, const std::vector<double>& wavelengths = {},
28  const std::vector<double>& bandwidths = {});
29 
30  ThermodynamicFunction GetRadiationPropertiesFunction(RadiationProperty property, const std::vector<domain::Field>& fields) const;
31  ThermodynamicTemperatureFunction GetRadiationPropertiesTemperatureFunction(RadiationProperty property, const std::vector<domain::Field>& fields) const;
32 
33  static PetscErrorCode SootAbsorptionTemperatureFunction(const PetscReal* conserved, PetscReal temperature, PetscReal* kappa, void* ctx);
34 
35  static PetscErrorCode SootEmissionTemperatureFunction(const PetscReal* conserved, PetscReal temperature, PetscReal* epsilon, void* ctx);
36 
38  static inline PetscReal GetRefractiveIndex(PetscReal lambda) {
39  lambda *= 1E6;
40  return 1.811 + (0.1263 * log(lambda)) + (0.027 * log(lambda) * log(lambda)) + (0.0417 * log(lambda) * log(lambda) * log(lambda));
41  }
42  static inline PetscReal GetAbsorptiveIndex(PetscReal lambda) {
43  lambda *= 1E6;
44  return 0.5821 + (0.1213 * log(lambda)) + (0.2309 * log(lambda) * log(lambda)) - (0.01 * log(lambda) * log(lambda) * log(lambda));
45  }
46 };
47 } // namespace ablate::eos::radiationProperties
48 
49 #endif // ABLATELIBRARY_SOOTSPECTRUMPROPERTIES_HPP
Definition: radiationProperties.hpp:11
Definition: sootSpectrumProperties.hpp:11
SootSpectrumProperties(std::shared_ptr< eos::EOS > eosIn, int num=0, double min=0.4E-6, double max=30E-6, const std::vector< double > &wavelengths={}, const std::vector< double > &bandwidths={})
Definition: sootSpectrumProperties.cpp:3
static PetscErrorCode SootAbsorptionTemperatureFunction(const PetscReal *conserved, PetscReal temperature, PetscReal *kappa, void *ctx)
Definition: sootSpectrumProperties.cpp:43
static PetscErrorCode SootEmissionTemperatureFunction(const PetscReal *conserved, PetscReal temperature, PetscReal *epsilon, void *ctx)
Definition: sootSpectrumProperties.cpp:27
static PetscReal GetRefractiveIndex(PetscReal lambda)
Polynomial fits to soot data for hydrocarbon combustion conditions (Modest, ch. 11 pg....
Definition: sootSpectrumProperties.hpp:38
ThermodynamicTemperatureFunction GetRadiationPropertiesTemperatureFunction(RadiationProperty property, const std::vector< domain::Field > &fields) const
Definition: sootSpectrumProperties.cpp:69