ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / RAS / qZeta / qZeta.H
blobba8dc1b413f76bf329855600808da3e17c2da946
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::qZeta
27 Description
28     Gibson and Dafa'Alla's q-zeta two-equation low-Re turbulence model
29     for incompressible flows
31     References:
32     \verbatim
33         "Calculation of oscillating boundary layers with the
34         q-[zeta] turbulence model"
35         A.A. Dafa'Alla, E. Juntasaro and M.M. Gibson
37         Engineering Turbulence Modelling and Experiments 3:
38         Proceedings of the Third International Symposium,
39         Crete, Greece, May 27-29, 1996, p. 141.
41         Editors: Wolfgang Rodi and G. Bergeles
42     \endverbatim
44 SourceFiles
45     qZeta.C
47 \*---------------------------------------------------------------------------*/
49 #ifndef qZeta_H
50 #define qZeta_H
52 #include "RASModel.H"
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 namespace Foam
58 namespace incompressible
60 namespace RASModels
63 /*---------------------------------------------------------------------------*\
64                            Class qZeta Declaration
65 \*---------------------------------------------------------------------------*/
67 class qZeta
69     public RASModel
72 protected:
74     // Protected data
76         // Model coefficients
78             dimensionedScalar Cmu_;
79             dimensionedScalar C1_;
80             dimensionedScalar C2_;
81             dimensionedScalar sigmaZeta_;
82             Switch anisotropic_;
84         //- Lower limit of q
85         dimensionedScalar qMin_;
87         //- Lower limit of zeta
88         dimensionedScalar zetaMin_;
90         // Fields
92             volScalarField k_;
93             volScalarField epsilon_;
95             volScalarField q_;
96             volScalarField zeta_;
98             volScalarField nut_;
101     // Protected Member Functions
103         tmp<volScalarField> fMu() const;
104         tmp<volScalarField> f2() const;
107 public:
109     //- Runtime type information
110     TypeName("qZeta");
112     // Constructors
114         //- Construct from components
115         qZeta
116         (
117             const volVectorField& U,
118             const surfaceScalarField& phi,
119             transportModel& transport,
120             const word& turbulenceModelName = turbulenceModel::typeName,
121             const word& modelName = typeName
122         );
125     //- Destructor
126     virtual ~qZeta()
127     {}
130     // Member Functions
132        // Access
134             //- Return the lower allowable limit for q (default: SMALL)
135             const dimensionedScalar& qMin() const
136             {
137                 return qMin_;
138             }
140             //- Return the lower allowable limit for zeta (default: SMALL)
141             const dimensionedScalar& zetaMin() const
142             {
143                 return zetaMin_;
144             }
146             //- Allow qMin to be changed
147             dimensionedScalar& qMin()
148             {
149                 return qMin_;
150             }
152             //- Allow zetaMin to be changed
153             dimensionedScalar& zetaMin()
154             {
155                 return zetaMin_;
156             }
159         //- Return the turbulence viscosity
160         virtual tmp<volScalarField> nut() const
161         {
162             return nut_;
163         }
165         //- Return the effective diffusivity for q
166         tmp<volScalarField> DqEff() const
167         {
168             return tmp<volScalarField>
169             (
170                 new volScalarField("DqEff", nut_ + nu())
171             );
172         }
174         //- Return the effective diffusivity for epsilon
175         tmp<volScalarField> DzetaEff() const
176         {
177             return tmp<volScalarField>
178             (
179                 new volScalarField("DzetaEff", nut_/sigmaZeta_ + nu())
180             );
181         }
183         //- Return the turbulence kinetic energy
184         virtual tmp<volScalarField> k() const
185         {
186             return k_;
187         }
189         //- Return the turbulence kinetic energy dissipation rate
190         virtual tmp<volScalarField> epsilon() const
191         {
192             return epsilon_;
193         }
195         //- Return the Reynolds stress tensor
196         virtual tmp<volSymmTensorField> R() const;
198         virtual const volScalarField& q() const
199         {
200             return q_;
201         }
203         virtual const volScalarField& zeta() const
204         {
205             return zeta_;
206         }
208         //- Return the effective stress tensor including the laminar stress
209         virtual tmp<volSymmTensorField> devReff() const;
211         //- Return the source term for the momentum equation
212         virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
214         //- Solve the turbulence equations and correct the turbulence viscosity
215         virtual void correct();
217         //- Read RASProperties dictionary
218         virtual bool read();
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 } // End namespace RASModels
225 } // End namespace incompressible
226 } // End namespace Foam
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 #endif
232 // ************************************************************************* //