initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / turbulenceModels / RAS / compressible / RNGkEpsilon / RNGkEpsilon.H
blob5063615ccbe5a19fe87b52f8dda97713d8fc03be
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::compressible::RASModels::RNGkEpsilon
28 Description
29     Renormalisation group k-epsilon turbulence model for compressible 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
49     RNGkEpsilonCorrect.C
51 \*---------------------------------------------------------------------------*/
53 #ifndef compressibleRNGkEpsilon_H
54 #define compressibleRNGkEpsilon_H
56 #include "RASModel.H"
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 namespace Foam
62 namespace compressible
64 namespace RASModels
67 /*---------------------------------------------------------------------------*\
68                            Class RNGkEpsilon Declaration
69 \*---------------------------------------------------------------------------*/
71 class RNGkEpsilon
73     public RASModel
75     // Private data
77         dimensionedScalar Cmu_;
78         dimensionedScalar C1_;
79         dimensionedScalar C2_;
80         dimensionedScalar C3_;
81         dimensionedScalar alphak_;
82         dimensionedScalar alphaEps_;
83         dimensionedScalar alphah_;
84         dimensionedScalar eta0_;
85         dimensionedScalar beta_;
87         volScalarField k_;
88         volScalarField epsilon_;
89         volScalarField mut_;
92 public:
94     //- Runtime type information
95     TypeName("RNGkEpsilon");
97     // Constructors
99         //- from components
100         RNGkEpsilon
101         (
102             const volScalarField& rho,
103             const volVectorField& U,
104             const surfaceScalarField& phi,
105             basicThermo& thermophysicalModel
106         );
109     // Destructor
111         ~RNGkEpsilon(){}
114     // Member Functions
116         //- Return the turbulence viscosity
117         tmp<volScalarField> mut() const
118         {
119             return mut_;
120         }
122         //- Return the effective diffusivity for k
123         tmp<volScalarField> DkEff() const
124         {
125             return tmp<volScalarField>
126             (
127                 new volScalarField("DkEff", alphak_*mut_ + mu())
128             );
129         }
131         //- Return the effective diffusivity for epsilon
132         tmp<volScalarField> DepsilonEff() const
133         {
134             return tmp<volScalarField>
135             (
136                 new volScalarField("DepsilonEff", alphaEps_*mut_ + mu())
137             );
138         }
140         //- Return the effective turbulent thermal diffusivity
141         tmp<volScalarField> alphaEff() const
142         {
143             return tmp<volScalarField>
144             (
145                 new volScalarField("alphaEff", alphah_*mut_ + alpha())
146             );
147         }
149         //- Return the turbulence kinetic energy
150         tmp<volScalarField> k() const
151         {
152             return k_;
153         }
155         //- Return the turbulence kinetic energy dissipation rate
156         tmp<volScalarField> epsilon() const
157         {
158             return epsilon_;
159         }
161         //- Return the Reynolds stress tensor
162         tmp<volSymmTensorField> R() const;
164         //- Return the effective stress tensor including the laminar stress
165         tmp<volSymmTensorField> devRhoReff() const;
167         //- Return the effective stress tensor including the laminar stress
168         tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U) const;
170         //- Solve the turbulence equations and correct the turbulence viscosity
171         void correct();
173         //- Read turbulenceProperties dictionary
174         bool read();
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 } // End namespace RASModels
181 } // End namespace compressible
182 } // End namespace Foam
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 #endif
188 // ************************************************************************* //