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 incompressible turbulence models. More...

#include <incompressibleRASModels/RASModel.H>


Detailed Description

Abstract base class for incompressible turbulence models.

Source files

Definition at line 68 of file RASModel.H.

Inheritance diagram for incompressible::RASModel:
Collaboration diagram for incompressible::RASModel:

List of all members.

Public Member Functions

 TypeName ("RASModel")
 Runtime type information.
 declareRunTimeSelectionTable (autoPtr, RASModel, dictionary,(const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel),(U, phi, lamTransportModel))
 RASModel (const word &type, const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel)
 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 >  nut () const =0
 Return the turbulence viscosity.
virtual tmp< volScalarField >  nuEff () const
 Return the effective viscosity.
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 >  devReff () const =0
 Return the effective stress tensor including the laminar stress.
virtual tmp< fvVectorMatrix >  divDevReff (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 volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel)
 Return a reference to the selected RAS 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.
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 volVectorField &   U,
const surfaceScalarField &   phi,
transportModel &   lamTransportModel  
)

Construct from components.

Definition at line 57 of file RASModel.C.

References dictionary::readIfPresent().

virtual ~RASModel (  ) [inline, virtual]

Destructor.

Definition at line 169 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 volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel)   ,
(U, phi, lamTransportModel)    
)
autoPtr< RASModel > New ( const volVectorField &   U,
const surfaceScalarField &   phi,
transportModel &   lamTransportModel  
) [static]
const dimensionedScalar& k0 (  ) const [inline]

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

Definition at line 178 of file RASModel.H.

References incompressible::RASModel::k0_.

const dimensionedScalar& epsilon0 (  ) const [inline]

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

less than

Definition at line 185 of file RASModel.H.

References incompressible::RASModel::epsilon0_.

const dimensionedScalar& epsilonSmall (  ) const [inline]

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

calculating nut

Definition at line 192 of file RASModel.H.

References incompressible::RASModel::epsilonSmall_.

const dimensionedScalar& omega0 (  ) const [inline]

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

less than

Definition at line 199 of file RASModel.H.

References incompressible::RASModel::omega0_.

const dimensionedScalar& omegaSmall (  ) const [inline]

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

calculating nut

Definition at line 206 of file RASModel.H.

References incompressible::RASModel::omegaSmall_.

dimensionedScalar& k0 (  ) [inline]

Allow k0 to be changed.

Definition at line 212 of file RASModel.H.

References incompressible::RASModel::k0_.

dimensionedScalar& epsilon0 (  ) [inline]

Allow epsilon0 to be changed.

Definition at line 218 of file RASModel.H.

References incompressible::RASModel::epsilon0_.

dimensionedScalar& epsilonSmall (  ) [inline]

Allow epsilonSmall to be changed.

Definition at line 224 of file RASModel.H.

References incompressible::RASModel::epsilonSmall_.

dimensionedScalar& omega0 (  ) [inline]

Allow omega0 to be changed.

Definition at line 230 of file RASModel.H.

References incompressible::RASModel::omega0_.

dimensionedScalar& omegaSmall (  ) [inline]

Allow omegaSmall to be changed.

Definition at line 236 of file RASModel.H.

References incompressible::RASModel::omegaSmall_.

const nearWallDist& y (  ) const [inline]

Return the near wall distances.

Definition at line 242 of file RASModel.H.

References incompressible::RASModel::y_.

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

Calculate y+ at the edge of the laminar sublayer.

Definition at line 157 of file RASModel.C.

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

virtual tmp<volScalarField> epsilon (  ) const [pure virtual]
virtual tmp<fvVectorMatrix> divDevReff ( volVectorField &   U  ) const [pure virtual]
tmp< scalarField > yPlus ( const label   patchI,
const scalar   Cmu  
) const [virtual]

Return yPlus for the given patch.

Definition at line 170 of file RASModel.C.

References boundaryField(), Foam::endl(), k(), Foam::nl, nu, Foam::pow(), List< T >::setSize(), fvPatch::size(), Foam::sqrt(), and WarningIn.

void correct (  ) [pure virtual]

Solve the turbulence equations and correct the turbulence viscosity.

Implements incompressible::turbulenceModel.

Implemented in kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, realizableKE, RNGkEpsilon, and SpalartAllmaras.

Definition at line 199 of file RASModel.C.

References correct().

bool read (  ) [pure virtual]

Member Data Documentation

Switch printCoeffs_ [protected]

Flag to print the model coeffs at run-time.

Definition at line 82 of file RASModel.H.

dictionary coeffDict_ [protected]

Model coefficients dictionary.

Definition at line 85 of file RASModel.H.

Referenced by incompressible::RASModel::coeffDict().

dimensionedScalar epsilonSmall_ [protected]

Small epsilon value used to avoid divide by zero.

Definition at line 94 of file RASModel.H.

Referenced by incompressible::RASModel::epsilonSmall().

dimensionedScalar omega0_ [protected]

Lower limit for omega.

Definition at line 97 of file RASModel.H.

Referenced by kOmega::correct(), kOmegaSST::correct(), and incompressible::RASModel::omega0().

dimensionedScalar omegaSmall_ [protected]

Small omega value used to avoid divide by zero.

Definition at line 100 of file RASModel.H.

Referenced by kOmega::correct(), and incompressible::RASModel::omegaSmall().

nearWallDist y_ [protected]

Near wall distance boundary field.

Definition at line 103 of file RASModel.H.

Referenced by incompressible::RASModel::y().


The documentation for this class was generated from the following files:
  • src/turbulenceModels/incompressible/RAS/RASModel/RASModel.H
  • src/turbulenceModels/incompressible/RAS/RASModel/RASModel.C