1 #ifndef ABLATELIBRARY_SOOTSPECTRUMPROPERTIES_HPP
2 #define ABLATELIBRARY_SOOTSPECTRUMPROPERTIES_HPP
4 #include "eos/tChemSoot.hpp"
5 #include "finiteVolume/compressibleFlowFields.hpp"
6 #include "radiation/radiation.hpp"
7 #include "radiationProperties.hpp"
8 #include "utilities/constants.hpp"
10 namespace ablate::eos::radiationProperties {
13 struct FunctionContext {
14 PetscInt densityYiCSolidCOffset;
17 const std::vector<PetscReal> wavelengths;
18 const std::vector<PetscReal> bandwidths;
20 const std::shared_ptr<eos::EOS> eos;
21 constexpr
static PetscReal rhoC = 2000;
23 std::vector<PetscReal> wavelengthsIn;
24 std::vector<PetscReal> bandwidthsIn;
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 = {});
30 ThermodynamicFunction GetRadiationPropertiesFunction(RadiationProperty property,
const std::vector<domain::Field>& fields)
const;
40 return 1.811 + (0.1263 * log(lambda)) + (0.027 * log(lambda) * log(lambda)) + (0.0417 * log(lambda) * log(lambda) * log(lambda));
42 static inline PetscReal GetAbsorptiveIndex(PetscReal lambda) {
44 return 0.5821 + (0.1213 * log(lambda)) + (0.2309 * log(lambda) * log(lambda)) - (0.01 * log(lambda) * log(lambda) * log(lambda));
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