initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / turbulenceModels / RAS / incompressible / RNGkEpsilon / RNGkEpsilon.H
blob960d0a344d740f510cfb03f00f6f9428e34f4b81
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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
25 Class
26     Foam::incompressible::RASModels::RNGkEpsilon
28 Description
29     Renormalisation group k-epsilon turbulence model for incompressible flows.
31     The default model coefficients correspond to the following:
32     @verbatim
33         RNGkEpsilonCoeffs
34         {
35             Cmu         0.0845;
36             C1          1.42;
37             C2          1.68;
38             C3          -0.33;  // only for compressible
39             alphah      1.0;    // only for compressible
40             alphak      1.39;
41             alphaEps    1.39;
42             eta0        4.38;
43             beta        0.012;
44         }
45     @endverbatim
47 SourceFiles
48     RNGkEpsilon.C
50 \*---------------------------------------------------------------------------*/
52 #ifndef RNGkEpsilon_H
53 #define RNGkEpsilon_H
55 #include "RASModel.H"
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 namespace Foam
61 namespace incompressible
63 namespace RASModels
66 /*---------------------------------------------------------------------------*\
67                            Class RNGkEpsilon Declaration
68 \*---------------------------------------------------------------------------*/
70 class RNGkEpsilon
72     public RASModel
74     // Private data
76         dimensionedScalar Cmu_;
77         dimensionedScalar C1_;
78         dimensionedScalar C2_;
79         dimensionedScalar alphak_;
80         dimensionedScalar alphaEps_;
81         dimensionedScalar eta0_;
82         dimensionedScalar beta_;
84         volScalarField k_;
85         volScalarField epsilon_;
86         volScalarField nut_;
89 public:
91     //- Runtime type information
92     TypeName("RNGkEpsilon");
94     // Constructors
96         //- from components
97         RNGkEpsilon
98         (
99             const volVectorField& U,
100             const surfaceScalarField& phi,
101             transportModel& transport
102         );
105     // Destructor
107         ~RNGkEpsilon()
108         {}
111     // Member Functions
113         //- Return the turbulence viscosity
114         tmp<volScalarField> nut() const
115         {
116             return nut_;
117         }
119         //- Return the effective diffusivity for k
120         tmp<volScalarField> DkEff() const
121         {
122             return tmp<volScalarField>
123             (
124                 new volScalarField("DkEff", alphak_*nut_ + nu())
125             );
126         }
128         //- Return the effective diffusivity for epsilon
129         tmp<volScalarField> DepsilonEff() const
130         {
131             return tmp<volScalarField>
132             (
133                 new volScalarField("DepsilonEff", alphaEps_*nut_ + nu())
134             );
135         }
137         //- Return the turbulence kinetic energy
138         tmp<volScalarField> k() const
139         {
140             return k_;
141         }
143         //- Return the turbulence kinetic energy dissipation rate
144         tmp<volScalarField> epsilon() const
145         {
146             return epsilon_;
147         }
149         //- Return the Reynolds stress tensor
150         tmp<volSymmTensorField> R() const;
152         //- Return the effective stress tensor including the laminar stress
153         tmp<volSymmTensorField> devReff() const;
155         //- Return the source term for the momentum equation
156         tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
158         //- Solve the turbulence equations and correct the turbulence viscosity
159         void correct();
161         //- Read turbulenceProperties dictionary
162         bool read();
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 } // End namespace RASModels
169 } // End namespace incompressible
170 } // End namespace Foam
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 #endif
176 // ************************************************************************* //