Radiation intensity for a ray in a given direction. More...
Public Member Functions | |
radiativeIntensityRay (const fvDOM &dom, const fvMesh &mesh, const scalar phi, const scalar theta, const scalar deltaPhi, const scalar deltaTheta, const label lambda, const absorptionEmissionModel &absEmmModel_, const blackBodyEmission &blackBody, const label rayId) | |
Construct form components. More... | |
~radiativeIntensityRay () | |
Destructor. More... | |
scalar | correct () |
Update radiative intensity on i direction. More... | |
void | init (const scalar phi, const scalar theta, const scalar deltaPhi, const scalar deltaTheta, const scalar lambda) |
Initialise the ray in i direction. More... | |
void | addIntensity () |
Add radiative intensities from all the bands. More... | |
const volScalarField & | I () const |
Return intensity. More... | |
const volScalarField & | qr () const |
Return const access to the boundary heat flux. More... | |
volScalarField & | qr () |
Return non-const access to the boundary heat flux. More... | |
volScalarField & | qin () |
Return non-const access to the boundary incident heat flux. More... | |
volScalarField & | qem () |
Return non-const access to the boundary emitted heat flux. More... | |
const volScalarField & | qin () const |
Return const access to the boundary incident heat flux. More... | |
const volScalarField & | qem () const |
Return const access to the boundary emitted heat flux. More... | |
const vector & | d () const |
Return direction. More... | |
const vector & | dAve () const |
Return the average vector inside the solid angle. More... | |
vector & | d () |
Return direction. More... | |
vector & | dAve () |
Return the average vector inside the solid angle. More... | |
scalar | nLambda () const |
Return the number of bands. More... | |
scalar | phi () const |
Return the phi angle. More... | |
scalar | theta () const |
Return the theta angle. More... | |
scalar | omega () const |
Return the solid angle. More... | |
const volScalarField & | ILambda (const label lambdaI) const |
Return the radiative intensity for a given wavelength. More... | |
Static Public Attributes | |
static const word | intensityPrefix |
Radiation intensity for a ray in a given direction.
Definition at line 57 of file radiativeIntensityRay.H.
radiativeIntensityRay | ( | const fvDOM & | dom, |
const fvMesh & | mesh, | ||
const scalar | phi, | ||
const scalar | theta, | ||
const scalar | deltaPhi, | ||
const scalar | deltaTheta, | ||
const label | lambda, | ||
const absorptionEmissionModel & | absEmmModel_, | ||
const blackBodyEmission & | blackBody, | ||
const label | rayId | ||
) |
Construct form components.
Definition at line 43 of file radiativeIntensityRay.C.
References IOobject::AUTO_WRITE, Foam::cos(), Foam::vectorTools::cosPhi(), forAll, IOobject::MUST_READ, Foam::name(), VectorSpace< Vector< scalar >, scalar, 3 >::nComponents, IOobject::NO_READ, IOobject::NO_WRITE, phi, IOobject::readOpt(), autoPtr< T >::reset(), Foam::rotationTensor(), Foam::sin(), VectorSpace< Vector< scalar >, scalar, 3 >::zero, and Foam::Zero.
Destructor.
Definition at line 263 of file radiativeIntensityRay.C.
Foam::scalar correct | ( | ) |
Update radiative intensity on i direction.
Definition at line 269 of file radiativeIntensityRay.C.
References Foam::fvm::div(), forAll, SolverPerformance< Type >::initialResidual(), k, Foam::max(), Foam::constant::mathematical::pi(), Foam::solve(), and Foam::fvm::Sp().
void init | ( | const scalar | phi, |
const scalar | theta, | ||
const scalar | deltaPhi, | ||
const scalar | deltaTheta, | ||
const scalar | lambda | ||
) |
Initialise the ray in i direction.
void addIntensity | ( | ) |
Add radiative intensities from all the bands.
Definition at line 316 of file radiativeIntensityRay.C.
References Foam::dimMass, Foam::dimTime, forAll, Foam::pow3(), and Foam::Zero.
|
inline |
Return intensity.
Definition at line 29 of file radiativeIntensityRayI.H.
|
inline |
Return const access to the boundary heat flux.
Definition at line 36 of file radiativeIntensityRayI.H.
Referenced by wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs(), and greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs().
|
inline |
Return non-const access to the boundary heat flux.
Definition at line 42 of file radiativeIntensityRayI.H.
|
inline |
Return non-const access to the boundary incident heat flux.
Definition at line 54 of file radiativeIntensityRayI.H.
Referenced by wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs(), and greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs().
|
inline |
Return non-const access to the boundary emitted heat flux.
Definition at line 67 of file radiativeIntensityRayI.H.
Referenced by wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs(), and greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs().
|
inline |
Return const access to the boundary incident heat flux.
Definition at line 48 of file radiativeIntensityRayI.H.
|
inline |
Return const access to the boundary emitted heat flux.
Definition at line 61 of file radiativeIntensityRayI.H.
|
inline |
Return direction.
Definition at line 73 of file radiativeIntensityRayI.H.
|
inline |
Return the average vector inside the solid angle.
Definition at line 79 of file radiativeIntensityRayI.H.
Referenced by wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs(), and greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs().
|
inline |
Return direction.
Definition at line 85 of file radiativeIntensityRayI.H.
|
inline |
Return the average vector inside the solid angle.
Definition at line 91 of file radiativeIntensityRayI.H.
|
inline |
Return the number of bands.
Definition at line 97 of file radiativeIntensityRayI.H.
|
inline |
Return the phi angle.
Definition at line 103 of file radiativeIntensityRayI.H.
|
inline |
Return the theta angle.
Definition at line 109 of file radiativeIntensityRayI.H.
|
inline |
Return the solid angle.
Definition at line 115 of file radiativeIntensityRayI.H.
|
inline |
Return the radiative intensity for a given wavelength.
Definition at line 123 of file radiativeIntensityRayI.H.
|
static |
Definition at line 61 of file radiativeIntensityRay.H.