ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / RAS / NonlinearKEShih / NonlinearKEShih.H
blob9c18ab71e6d395264d44a236e9242fa0a52458b7
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::NonlinearKEShih
27 Description
28     Shih's quadratic non-linear k-epsilon turbulence model for
29     incompressible flows
31 SourceFiles
32     NonlinearKEShih.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef NonlinearKEShih_H
37 #define NonlinearKEShih_H
39 #include "RASModel.H"
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 namespace Foam
45 namespace incompressible
47 namespace RASModels
50 /*---------------------------------------------------------------------------*\
51                       Class NonlinearKEShih Declaration
52 \*---------------------------------------------------------------------------*/
54 class NonlinearKEShih
56     public RASModel
59 protected:
61     // Protected data
63         // Model coefficients
65             dimensionedScalar C1_;
66             dimensionedScalar C2_;
67             dimensionedScalar sigmak_;
68             dimensionedScalar sigmaEps_;
69             dimensionedScalar A1_;
70             dimensionedScalar A2_;
71             dimensionedScalar Ctau1_;
72             dimensionedScalar Ctau2_;
73             dimensionedScalar Ctau3_;
74             dimensionedScalar alphaKsi_;
76             dimensionedScalar kappa_;
77             dimensionedScalar E_;
80         // Fields
82             volScalarField k_;
83             volScalarField epsilon_;
85             volTensorField gradU_;
86             volScalarField eta_;
87             volScalarField ksi_;
88             volScalarField Cmu_;
89             volScalarField fEta_;
91             volScalarField nut_;
93             volSymmTensorField nonlinearStress_;
96 public:
98     //- Runtime type information
99     TypeName("NonlinearKEShih");
101     // Constructors
103         //- Construct from components
104         NonlinearKEShih
105         (
106             const volVectorField& U,
107             const surfaceScalarField& phi,
108             transportModel& transport,
109             const word& turbulenceModelName = turbulenceModel::typeName,
110             const word& modelName = typeName
111         );
114     //- Destructor
115     virtual ~NonlinearKEShih()
116     {}
119     // Member Functions
121         //- Return the turbulence viscosity
122         virtual tmp<volScalarField> nut() const
123         {
124             return nut_;
125         }
127         //- Return the effective diffusivity for k
128         tmp<volScalarField> DkEff() const
129         {
130             return tmp<volScalarField>
131             (
132                 new volScalarField("DkEff", nut_/sigmak_ + nu())
133             );
134         }
136         //- Return the effective diffusivity for epsilon
137         tmp<volScalarField> DepsilonEff() const
138         {
139             return tmp<volScalarField>
140             (
141                 new volScalarField("DepsilonEff", nut_/sigmaEps_ + nu())
142             );
143         }
145         //- Return the turbulence kinetic energy
146         virtual tmp<volScalarField> k() const
147         {
148             return k_;
149         }
151         //- Return the turbulence kinetic energy dissipation rate
152         virtual tmp<volScalarField> epsilon() const
153         {
154             return epsilon_;
155         }
157         //- Return the Reynolds stress tensor
158         virtual tmp<volSymmTensorField> R() const;
160         //- Return the effective stress tensor including the laminar stress
161         virtual tmp<volSymmTensorField> devReff() const;
163         //- Return the source term for the momentum equation
164         virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
166         //- Solve the turbulence equations and correct the turbulence viscosity
167         virtual void correct();
169         //- Read RASProperties dictionary
170         virtual bool read();
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 } // End namespace RASModels
177 } // End namespace incompressible
178 } // End namespace Foam
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 #endif
184 // ************************************************************************* //