ABLATE Source Documentation  0.12.34
radiationLoss.hpp
1 #ifndef ABLATELIBRARY_RADIATIONLOSS_HPP
2 #define ABLATELIBRARY_RADIATIONLOSS_HPP
3 
4 #include "eos/radiationProperties/radiationProperties.hpp"
5 #include "finiteVolume/compressibleFlowFields.hpp"
6 #include "process.hpp"
7 #include "utilities/constants.hpp"
8 
9 namespace ablate::finiteVolume::processes {
13 class RadiationLoss : public Process {
15  const std::map<std::string, std::shared_ptr<ablate::mathFunctions::MathFunction>> functions;
16 
21  static PetscErrorCode ComputeRadiationLoss(PetscInt dim, PetscReal time, const PetscFVCellGeom* cg, const PetscInt uOff[], const PetscScalar u[], const PetscInt aOff[], const PetscScalar a[],
22  PetscScalar f[], void* ctx);
23 
24  public:
25  explicit RadiationLoss(std::shared_ptr<eos::radiationProperties::RadiationModel> radiationModelIn, double tInfinityIn = 300);
26 
31  void Setup(ablate::finiteVolume::FiniteVolumeSolver& fvmSolver) override;
32 
33  inline std::shared_ptr<eos::radiationProperties::RadiationModel> GetRadiationModel() { return radiationModel; }
34 
41  static inline PetscReal GetIntensity(PetscReal tInfinity, PetscReal temperature, PetscReal kappa) {
42  tInfinity = PetscMax(tInfinity, (temperature - 500));
43 
44  // Compute the losses
45  PetscReal netIntensity = -4.0 * ablate::utilities::Constants::sbc * temperature * temperature * temperature * temperature;
46 
47  netIntensity += 4.0 * ablate::utilities::Constants::sbc * tInfinity * tInfinity * tInfinity * tInfinity;
48 
49  // scale by kappa
50  netIntensity *= kappa;
51 
52  return PetscAbsReal(netIntensity) > ablate::utilities::Constants::large ? ablate::utilities::Constants::large * PetscSignReal(netIntensity) : netIntensity;
53  }
54 
56  const std::shared_ptr<eos::radiationProperties::RadiationModel> radiationModel;
57 
60 
61  PetscReal tInfinity;
62 };
63 
64 } // namespace ablate::finiteVolume::processes
65 
66 #endif // ABLATELIBRARY_RADIATIONLOSS_HPP
Definition: finiteVolumeSolver.hpp:28
Definition: radiationLoss.hpp:13
eos::ThermodynamicTemperatureFunction absorptivityFunction
hold a pointer to the absorptivity function
Definition: radiationLoss.hpp:59
void Setup(ablate::finiteVolume::FiniteVolumeSolver &fvmSolver) override
Definition: radiationLoss.cpp:6
static PetscReal GetIntensity(PetscReal tInfinity, PetscReal temperature, PetscReal kappa)
Definition: radiationLoss.hpp:41
const std::shared_ptr< eos::radiationProperties::RadiationModel > radiationModel
model used to provided the absorptivity function
Definition: radiationLoss.hpp:56
constexpr static PetscReal sbc
Stefan-Boltzman Constant (J/K)
Definition: constants.hpp:10
constexpr static PetscReal large
A somewhat large number.
Definition: constants.hpp:31