Forward compatibility: flex
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / RAS / LienCubicKE / LienCubicKE.H
blob95467863dd603cf343e4edae3e26cef824c5bacb
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::incompressible::RASModels::LienCubicKE
27 Description
28     Lien cubic non-linear k-epsilon turbulence model for incompressible flows.
30 SourceFiles
31     LienCubicKE.C
33 \*---------------------------------------------------------------------------*/
35 #ifndef LienCubicKE_H
36 #define LienCubicKE_H
38 #include "RASModel.H"
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 namespace Foam
44 namespace incompressible
46 namespace RASModels
49 /*---------------------------------------------------------------------------*\
50                            Class LienCubicKE Declaration
51 \*---------------------------------------------------------------------------*/
53 class LienCubicKE
55     public RASModel
57     // Private data
59         // Model coefficients
61             dimensionedScalar C1_;
62             dimensionedScalar C2_;
63             dimensionedScalar sigmak_;
64             dimensionedScalar sigmaEps_;
65             dimensionedScalar A1_;
66             dimensionedScalar A2_;
67             dimensionedScalar Ctau1_;
68             dimensionedScalar Ctau2_;
69             dimensionedScalar Ctau3_;
70             dimensionedScalar alphaKsi_;
73         // Fields
75             volScalarField k_;
76             volScalarField epsilon_;
78             volTensorField gradU_;
79             volScalarField eta_;
80             volScalarField ksi_;
81             volScalarField Cmu_;
82             volScalarField fEta_;
83             volScalarField C5viscosity_;
85             volScalarField nut_;
87             volSymmTensorField nonlinearStress_;
90 public:
92     //- Runtime type information
93     TypeName("LienCubicKE");
95     // Constructors
97         //- Construct from components
98         LienCubicKE
99         (
100             const volVectorField& U,
101             const surfaceScalarField& phi,
102             transportModel& transport
103         );
106     // Destructor
108         virtual ~LienCubicKE()
109         {}
112     // Member Functions
114         //- Return the turbulence viscosity
115         virtual tmp<volScalarField> nut() const
116         {
117             return nut_;
118         }
120         //- Return the effective diffusivity for k
121         tmp<volScalarField> DkEff() const
122         {
123             return tmp<volScalarField>
124             (
125                 new volScalarField("DkEff", nut_/sigmak_ + nu())
126             );
127         }
129         //- Return the effective diffusivity for epsilon
130         tmp<volScalarField> DepsilonEff() const
131         {
132             return tmp<volScalarField>
133             (
134                 new volScalarField("DepsilonEff", nut_/sigmaEps_ + nu())
135             );
136         }
138         //- Return the turbulence kinetic energy
139         virtual tmp<volScalarField> k() const
140         {
141             return k_;
142         }
144         //- Return the turbulence kinetic energy dissipation rate
145         virtual tmp<volScalarField> epsilon() const
146         {
147             return epsilon_;
148         }
150         //- Return the Reynolds stress tensor
151         virtual tmp<volSymmTensorField> R() const;
153         //- Return the effective stress tensor including the laminar stress
154         virtual tmp<volSymmTensorField> devReff() const;
156         //- Return the source term for the momentum equation
157         virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
159         //- Solve the turbulence equations and correct the turbulence viscosity
160         virtual void correct();
162         //- Read RASProperties dictionary
163         virtual bool read();
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 } // End namespace RASModels
170 } // End namespace incompressible
171 } // End namespace Foam
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 #endif
177 // ************************************************************************* //