Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | Static Private Attributes
radiativeIntensityRay Class Reference

Radiation intensity for a ray in a given direction. More...

Collaboration diagram for radiativeIntensityRay:
Collaboration graph
[legend]

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 volScalarFieldI () const
 Return intensity. More...
 
const volScalarFieldQr () const
 Return const access to the boundary heat flux. More...
 
volScalarFieldQr ()
 Return non-const access to the boundary heat flux. More...
 
volScalarFieldQin ()
 Return non-const access to the boundary incident heat flux. More...
 
volScalarFieldQem ()
 Return non-const access to the boundary emmited heat flux. More...
 
const volScalarFieldQin () const
 Return const access to the boundary incident heat flux. More...
 
const volScalarFieldQem () const
 Return const access to the boundary emmited heat flux. More...
 
const vectord () const
 Return direction. More...
 
const vectordAve () 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 volScalarFieldILambda (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 fvDOMdom_
 Refence to the owner fvDOM object. More...
 
const fvMeshmesh_
 Reference to the mesh. More...
 
const absorptionEmissionModelabsorptionEmission_
 Absorption/emission model. More...
 
const blackBodyEmissionblackBody_
 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< volScalarFieldILambda_
 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...
 

Detailed Description

Radiation intensity for a ray in a given direction.

Source files

Definition at line 55 of file radiativeIntensityRay.H.

Constructor & Destructor Documentation

◆ radiativeIntensityRay() [1/2]

Disallow default bitwise copy construct.

◆ radiativeIntensityRay() [2/2]

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 
)

◆ ~radiativeIntensityRay()

Destructor.

Definition at line 245 of file radiativeIntensityRay.C.

Member Function Documentation

◆ operator=()

void operator= ( const radiativeIntensityRay )
private

Disallow default bitwise assignment.

◆ correct()

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().

Here is the call graph for this function:

◆ init()

void init ( const scalar  phi,
const scalar  theta,
const scalar  deltaPhi,
const scalar  deltaTheta,
const scalar  lambda 
)

Initialise the ray in i direction.

◆ addIntensity()

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().

Here is the call graph for this function:

◆ I()

const Foam::volScalarField & I ( ) const
inline

Return intensity.

Definition at line 27 of file radiativeIntensityRayI.H.

References radiativeIntensityRay::I_.

◆ Qr() [1/2]

Foam::volScalarField & Qr ( ) const
inline

Return const access to the boundary heat flux.

Definition at line 34 of file radiativeIntensityRayI.H.

Referenced by wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs(), and greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs().

Here is the caller graph for this function:

◆ Qr() [2/2]

volScalarField& Qr ( )
inline

Return non-const access to the boundary heat flux.

◆ Qin() [1/2]

Foam::volScalarField & Qin ( )
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().

Here is the caller graph for this function:

◆ Qem() [1/2]

Foam::volScalarField & Qem ( )
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().

Here is the caller graph for this function:

◆ Qin() [2/2]

const volScalarField& Qin ( ) const
inline

Return const access to the boundary incident heat flux.

◆ Qem() [2/2]

const volScalarField& Qem ( ) const
inline

Return const access to the boundary emmited heat flux.

◆ d()

const Foam::vector & d ( ) const
inline

Return direction.

Definition at line 71 of file radiativeIntensityRayI.H.

◆ dAve()

const Foam::vector & dAve ( ) const
inline

Return the average vector inside the solid angle.

Definition at line 77 of file radiativeIntensityRayI.H.

Referenced by wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs(), and greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs().

Here is the caller graph for this function:

◆ nLambda()

Foam::scalar nLambda ( ) const
inline

Return the number of bands.

Definition at line 83 of file radiativeIntensityRayI.H.

◆ phi()

Foam::scalar phi ( ) const
inline

Return the phi angle.

Definition at line 89 of file radiativeIntensityRayI.H.

◆ theta()

Foam::scalar theta ( ) const
inline

Return the theta angle.

Definition at line 95 of file radiativeIntensityRayI.H.

◆ omega()

Foam::scalar omega ( ) const
inline

Return the solid angle.

Definition at line 101 of file radiativeIntensityRayI.H.

◆ ILambda()

const Foam::volScalarField & ILambda ( const label  lambdaI) const
inline

Return the radiative intensity for a given wavelength.

Definition at line 109 of file radiativeIntensityRayI.H.

Field Documentation

◆ intensityPrefix

const Foam::word intensityPrefix
static

Definition at line 59 of file radiativeIntensityRay.H.

◆ dom_

const fvDOM& dom_
private

Refence to the owner fvDOM object.

Definition at line 67 of file radiativeIntensityRay.H.

◆ mesh_

const fvMesh& mesh_
private

Reference to the mesh.

Definition at line 70 of file radiativeIntensityRay.H.

◆ absorptionEmission_

const absorptionEmissionModel& absorptionEmission_
private

Absorption/emission model.

Definition at line 73 of file radiativeIntensityRay.H.

◆ blackBody_

const blackBodyEmission& blackBody_
private

Black body.

Definition at line 76 of file radiativeIntensityRay.H.

◆ I_

volScalarField I_
private

Total radiative intensity / [W/m2].

Definition at line 79 of file radiativeIntensityRay.H.

Referenced by radiativeIntensityRay::I().

◆ Qr_

volScalarField Qr_
private

Total radiative heat flux on boundary.

Definition at line 82 of file radiativeIntensityRay.H.

◆ Qin_

volScalarField Qin_
private

Incident radiative heat flux on boundary.

Definition at line 85 of file radiativeIntensityRay.H.

◆ Qem_

volScalarField Qem_
private

Emitted radiative heat flux on boundary.

Definition at line 88 of file radiativeIntensityRay.H.

◆ d_

vector d_
private

Direction.

Definition at line 91 of file radiativeIntensityRay.H.

◆ dAve_

vector dAve_
private

Average direction vector inside the solid angle.

Definition at line 94 of file radiativeIntensityRay.H.

◆ theta_

scalar theta_
private

Theta angle.

Definition at line 97 of file radiativeIntensityRay.H.

◆ phi_

scalar phi_
private

Phi angle.

Definition at line 100 of file radiativeIntensityRay.H.

◆ omega_

scalar omega_
private

Solid angle.

Definition at line 103 of file radiativeIntensityRay.H.

◆ nLambda_

label nLambda_
private

Number of wavelengths/bands.

Definition at line 106 of file radiativeIntensityRay.H.

◆ ILambda_

PtrList<volScalarField> ILambda_
private

List of pointers to radiative intensity fields for given wavelengths.

Definition at line 109 of file radiativeIntensityRay.H.

◆ rayId

label rayId
staticprivate

Global ray id - incremented in constructor.

Definition at line 112 of file radiativeIntensityRay.H.

◆ myRayId_

label myRayId_
private

My ray Id.

Definition at line 115 of file radiativeIntensityRay.H.


The documentation for this class was generated from the following files: