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

realizableKE Class Reference

Realizable k-epsilon turbulence model for incompressible flows. More...

#include <incompressibleRASModels/realizableKE.H>


Detailed Description

Realizable k-epsilon turbulence model for incompressible flows.

Model described in the paper:

"A New k-epsilon Eddy Viscosity Model for High Reynolds Number
Turbulent Flows"

Tsan-Hsing Shih, William W. Liou, Aamir Shabbir, Zhigang Tang and
Jiang Zhu

Computers and Fluids Vol. 24, No. 3, pp. 227-238, 1995

The default model coefficients correspond to the following:

realizableKECoeffs
{
Cmu 0.09;
A0  4.0;
C2  1.9;
sigmak  1.0;
sigmaEps1.2;
}
Source files

Definition at line 76 of file realizableKE.H.

Inheritance diagram for realizableKE:
Collaboration diagram for realizableKE:

List of all members.

Public Member Functions

 TypeName ("realizableKE")
 Runtime type information.
 realizableKE (const volVectorField &U, const surfaceScalarField &phi, transportModel &transport)
 Construct from components.
virtual  ~realizableKE ()
 Destructor.
virtual tmp< volScalarField >  nut () const
 Return the turbulence viscosity.
tmp< volScalarField >  DkEff () const
 Return the effective diffusivity for k.
tmp< volScalarField >  DepsilonEff () const
 Return the effective diffusivity for epsilon.
virtual tmp< volScalarField >  k () const
 Return the turbulence kinetic energy.
virtual tmp< volScalarField >  epsilon () const
 Return the turbulence kinetic energy dissipation rate.
virtual tmp< volSymmTensorField >  R () const
 Return the Reynolds stress tensor.
virtual tmp< volSymmTensorField >  devReff () const
 Return the effective stress tensor including the laminar stress.
virtual tmp< fvVectorMatrix >  divDevReff (volVectorField &U) const
 Return the source term for the momentum equation.
virtual void  correct ()
 Solve the turbulence equations and correct the turbulence viscosity.
virtual bool  read ()
 Read RASProperties dictionary.

Constructor & Destructor Documentation

realizableKE ( const volVectorField &   U,
const surfaceScalarField &   phi,
transportModel &   transport  
)

Construct from components.

Definition at line 89 of file realizableKE.C.

References Foam::bound(), Foam::fvc::grad(), and Foam::sqr().

virtual ~realizableKE (  ) [inline, virtual]

Destructor.

Definition at line 127 of file realizableKE.H.


Member Function Documentation

TypeName ( "realizableKE"    )

Runtime type information.

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

Return the turbulence viscosity.

Implements incompressible::RASModel.

Definition at line 134 of file realizableKE.H.

tmp<volScalarField> DkEff (  ) const [inline]

Return the effective diffusivity for k.

Definition at line 140 of file realizableKE.H.

References incompressible::turbulenceModel::nu().

Referenced by realizableKE::correct().

tmp<volScalarField> DepsilonEff (  ) const [inline]

Return the effective diffusivity for epsilon.

Definition at line 149 of file realizableKE.H.

References incompressible::turbulenceModel::nu().

Referenced by realizableKE::correct().

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

Return the turbulence kinetic energy.

Implements incompressible::RASModel.

Definition at line 158 of file realizableKE.H.

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

Return the turbulence kinetic energy dissipation rate.

Implements incompressible::RASModel.

Definition at line 164 of file realizableKE.H.

tmp< fvVectorMatrix > divDevReff ( volVectorField &   U  ) const [virtual]

Return the source term for the momentum equation.

Implements incompressible::RASModel.

Definition at line 233 of file realizableKE.C.

References Foam::dev(), Foam::fvc::div(), Foam::fvc::grad(), Foam::fvm::laplacian(), incompressible::RASModel::nuEff(), and Foam::T().

bool read (  ) [virtual]

Read RASProperties dictionary.

Implements incompressible::RASModel.

Definition at line 243 of file realizableKE.C.

References incompressible::RASModel::coeffDict(), and dimensioned< Type >::readIfPresent().


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