ABLATE Source Documentation  0.12.34
surfaceRadiation.hpp
1 #ifndef ABLATELIBRARY_SURFACERADIATION_HPP
2 #define ABLATELIBRARY_SURFACERADIATION_HPP
3 
4 #include "domain/range.hpp"
5 #include "domain/reverseRange.hpp"
6 #include "radiation.hpp"
7 #include "utilities/constants.hpp"
8 
9 namespace ablate::radiation {
10 
12  protected:
15 
16  public:
17  SurfaceRadiation(const std::string& solverId, const std::shared_ptr<domain::Region>& region, const PetscInt raynumber, std::shared_ptr<eos::radiationProperties::RadiationModel> radiationModelIn,
18  std::shared_ptr<ablate::monitors::logs::Log> = {});
20 
21  void Initialize(const ablate::domain::Range& cellRange, ablate::domain::SubDomain& subDomain) override;
30  PetscReal SurfaceComponent(const PetscReal normal[], PetscInt iCell, PetscInt nphi, PetscInt ntheta) override;
31 
36  static inline std::string GetClassType() { return "SurfaceRadiation"; }
37 
45  virtual inline void GetSurfaceIntensity(PetscReal* intensity, PetscInt faceId, PetscReal temperature, PetscReal emissivity = 1.0, PetscReal absorptivity = 1.0) {
46  // add in precomputed gains
47  for (int i = 0; i < (int)absorptivityFunction.propertySize; ++i) { // Compute the losses
48  // Emitted From Surface = \epsilon \sigma T^4
49  PetscReal netIntensity = -emissivity * ablate::utilities::Constants::sbc * temperature * temperature * temperature * temperature;
50 
51  // Absorbed at the surface = q''_Rad * \alpha
52  netIntensity += absorptivity * evaluatedGains[absorptivityFunction.propertySize * indexLookup.GetAbsoluteIndex(faceId) + i];
53 
54  // This line just ensures that the netIntenisty isn't going to extreme values, If this clipping is really needed, it probably should be fixed somewhere else -klb
55  intensity[i] = abs(netIntensity) > ablate::utilities::Constants::large ? ablate::utilities::Constants::large * PetscSignReal(netIntensity) : netIntensity;
56  }
57  }
58 };
59 } // namespace ablate::radiation
60 #endif // ABLATELIBRARY_SURFACERADIATION_HPP
Definition: subDomain.hpp:19
Definition: radiation.hpp:25
std::vector< PetscScalar > evaluatedGains
size up the evaluated gains, this index is based upon order of the requested cells
Definition: radiation.hpp:252
eos::ThermodynamicTemperatureFunction absorptivityFunction
hold a pointer to the absorptivity function
Definition: radiation.hpp:267
std::string solverId
the name of this solver
Definition: radiation.hpp:258
const std::shared_ptr< domain::Region > region
the region for which this solver applies
Definition: radiation.hpp:261
Definition: surfaceRadiation.hpp:11
static std::string GetClassType()
Definition: surfaceRadiation.hpp:36
ablate::domain::ReverseRange indexLookup
used to look up from the face id to range index
Definition: surfaceRadiation.hpp:14
void Initialize(const ablate::domain::Range &cellRange, ablate::domain::SubDomain &subDomain) override
Definition: surfaceRadiation.cpp:9
virtual void GetSurfaceIntensity(PetscReal *intensity, PetscInt faceId, PetscReal temperature, PetscReal emissivity=1.0, PetscReal absorptivity=1.0)
Definition: surfaceRadiation.hpp:45
PetscReal SurfaceComponent(const PetscReal normal[], PetscInt iCell, PetscInt nphi, PetscInt ntheta) override
Definition: surfaceRadiation.cpp:19
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
Definition: range.hpp:11
Definition: reverseRange.hpp:11
PetscInt GetAbsoluteIndex(PetscInt point) const
Definition: reverseRange.hpp:52
PetscInt propertySize
the property size being set
Definition: eos.hpp:48