initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / turbulenceModels / RAS / incompressible / LienCubicKE / LienCubicKE.H
blob005c8257751e255d007557c6210cd3b5abf3ee8f
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::LienCubicKE
28 Description
29     Lien cubic non-linear k-epsilon turbulence model for incompressible flows.
31 SourceFiles
32     LienCubicKE.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef LienCubicKE_H
37 #define LienCubicKE_H
39 #include "RASModel.H"
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 namespace Foam
45 namespace incompressible
47 namespace RASModels
50 /*---------------------------------------------------------------------------*\
51                            Class LienCubicKE Declaration
52 \*---------------------------------------------------------------------------*/
54 class LienCubicKE
56     public RASModel
58     // Private data
60         dimensionedScalar C1_;
61         dimensionedScalar C2_;
62         dimensionedScalar alphak_;
63         dimensionedScalar alphaEps_;
64         dimensionedScalar A1_;
65         dimensionedScalar A2_;
66         dimensionedScalar Ctau1_;
67         dimensionedScalar Ctau2_;
68         dimensionedScalar Ctau3_;
69         dimensionedScalar alphaKsi_;
71         volScalarField k_;
72         volScalarField epsilon_;
74         volTensorField gradU_;
75         volScalarField eta_;
76         volScalarField ksi_;
77         volScalarField Cmu_;
78         volScalarField fEta_;
79         volScalarField C5viscosity_;
81         volScalarField nut_;
83         volSymmTensorField nonlinearStress_;
86 public:
88     //- Runtime type information
89     TypeName("LienCubicKE");
91     // Constructors
93         //- from components
94         LienCubicKE
95         (
96             const volVectorField& U,
97             const surfaceScalarField& phi,
98             transportModel& transport
99         );
102     // Destructor
104         ~LienCubicKE()
105         {}
108     // Member Functions
110         //- Return the turbulence viscosity
111         tmp<volScalarField> nut() const
112         {
113             return nut_;
114         }
116         //- Return the effective diffusivity for k
117         tmp<volScalarField> DkEff() const
118         {
119             return tmp<volScalarField>
120             (
121                 new volScalarField("DkEff", alphak_*nut_ + nu())
122             );
123         }
125         //- Return the effective diffusivity for epsilon
126         tmp<volScalarField> DepsilonEff() const
127         {
128             return tmp<volScalarField>
129             (
130                 new volScalarField("DepsilonEff", alphaEps_*nut_ + nu())
131             );
132         }
134         //- Return the turbulence kinetic energy
135         tmp<volScalarField> k() const
136         {
137             return k_;
138         }
140         //- Return the turbulence kinetic energy dissipation rate
141         tmp<volScalarField> epsilon() const
142         {
143             return epsilon_;
144         }
146         //- Return the Reynolds stress tensor
147         tmp<volSymmTensorField> R() const;
149         //- Return the effective stress tensor including the laminar stress
150         tmp<volSymmTensorField> devReff() const;
152         //- Return the source term for the momentum equation
153         tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
155         //- Solve the turbulence equations and correct the turbulence viscosity
156         void correct();
158         //- Read turbulenceProperties dictionary
159         bool read();
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 } // End namespace RASModels
166 } // End namespace incompressible
167 } // End namespace Foam
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 #endif
173 // ************************************************************************* //