OpenFOAM logo
The Open Source CFD Toolbox
  Source Guide OpenCFD Solutions Contact OpenFOAM

fvDOM Class Reference

Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating media, not including scatter. More...

Inheritance diagram for fvDOM:
Collaboration diagram for fvDOM:

List of all members.


Public Member Functions

 TypeName ("fvDOM")
 Runtime type information.
 fvDOM (const volScalarField &T)
 Construct from components.
virtual ~fvDOM ()
 Destructor.
void calculate ()
 Solve radiation equation(s).
bool read ()
 Read radiation properties dictionary.
void updateG ()
 Update G and calculate total heat flux on boundary.
void setRayIdLambdaId (const word &name, label &rayId, label &lambdaId) const
 Set the rayId and lambdaId from by decomposing an intensity.
virtual tmp< volScalarFieldRp () const
 Source term component (for power of T^4).
virtual tmp< DimensionedField
< scalar, volMesh > > 
Ru () const
 Source term component (constant).
const radiativeIntensityRayIRay (const label rayI) const
 Ray intensity for rayI.
const volScalarFieldIRayLambda (const label rayI, const label lambdaI) const
 Ray intensity for rayI and lambda bandwidth.
label nTheta () const
 Number of angles in theta.
label nPhi () const
 Number of angles in phi.
label nRay () const
 Number of rays.
label nLambda () const
 Number of wavelengths.
const volScalarFielda () const
 Const access to total absorption coefficient.
const volScalarFieldaLambda (const label lambdaI) const
 Const access to wavelength total absorption coefficient.
const volScalarFieldG () const
 Const access to incident radiation field.
const volScalarFieldQr () const
 Const access to total radiative heat flux field.
const blackBodyEmissionblackBody () const
 Const access to black body.

Detailed Description

Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating media, not including scatter.

Available absorption models: greyMeanAbsoprtionEmission wideBandAbsorptionEmission

i.e. dictionary fvDOMCoeffs { nPhi 1; // azimuthal angles in PI/2 on X-Y.(from Y to X) nTheta 2; // polar angles in PI (from Z to X-Y plane) convergence 1e-4; // convergence criteria for radiation iteration }

solverFreq 1; // Number of flow iterations per radiation iteration

The total number of solid angles is 4*nPhi*nTheta.

In 1D the direction of the rays is X (nPhi and nTheta are ignored) In 2D the direction of the rays is on X-Y plane (only nPhi is considered) In 3D (nPhi and nTheta are considered)

Source files

Definition at line 69 of file fvDOM.H.


Constructor & Destructor Documentation

fvDOM ( const volScalarField T  ) 

Construct from components.

Definition at line 47 of file fvDOM.C.

~fvDOM (  )  [virtual]

Destructor.

Definition at line 258 of file fvDOM.C.


Member Function Documentation

TypeName ( "fvDOM"   ) 

Runtime type information.

void calculate (  )  [virtual]

Solve radiation equation(s).

Implements radiationModel.

Definition at line 282 of file fvDOM.C.

bool read (  )  [virtual]

Read radiation properties dictionary.

Implements radiationModel.

Definition at line 264 of file fvDOM.C.

void updateG (  ) 

Update G and calculate total heat flux on boundary.

Definition at line 353 of file fvDOM.C.

References radiationModel::absorptionEmission_, and blackBodyEmission::correct().

Here is the call graph for this function:

void setRayIdLambdaId ( const word name,
label rayId,
label lambdaId 
) const

Set the rayId and lambdaId from by decomposing an intensity.

field name

Definition at line 369 of file fvDOM.C.

Foam::tmp< Foam::volScalarField > Rp (  )  const [virtual]

Source term component (for power of T^4).

Implements radiationModel.

Definition at line 308 of file fvDOM.C.

Foam::tmp< Foam::DimensionedField< Foam::scalar, Foam::volMesh > > Ru (  )  const [virtual]

Source term component (constant).

Implements radiationModel.

Definition at line 330 of file fvDOM.C.

const Foam::radiation::radiativeIntensityRay & IRay ( const label  rayI  )  const [inline]

Ray intensity for rayI.

Definition at line 19 of file fvDOMI.H.

const Foam::volScalarField & IRayLambda ( const label  rayI,
const label  lambdaI 
) const [inline]

Ray intensity for rayI and lambda bandwidth.

Definition at line 27 of file fvDOMI.H.

Foam::label nTheta (  )  const [inline]

Number of angles in theta.

Definition at line 36 of file fvDOMI.H.

Foam::label nPhi (  )  const [inline]

Number of angles in phi.

Definition at line 42 of file fvDOMI.H.

Foam::label nRay (  )  const [inline]

Number of rays.

Definition at line 48 of file fvDOMI.H.

Foam::label nLambda (  )  const [inline]

Number of wavelengths.

Definition at line 54 of file fvDOMI.H.

const Foam::volScalarField & a (  )  const [inline]

Const access to total absorption coefficient.

Definition at line 60 of file fvDOMI.H.

const Foam::volScalarField & aLambda ( const label  lambdaI  )  const [inline]

Const access to wavelength total absorption coefficient.

Definition at line 67 of file fvDOMI.H.

const Foam::volScalarField & G (  )  const [inline]

Const access to incident radiation field.

Definition at line 75 of file fvDOMI.H.

const Foam::volScalarField & Qr (  )  const [inline]

Const access to total radiative heat flux field.

Definition at line 81 of file fvDOMI.H.

const Foam::radiation::blackBodyEmission & blackBody (  )  const [inline]

Const access to black body.

Definition at line 88 of file fvDOMI.H.


The documentation for this class was generated from the following files:
  • src/thermophysicalModels/radiation/radiationModel/fvDOM/fvDOM/fvDOM.H
  • src/thermophysicalModels/radiation/radiationModel/fvDOM/fvDOM/fvDOM.C
  • src/thermophysicalModels/radiation/radiationModel/fvDOM/fvDOM/fvDOMI.H
Copyright © 2000-2009 OpenCFD Ltd