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 emmited 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 emmited heat flux. More... | |
const vector & | d () const |
Return direction. More... | |
const vector & | dAve () const |
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 |
Private Member Functions | |
radiativeIntensityRay (const radiativeIntensityRay &) | |
Disallow default bitwise copy construct. More... | |
void | operator= (const radiativeIntensityRay &) |
Disallow default bitwise assignment. More... | |
Private Attributes | |
const fvDOM & | dom_ |
Refence to the owner fvDOM object. More... | |
const fvMesh & | mesh_ |
Reference to the mesh. More... | |
const absorptionEmissionModel & | absorptionEmission_ |
Absorption/emission model. More... | |
const blackBodyEmission & | blackBody_ |
Black body. More... | |
volScalarField | I_ |
Total radiative intensity / [W/m2]. More... | |
volScalarField | Qr_ |
Total radiative heat flux on boundary. More... | |
volScalarField | Qin_ |
Incident radiative heat flux on boundary. More... | |
volScalarField | Qem_ |
Emitted radiative heat flux on boundary. More... | |
vector | d_ |
Direction. More... | |
vector | dAve_ |
Average direction vector inside the solid angle. More... | |
scalar | theta_ |
Theta angle. More... | |
scalar | phi_ |
Phi angle. More... | |
scalar | omega_ |
Solid angle. More... | |
label | nLambda_ |
Number of wavelengths/bands. More... | |
PtrList< volScalarField > | ILambda_ |
List of pointers to radiative intensity fields for given wavelengths. More... | |
label | myRayId_ |
My ray Id. More... | |
Static Private Attributes | |
static label | rayId |
Global ray id - incremented in constructor. More... | |
Radiation intensity for a ray in a given direction.
Definition at line 55 of file radiativeIntensityRay.H.
|
private |
Disallow default bitwise copy construct.
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 40 of file radiativeIntensityRay.C.
References IOobject::AUTO_WRITE, Foam::cos(), Foam::vectorTools::cosPhi(), forAll, IOobject::headerOk(), IOobject::MUST_READ, Foam::name(), VectorSpace< Vector< scalar >, scalar, 3 >::nComponents, IOobject::NO_READ, IOobject::NO_WRITE, phi, IOobject::readOpt(), autoPtr::reset(), Foam::rotationTensor(), Foam::sin(), autoPtr::valid(), and Vector< scalar >::zero.
Destructor.
Definition at line 245 of file radiativeIntensityRay.C.
|
private |
Disallow default bitwise assignment.
Foam::scalar correct | ( | ) |
Update radiative intensity on i direction.
Definition at line 251 of file radiativeIntensityRay.C.
References Foam::fvm::div(), forAll, SolverPerformance::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, and Foam::pow3().
|
inline |
Return intensity.
Definition at line 27 of file radiativeIntensityRayI.H.
References radiativeIntensityRay::I_.
|
inline |
Return const access to the boundary heat flux.
Definition at line 34 of file radiativeIntensityRayI.H.
Referenced by wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs(), and greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs().
|
inline |
Return non-const access to the boundary heat flux.
|
inline |
Return non-const access to the boundary incident heat flux.
Definition at line 46 of file radiativeIntensityRayI.H.
Referenced by wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs(), and greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs().
|
inline |
Return non-const access to the boundary emmited heat flux.
Definition at line 59 of file radiativeIntensityRayI.H.
Referenced by wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs(), and greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs().
|
inline |
Return const access to the boundary incident heat flux.
|
inline |
Return const access to the boundary emmited heat flux.
|
inline |
Return direction.
Definition at line 71 of file radiativeIntensityRayI.H.
|
inline |
Return the average vector inside the solid angle.
Definition at line 77 of file radiativeIntensityRayI.H.
Referenced by wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs(), and greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs().
|
inline |
Return the number of bands.
Definition at line 83 of file radiativeIntensityRayI.H.
|
inline |
Return the phi angle.
Definition at line 89 of file radiativeIntensityRayI.H.
|
inline |
Return the theta angle.
Definition at line 95 of file radiativeIntensityRayI.H.
|
inline |
Return the solid angle.
Definition at line 101 of file radiativeIntensityRayI.H.
|
inline |
Return the radiative intensity for a given wavelength.
Definition at line 109 of file radiativeIntensityRayI.H.
|
static |
Definition at line 59 of file radiativeIntensityRay.H.
|
private |
Refence to the owner fvDOM object.
Definition at line 67 of file radiativeIntensityRay.H.
|
private |
Reference to the mesh.
Definition at line 70 of file radiativeIntensityRay.H.
|
private |
Absorption/emission model.
Definition at line 73 of file radiativeIntensityRay.H.
|
private |
Black body.
Definition at line 76 of file radiativeIntensityRay.H.
|
private |
Total radiative intensity / [W/m2].
Definition at line 79 of file radiativeIntensityRay.H.
Referenced by radiativeIntensityRay::I().
|
private |
Total radiative heat flux on boundary.
Definition at line 82 of file radiativeIntensityRay.H.
|
private |
Incident radiative heat flux on boundary.
Definition at line 85 of file radiativeIntensityRay.H.
|
private |
Emitted radiative heat flux on boundary.
Definition at line 88 of file radiativeIntensityRay.H.
|
private |
Direction.
Definition at line 91 of file radiativeIntensityRay.H.
|
private |
Average direction vector inside the solid angle.
Definition at line 94 of file radiativeIntensityRay.H.
|
private |
Theta angle.
Definition at line 97 of file radiativeIntensityRay.H.
|
private |
Phi angle.
Definition at line 100 of file radiativeIntensityRay.H.
|
private |
Solid angle.
Definition at line 103 of file radiativeIntensityRay.H.
|
private |
Number of wavelengths/bands.
Definition at line 106 of file radiativeIntensityRay.H.
|
private |
List of pointers to radiative intensity fields for given wavelengths.
Definition at line 109 of file radiativeIntensityRay.H.
|
staticprivate |
Global ray id - incremented in constructor.
Definition at line 112 of file radiativeIntensityRay.H.
|
private |
My ray Id.
Definition at line 115 of file radiativeIntensityRay.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.