1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software; you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
26 Foam::incompressible::RASModels::kEpsilon
29 Standard k-epsilon turbulence model for incompressible flows.
31 The default model coefficients correspond to the following:
45 \*---------------------------------------------------------------------------*/
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 namespace incompressible
61 /*---------------------------------------------------------------------------*\
62 Class kEpsilon Declaration
63 \*---------------------------------------------------------------------------*/
73 dimensionedScalar Cmu_;
74 dimensionedScalar C1_;
75 dimensionedScalar C2_;
76 dimensionedScalar sigmaEps_;
82 volScalarField epsilon_;
88 //- Runtime type information
93 //- Construct from components
96 const volVectorField& U,
97 const surfaceScalarField& phi,
98 transportModel& transport
109 //- Return the turbulence viscosity
110 virtual tmp<volScalarField> nut() const
115 //- Return the effective diffusivity for k
116 tmp<volScalarField> DkEff() const
118 return tmp<volScalarField>
120 new volScalarField("DkEff", nut_ + nu())
124 //- Return the effective diffusivity for epsilon
125 tmp<volScalarField> DepsilonEff() const
127 return tmp<volScalarField>
129 new volScalarField("DepsilonEff", nut_/sigmaEps_ + nu())
133 //- Return the turbulence kinetic energy
134 virtual tmp<volScalarField> k() const
139 //- Return the turbulence kinetic energy dissipation rate
140 virtual tmp<volScalarField> epsilon() const
145 //- Return the Reynolds stress tensor
146 virtual tmp<volSymmTensorField> R() const;
148 //- Return the effective stress tensor including the laminar stress
149 virtual tmp<volSymmTensorField> devReff() const;
151 //- Return the source term for the momentum equation
152 virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
154 //- Solve the turbulence equations and correct the turbulence viscosity
155 virtual void correct();
157 //- Read RASProperties dictionary
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 } // End namespace RASModels
165 } // End namespace incompressible
166 } // End namespace Foam
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 // ************************************************************************* //