ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / RAS / LamBremhorstKE / LamBremhorstKE.H
blob449055f4c5dde75b5f4e63ab52eda72b9cc89ee7
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::incompressible::RASModels::LamBremhorstKE
27 Description
28     Lam and Bremhorst low-Reynolds number k-epsilon turbulence model
29     for incompressible flows
31 SourceFiles
32     LamBremhorstKE.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef LamBremhorstKE_H
37 #define LamBremhorstKE_H
39 #include "RASModel.H"
40 #include "wallDist.H"
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 namespace Foam
46 namespace incompressible
48 namespace RASModels
51 /*---------------------------------------------------------------------------*\
52                       Class LamBremhorstKE Declaration
53 \*---------------------------------------------------------------------------*/
55 class LamBremhorstKE
57     public RASModel
60 protected:
62     // Protected data
64         dimensionedScalar Cmu_;
65         dimensionedScalar C1_;
66         dimensionedScalar C2_;
67         dimensionedScalar sigmaEps_;
69         volScalarField k_;
70         volScalarField epsilon_;
72         wallDist y_;
73         volScalarField Rt_;
75         volScalarField fMu_;
76         volScalarField nut_;
79 public:
81     //- Runtime type information
82     TypeName("LamBremhorstKE");
85     // Constructors
87         //- Construct from components
88         LamBremhorstKE
89         (
90             const volVectorField& U,
91             const surfaceScalarField& phi,
92             transportModel& transport,
93             const word& turbulenceModelName = turbulenceModel::typeName,
94             const word& modelName = typeName
95         );
98     //- Destructor
99     virtual ~LamBremhorstKE()
100     {}
103     // Member Functions
105         //- Return the turbulence viscosity
106         virtual tmp<volScalarField> nut() const
107         {
108             return nut_;
109         }
111         //- Return the effective diffusivity for k
112         tmp<volScalarField> DkEff() const
113         {
114             return tmp<volScalarField>
115             (
116                 new volScalarField("DkEff", nut_ + nu())
117             );
118         }
120         //- Return the effective diffusivity for epsilon
121         tmp<volScalarField> DepsilonEff() const
122         {
123             return tmp<volScalarField>
124             (
125                 new volScalarField("DepsilonEff", nut_/sigmaEps_ + nu())
126             );
127         }
129         //- Return the turbulence kinetic energy
130         virtual tmp<volScalarField> k() const
131         {
132             return k_;
133         }
135         //- Return the turbulence kinetic energy dissipation rate
136         virtual tmp<volScalarField> epsilon() const
137         {
138             return epsilon_;
139         }
141         //- Return the Reynolds stress tensor
142         virtual tmp<volSymmTensorField> R() const;
144         //- Return the effective stress tensor including the laminar stress
145         virtual tmp<volSymmTensorField> devReff() const;
147         //- Return the source term for the momentum equation
148         virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
150         //- Solve the turbulence equations and correct the turbulence viscosity
151         virtual void correct();
153         //- Read RASProperties dictionary
154         virtual bool read();
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 } // End namespace RASModels
161 } // End namespace incompressible
162 } // End namespace Foam
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 #endif
168 // ************************************************************************* //