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

RASModel Class Reference

Abstract base class for turbulence models for compressible and combusting flows. More...

Inheritance diagram for RASModel:
Collaboration diagram for RASModel:

List of all members.


Public Member Functions

 TypeName ("RASModel")
 Runtime type information.
 declareRunTimeSelectionTable (autoPtr, RASModel, dictionary,(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel),(rho, U, phi, thermoPhysicalModel))
 RASModel (const word &type, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel)
 Construct from components.
virtual ~RASModel ()
 Destructor.
const dimensionedScalark0 () const
 Return the value of k0 which k is not allowed to be less than.
const dimensionedScalarepsilon0 () const
 Return the value of epsilon0 which epsilon is not allowed to be.
const dimensionedScalarepsilonSmall () const
 Return the value of epsilonSmall which is added to epsilon when.
const dimensionedScalaromega0 () const
 Return the value of omega0 which epsilon is not allowed to be.
const dimensionedScalaromegaSmall () const
 Return the value of omegaSmall which is added to epsilon when.
dimensionedScalark0 ()
 Allow k0 to be changed.
dimensionedScalarepsilon0 ()
 Allow epsilon0 to be changed.
dimensionedScalarepsilonSmall ()
 Allow epsilonSmall to be changed.
dimensionedScalaromega0 ()
 Allow omega0 to be changed.
dimensionedScalaromegaSmall ()
 Allow omegaSmall to be changed.
const nearWallDisty () const
 Return the near wall distances.
scalar yPlusLam (const scalar kappa, const scalar E) const
 Calculate y+ at the edge of the laminar sublayer.
const dictionarycoeffDict () const
 Const access to the coefficients dictionary.
virtual tmp< volScalarFieldmut () const =0
 Return the turbulence viscosity.
virtual tmp< volScalarFieldmuEff () const
 Return the effective viscosity.
virtual tmp< volScalarFieldalphaEff () const =0
 Return the effective turbulent thermal diffusivity.
virtual tmp< volScalarFieldk () const =0
 Return the turbulence kinetic energy.
virtual tmp< volScalarFieldepsilon () const =0
 Return the turbulence kinetic energy dissipation rate.
virtual tmp< volSymmTensorFieldR () const =0
 Return the Reynolds stress tensor.
virtual tmp< volSymmTensorFielddevRhoReff () const =0
 Return the effective stress tensor including the laminar stress.
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const =0
 Return the source term for the momentum equation.
virtual tmp< scalarFieldyPlus (const label patchI, const scalar Cmu) const
 Return yPlus for the given patch.
virtual void correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity.
virtual bool read ()=0
 Read RASProperties dictionary.

Static Public Member Functions

static autoPtr< RASModelNew (const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel)
 Return a reference to the selected turbulence model.

Protected Member Functions

virtual void printCoeffs ()
 Print model coefficients.

Protected Attributes

Switch turbulence_
 Turbulence on/off flag.
Switch printCoeffs_
 Flag to print the model coeffs at run-time.
dictionary coeffDict_
 Model coefficients dictionary.
scalar yPlusLam_
 Value of y+ at the edge of the laminar sublayer.
dimensionedScalar k0_
 Lower limit of k.
dimensionedScalar epsilon0_
 Lower limit of epsilon.
dimensionedScalar epsilonSmall_
 Small epsilon value used to avoid divide by zero.
dimensionedScalar omega0_
 Lower limit for omega.
dimensionedScalar omegaSmall_
 Small omega value used to avoid divide by zero.
nearWallDist y_
 Near wall distance boundary field.

Detailed Description

Abstract base class for turbulence models for compressible and combusting flows.

Source files

Definition at line 67 of file RASModel.H.


Constructor & Destructor Documentation

RASModel ( const word type,
const volScalarField rho,
const volVectorField U,
const surfaceScalarField phi,
const basicThermo thermoPhysicalModel 
)

Construct from components.

Definition at line 49 of file RASModel.C.

References RASModel::coeffDict_, Foam::endl(), Foam::Info, and Foam::type().

Here is the call graph for this function:

virtual ~RASModel (  )  [inline, virtual]

Destructor.

Definition at line 174 of file RASModel.H.


Member Function Documentation

void printCoeffs (  )  [protected, virtual]

Print model coefficients.

Definition at line 37 of file RASModel.C.

Referenced by LaunderGibsonRSTM::devRhoReff(), LRR::devRhoReff(), RNGkEpsilon::R(), kOmegaSST::R(), realizableKE::R(), LaunderSharmaKE::R(), SpalartAllmaras::R(), and kEpsilon::R().

Here is the caller graph for this function:

TypeName ( "RASModel"   ) 

Runtime type information.

declareRunTimeSelectionTable ( autoPtr  ,
RASModel  ,
dictionary  ,
(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel)  ,
(rho, U, phi, thermoPhysicalModel)   
)

autoPtr< RASModel > New ( const volScalarField rho,
const volVectorField U,
const surfaceScalarField phi,
const basicThermo thermoPhysicalModel 
) [static]

Return a reference to the selected turbulence model.

Definition at line 92 of file RASModel.C.

References turbulenceModel::mesh_.

const dimensionedScalar& k0 (  )  const [inline]

Return the value of k0 which k is not allowed to be less than.

Definition at line 183 of file RASModel.H.

References RASModel::k0(), and RASModel::k0_.

Referenced by RASModel::k0().

Here is the call graph for this function:

Here is the caller graph for this function:

const dimensionedScalar& epsilon0 (  )  const [inline]

Return the value of epsilon0 which epsilon is not allowed to be.

less than

Definition at line 192 of file RASModel.H.

const dimensionedScalar& epsilonSmall (  )  const [inline]

Return the value of epsilonSmall which is added to epsilon when.

calculating nut

Definition at line 201 of file RASModel.H.

References RASModel::epsilonSmall(), and RASModel::epsilonSmall_.

Referenced by RASModel::epsilonSmall().

Here is the call graph for this function:

Here is the caller graph for this function:

const dimensionedScalar& omega0 (  )  const [inline]

Return the value of omega0 which epsilon is not allowed to be.

less than

Definition at line 210 of file RASModel.H.

References RASModel::omega0(), and RASModel::omega0_.

Referenced by RASModel::omega0().

Here is the call graph for this function:

Here is the caller graph for this function:

const dimensionedScalar& omegaSmall (  )  const [inline]

Return the value of omegaSmall which is added to epsilon when.

calculating nut

Definition at line 219 of file RASModel.H.

References RASModel::omegaSmall_.

dimensionedScalar& k0 (  )  [inline]

Allow k0 to be changed.

Definition at line 225 of file RASModel.H.

References RASModel::k0_.

dimensionedScalar& epsilon0 (  )  [inline]

Allow epsilon0 to be changed.

Definition at line 231 of file RASModel.H.

References RASModel::epsilon0_.

dimensionedScalar& epsilonSmall (  )  [inline]

Allow epsilonSmall to be changed.

Definition at line 237 of file RASModel.H.

References RASModel::epsilonSmall_.

dimensionedScalar& omega0 (  )  [inline]

Allow omega0 to be changed.

Definition at line 243 of file RASModel.H.

References RASModel::omega0_.

dimensionedScalar& omegaSmall (  )  [inline]

Allow omegaSmall to be changed.

Definition at line 249 of file RASModel.H.

References RASModel::omegaSmall_.

const nearWallDist& y (  )  const [inline]

scalar yPlusLam ( const scalar  kappa,
const scalar  E 
) const

Calculate y+ at the edge of the laminar sublayer.

Definition at line 148 of file RASModel.C.

const dictionary& coeffDict (  )  const [inline]

Const access to the coefficients dictionary.

Definition at line 264 of file RASModel.H.

References RASModel::coeffDict_.

virtual tmp<volScalarField> mut (  )  const [pure virtual]

Return the turbulence viscosity.

Implements turbulenceModel.

Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.

Referenced by RASModel::muEff().

Here is the caller graph for this function:

virtual tmp<volScalarField> muEff (  )  const [inline, virtual]

Return the effective viscosity.

Implements turbulenceModel.

Reimplemented in laminar.

Definition at line 274 of file RASModel.H.

References RASModel::mut().

Here is the call graph for this function:

virtual tmp<volScalarField> alphaEff (  )  const [pure virtual]

Return the effective turbulent thermal diffusivity.

Implements turbulenceModel.

Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.

virtual tmp<volScalarField> k (  )  const [pure virtual]

Return the turbulence kinetic energy.

Implements turbulenceModel.

Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.

virtual tmp<volScalarField> epsilon (  )  const [pure virtual]

Return the turbulence kinetic energy dissipation rate.

Implements turbulenceModel.

Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.

virtual tmp<volSymmTensorField> R (  )  const [pure virtual]

virtual tmp<volSymmTensorField> devRhoReff (  )  const [pure virtual]

Return the effective stress tensor including the laminar stress.

Implements turbulenceModel.

Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.

virtual tmp<fvVectorMatrix> divDevRhoReff ( volVectorField U  )  const [pure virtual]

Return the source term for the momentum equation.

Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.

tmp< scalarField > yPlus ( const label  patchI,
const scalar  Cmu 
) const [virtual]

Return yPlus for the given patch.

Definition at line 161 of file RASModel.C.

References Foam::log().

Here is the call graph for this function:

void correct (  )  [pure virtual]

Solve the turbulence equations and correct the turbulence viscosity.

Implements turbulenceModel.

Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.

Definition at line 193 of file RASModel.C.

bool read (  )  [pure virtual]

Read RASProperties dictionary.

Implements turbulenceModel.

Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.

Definition at line 202 of file RASModel.C.

References polyMesh::changing(), nearWallDist::correct(), turbulenceModel::mesh_, and RASModel::y_.

Here is the call graph for this function:


Member Data Documentation

Switch turbulence_ [protected]

Turbulence on/off flag.

Definition at line 78 of file RASModel.H.

Switch printCoeffs_ [protected]

Flag to print the model coeffs at run-time.

Definition at line 81 of file RASModel.H.

dictionary coeffDict_ [protected]

Model coefficients dictionary.

Definition at line 84 of file RASModel.H.

Referenced by RASModel::coeffDict(), and RASModel::RASModel().

scalar yPlusLam_ [protected]

Value of y+ at the edge of the laminar sublayer.

Definition at line 87 of file RASModel.H.

dimensionedScalar k0_ [protected]

Lower limit of k.

Definition at line 90 of file RASModel.H.

Referenced by RASModel::k0().

dimensionedScalar epsilon0_ [protected]

Lower limit of epsilon.

Definition at line 93 of file RASModel.H.

Referenced by RASModel::epsilon0().

dimensionedScalar epsilonSmall_ [protected]

Small epsilon value used to avoid divide by zero.

Definition at line 96 of file RASModel.H.

Referenced by RASModel::epsilonSmall(), kEpsilon::kEpsilon(), and RNGkEpsilon::RNGkEpsilon().

dimensionedScalar omega0_ [protected]

Lower limit for omega.

Definition at line 99 of file RASModel.H.

Referenced by RASModel::omega0().

dimensionedScalar omegaSmall_ [protected]

Small omega value used to avoid divide by zero.

Definition at line 102 of file RASModel.H.

Referenced by RASModel::omegaSmall().

nearWallDist y_ [protected]

Near wall distance boundary field.

Definition at line 105 of file RASModel.H.

Referenced by RASModel::read(), and RASModel::y().


The documentation for this class was generated from the following files:
  • src/turbulenceModels/compressible/RAS/RASModel/RASModel.H
  • src/turbulenceModels/compressible/RAS/RASModel/RASModel.C
Copyright © 2000-2009 OpenCFD Ltd