FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

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

#include <compressibleRASModels/RASModel.H>


Detailed Description

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

Source files

Definition at line 71 of file RASModel.H.

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 dimensionedScalar &  k0 () const
 Return the value of k0 which k is not allowed to be less than.
const dimensionedScalar &  epsilon0 () const
 Return the value of epsilon0 which epsilon is not allowed to be.
const dimensionedScalar &  epsilonSmall () const
 Return the value of epsilonSmall which is added to epsilon when.
const dimensionedScalar &  omega0 () const
 Return the value of omega0 which epsilon is not allowed to be.
const dimensionedScalar &  omegaSmall () const
 Return the value of omegaSmall which is added to epsilon when.
dimensionedScalar &  k0 ()
 Allow k0 to be changed.
dimensionedScalar &  epsilon0 ()
 Allow epsilon0 to be changed.
dimensionedScalar &  epsilonSmall ()
 Allow epsilonSmall to be changed.
dimensionedScalar &  omega0 ()
 Allow omega0 to be changed.
dimensionedScalar &  omegaSmall ()
 Allow omegaSmall to be changed.
const nearWallDist &  y () 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 dictionary &  coeffDict () const
 Const access to the coefficients dictionary.
virtual tmp< volScalarField >  mut () const =0
 Return the turbulence viscosity.
virtual tmp< volScalarField >  muEff () const
 Return the effective viscosity.
virtual tmp< volScalarField >  alphaEff () const =0
 Return the effective turbulent thermal diffusivity.
virtual tmp< volScalarField >  k () const =0
 Return the turbulence kinetic energy.
virtual tmp< volScalarField >  epsilon () const =0
 Return the turbulence kinetic energy dissipation rate.
virtual tmp< volSymmTensorField >  R () const =0
 Return the Reynolds stress tensor.
virtual tmp< volSymmTensorField >  devRhoReff () const =0
 Return the effective stress tensor including the laminar stress.
virtual tmp< fvVectorMatrix >  divDevRhoReff (volVectorField &U) const =0
 Return the source term for the momentum equation.
virtual tmp< scalarField >  yPlus (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< RASModel >  New (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.

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 57 of file RASModel.C.

References dictionary::readIfPresent().

virtual ~RASModel (  ) [inline, virtual]

Destructor.

Definition at line 178 of file RASModel.H.


Member Function Documentation

void printCoeffs (  ) [protected, virtual]

Print model coefficients.

Definition at line 45 of file RASModel.C.

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

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 106 of file RASModel.C.

References TimePaths::constant(), IOobject::db(), Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorIn, Foam::Info, dictionary::lookup(), IOobject::MUST_READ, IOobject::NO_WRITE, phi, rho, IOobject::time(), and U.

const dimensionedScalar& k0 (  ) const [inline]

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

Definition at line 187 of file RASModel.H.

References RASModel::k0_.

const dimensionedScalar& epsilon0 (  ) const [inline]

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

less than

Definition at line 194 of file RASModel.H.

References RASModel::epsilon0_.

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_.

const dimensionedScalar& omega0 (  ) const [inline]

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

less than

Definition at line 208 of file RASModel.H.

References RASModel::omega0_.

const dimensionedScalar& omegaSmall (  ) const [inline]

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

calculating nut

Definition at line 215 of file RASModel.H.

References RASModel::omegaSmall_.

dimensionedScalar& k0 (  ) [inline]

Allow k0 to be changed.

Definition at line 221 of file RASModel.H.

References RASModel::k0_.

dimensionedScalar& epsilon0 (  ) [inline]

Allow epsilon0 to be changed.

Definition at line 227 of file RASModel.H.

References RASModel::epsilon0_.

dimensionedScalar& epsilonSmall (  ) [inline]

Allow epsilonSmall to be changed.

Definition at line 233 of file RASModel.H.

References RASModel::epsilonSmall_.

dimensionedScalar& omega0 (  ) [inline]

Allow omega0 to be changed.

Definition at line 239 of file RASModel.H.

References RASModel::omega0_.

dimensionedScalar& omegaSmall (  ) [inline]

Allow omegaSmall to be changed.

Definition at line 245 of file RASModel.H.

References RASModel::omegaSmall_.

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

Calculate y+ at the edge of the laminar sublayer.

Definition at line 162 of file RASModel.C.

References kappa(), Foam::log(), and Foam::max().

Referenced by mutRoughWallFunctionFvPatchScalarField::calcMut().

const dictionary& coeffDict (  ) const [inline]
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().

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

Return the effective viscosity.

Implements turbulenceModel.

Reimplemented in laminar.

Definition at line 270 of file RASModel.H.

References mu, and RASModel::mut().

virtual tmp<volScalarField> alphaEff (  ) const [pure virtual]
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]

Return the Reynolds stress tensor.

Implements turbulenceModel.

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

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]
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 207 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 216 of file RASModel.C.

References regIOobject::read(), and Foam::type().


Member Data Documentation

Switch turbulence_ [protected]

Turbulence on/off flag.

Definition at line 82 of file RASModel.H.

Switch printCoeffs_ [protected]

Flag to print the model coeffs at run-time.

Definition at line 85 of file RASModel.H.

dictionary coeffDict_ [protected]

Model coefficients dictionary.

Definition at line 88 of file RASModel.H.

Referenced by RASModel::coeffDict().

scalar yPlusLam_ [protected]

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

Definition at line 91 of file RASModel.H.

dimensionedScalar k0_ [protected]

Lower limit of k.

Definition at line 94 of file RASModel.H.

Referenced by RASModel::k0().

dimensionedScalar epsilon0_ [protected]

Lower limit of epsilon.

Definition at line 97 of file RASModel.H.

Referenced by RASModel::epsilon0().

dimensionedScalar epsilonSmall_ [protected]

Small epsilon value used to avoid divide by zero.

Definition at line 100 of file RASModel.H.

Referenced by RASModel::epsilonSmall().

dimensionedScalar omega0_ [protected]

Lower limit for omega.

Definition at line 103 of file RASModel.H.

Referenced by RASModel::omega0().

dimensionedScalar omegaSmall_ [protected]

Small omega value used to avoid divide by zero.

Definition at line 106 of file RASModel.H.

Referenced by RASModel::omegaSmall().

nearWallDist y_ [protected]

Near wall distance boundary field.

Definition at line 109 of file RASModel.H.

Referenced by 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