ABLATE Source Documentation  0.12.34
orthogonalRadiation.hpp
1 #ifndef ABLATELIBRARY_ORTHOGONALRADIATION_HPP
2 #define ABLATELIBRARY_ORTHOGONALRADIATION_HPP
3 
4 #include "domain/reverseRange.hpp"
5 #include "radiation.hpp"
6 #include "surfaceRadiation.hpp"
7 
8 namespace ablate::radiation {
9 
11  public:
12  OrthogonalRadiation(const std::string& solverId, const std::shared_ptr<domain::Region>& region, std::shared_ptr<eos::radiationProperties::RadiationModel> radiationModelIn,
13  std::shared_ptr<ablate::monitors::logs::Log> = {});
15 
16  void Setup(const ablate::domain::Range& cellRange, ablate::domain::SubDomain& subDomain) override;
17 
22  static inline std::string GetClassType() { return "OrthogonalRadiation"; }
23 
24  inline void GetSurfaceIntensity(PetscReal* intensityReturn, PetscInt faceId, PetscReal temperature, PetscReal emissivity = 1.0, PetscReal absorptivity = 1.0) override {
25  for (int i = 0; i < (int)absorptivityFunction.propertySize; ++i) { // Compute the losses
26  PetscReal netIntensity = evaluatedGains[absorptivityFunction.propertySize * indexLookup.GetAbsoluteIndex(faceId) + i];
27  intensityReturn[i] = abs(netIntensity) > ablate::utilities::Constants::large ? ablate::utilities::Constants::large * PetscSignReal(netIntensity) : netIntensity;
28  }
29  }
30 };
31 } // namespace ablate::radiation
32 #endif // ABLATELIBRARY_ORTHOGONALRADIATION_HPP
Definition: subDomain.hpp:19
Definition: orthogonalRadiation.hpp:10
void GetSurfaceIntensity(PetscReal *intensityReturn, PetscInt faceId, PetscReal temperature, PetscReal emissivity=1.0, PetscReal absorptivity=1.0) override
Definition: orthogonalRadiation.hpp:24
void Setup(const ablate::domain::Range &cellRange, ablate::domain::SubDomain &subDomain) override
Definition: orthogonalRadiation.cpp:9
~OrthogonalRadiation()
The ray number should never be used because there is only one ray emanating from every boundary face.
Definition: orthogonalRadiation.cpp:7
static std::string GetClassType()
Definition: orthogonalRadiation.hpp:22
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
ablate::domain::ReverseRange indexLookup
used to look up from the face id to range index
Definition: surfaceRadiation.hpp:14
constexpr static PetscReal large
A somewhat large number.
Definition: constants.hpp:31
Definition: range.hpp:11
PetscInt GetAbsoluteIndex(PetscInt point) const
Definition: reverseRange.hpp:52
PetscInt propertySize
the property size being set
Definition: eos.hpp:48