radiativeIntensityRay Class Reference

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

Collaboration diagram for radiativeIntensityRay:

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 (const radiativeIntensityRay &)=delete
 Disallow default bitwise copy construction. 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 emitted 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 emitted 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...
 
void operator= (const radiativeIntensityRay &)=delete
 Disallow default bitwise assignment. More...
 

Static Public Attributes

static const word intensityPrefix
 

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]

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() [2/2]

Disallow default bitwise copy construction.

◆ ~radiativeIntensityRay()

Destructor.

Definition at line 236 of file radiativeIntensityRay.C.

Member Function Documentation

◆ correct()

Foam::scalar correct ( )

Update radiative intensity on i direction.

Definition at line 242 of file radiativeIntensityRay.C.

References Foam::fvm::div(), forAll, SolverPerformance< Type >::initialResidual(), k, Foam::max(), Foam::constant::mathematical::pi(), fvMatrix< Type >::relax(), 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 284 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.

◆ qr() [1/2]

const 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]

Foam::volScalarField & qr ( )
inline

Return non-const access to the boundary heat flux.

Definition at line 40 of file radiativeIntensityRayI.H.

References radiativeIntensityRay::qin().

Here is the call graph for this function:

◆ qin() [1/2]

Foam::volScalarField & qin ( )
inline

Return non-const access to the boundary incident heat flux.

Definition at line 52 of file radiativeIntensityRayI.H.

References radiativeIntensityRay::qem().

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

Here is the call graph for this function:
Here is the caller graph for this function:

◆ qem() [1/2]

Foam::volScalarField & qem ( )
inline

Return non-const access to the boundary emitted heat flux.

Definition at line 65 of file radiativeIntensityRayI.H.

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

Here is the caller graph for this function:

◆ qin() [2/2]

const Foam::volScalarField & qin ( ) const
inline

Return const access to the boundary incident heat flux.

Definition at line 46 of file radiativeIntensityRayI.H.

◆ qem() [2/2]

const Foam::volScalarField & qem ( ) const
inline

Return const access to the boundary emitted heat flux.

Definition at line 59 of file radiativeIntensityRayI.H.

◆ d()

const Foam::vector & d ( ) const
inline

Return direction.

Definition at line 72 of file radiativeIntensityRayI.H.

Referenced by wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs().

Here is the caller graph for this function:

◆ dAve()

const Foam::vector & dAve ( ) const
inline

Return the average vector inside the solid angle.

Definition at line 79 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 86 of file radiativeIntensityRayI.H.

◆ phi()

Foam::scalar phi ( ) const
inline

Return the phi angle.

Definition at line 92 of file radiativeIntensityRayI.H.

◆ theta()

Foam::scalar theta ( ) const
inline

Return the theta angle.

Definition at line 98 of file radiativeIntensityRayI.H.

◆ omega()

Foam::scalar omega ( ) const
inline

Return the solid angle.

Definition at line 104 of file radiativeIntensityRayI.H.

References radiativeIntensityRay::ILambda().

Here is the call graph for this function:

◆ ILambda()

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

Return the radiative intensity for a given wavelength.

Definition at line 112 of file radiativeIntensityRayI.H.

Referenced by radiativeIntensityRay::omega().

Here is the caller graph for this function:

◆ operator=()

void operator= ( const radiativeIntensityRay )
delete

Disallow default bitwise assignment.

Member Data Documentation

◆ intensityPrefix

const Foam::word intensityPrefix
static

Definition at line 59 of file radiativeIntensityRay.H.


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