ABLATE Source Documentation  0.12.34
sootMeanProperties.hpp
1 #ifndef ABLATELIBRARY_SOOTMEANPROPERTIES_HPP
2 #define ABLATELIBRARY_SOOTMEANPROPERTIES_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 {
13  private:
14  struct FunctionContext {
15  PetscInt densityYiCSolidCOffset;
16  const ThermodynamicFunction temperatureFunction;
17  const ThermodynamicTemperatureFunction densityFunction;
18  };
19  const std::shared_ptr<eos::EOS> eos;
20  constexpr static PetscReal C_2 = (utilities::Constants::h * utilities::Constants::c) / (utilities::Constants::k);
21  constexpr static PetscReal C_0 = 7.0;
22  constexpr static PetscReal rhoC = 2000; // kg/m^3
23  public:
24  SootMeanProperties(std::shared_ptr<EOS> eosIn);
25  ThermodynamicFunction GetRadiationPropertiesFunction(RadiationProperty property, const std::vector<domain::Field>& fields) const;
26  ThermodynamicTemperatureFunction GetRadiationPropertiesTemperatureFunction(RadiationProperty property, const std::vector<domain::Field>& fields) const;
27 
28  static PetscErrorCode SootAbsorptionTemperatureFunction(const PetscReal* conserved, PetscReal temperature, PetscReal* kappa, void* ctx);
29 
30  static PetscErrorCode SootEmissionTemperatureFunction(const PetscReal* conserved, PetscReal temperature, PetscReal* epsilon, void* ctx);
31 
32  static inline PetscReal GetRefractiveIndex() { return 1; }
33 };
34 } // namespace ablate::eos::radiationProperties
35 #endif // ABLATELIBRARY_SOOTMEANPROPERTIES_HPP
Definition: radiationProperties.hpp:11
Definition: sootMeanProperties.hpp:12
static PetscErrorCode SootAbsorptionTemperatureFunction(const PetscReal *conserved, PetscReal temperature, PetscReal *kappa, void *ctx)
Definition: sootMeanProperties.cpp:14
ThermodynamicTemperatureFunction GetRadiationPropertiesTemperatureFunction(RadiationProperty property, const std::vector< domain::Field > &fields) const
Definition: sootMeanProperties.cpp:31
constexpr static PetscReal k
Boltzmann Constant.
Definition: constants.hpp:19
constexpr static PetscReal c
Speed of light.
Definition: constants.hpp:16
constexpr static PetscReal h
Planck Constant.
Definition: constants.hpp:13