qZeta.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2015 OpenFOAM Foundation
9 Copyright (C) 2019-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::incompressible::RASModels::qZeta
29
30Group
31 grpIcoRASTurbulence
32
33Description
34 Gibson and Dafa'Alla's q-zeta two-equation low-Re turbulence model
35 for incompressible flows
36
37 This turbulence model is described in:
38 \verbatim
39 Dafa'Alla, A.A., Juntasaro, E. & Gibson, M.M. (1996).
40 Calculation of oscillating boundary layers with the
41 q-zeta turbulence model.
42 Engineering Turbulence Modelling and Experiments 3:
43 Proceedings of the Third International Symposium,
44 Crete, Greece, May 27-29, 141.
45 \endverbatim
46 which is a development of the original q-zeta model described in:
47 \verbatim
48 Gibson, M. M., & Dafa'Alla, A. A. (1995).
49 Two-equation model for turbulent wall flow.
50 AIAA journal, 33(8), 1514-1518.
51 \endverbatim
52
53SourceFiles
54 qZeta.C
55
56\*---------------------------------------------------------------------------*/
57
58#ifndef qZeta_H
59#define qZeta_H
60
62#include "eddyViscosity.H"
63
64// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65
66namespace Foam
67{
68namespace incompressible
69{
70namespace RASModels
71{
72
73/*---------------------------------------------------------------------------*\
74 Class qZeta Declaration
75\*---------------------------------------------------------------------------*/
77class qZeta
78:
79 public eddyViscosity<incompressible::RASModel>
80{
81
82protected:
83
84 // Protected data
85
86 // Model coefficients
93
94 //- Lower limit of q
96
97 //- Lower limit of zeta
99
100 // Fields
107
108
109 // Protected Member Functions
110
111 tmp<volScalarField> fMu() const;
112 tmp<volScalarField> f2() const;
113 virtual void correctNut();
114
115
116public:
117
118 //- Runtime type information
119 TypeName("qZeta");
120
121 // Constructors
122
123 //- Construct from components
124 qZeta
125 (
127 const geometricOneField& rho,
128 const volVectorField& U,
129 const surfaceScalarField& alphaRhoPhi,
130 const surfaceScalarField& phi,
131 const transportModel& transport,
132 const word& propertiesName = turbulenceModel::propertiesName,
133 const word& type = typeName
134 );
135
136
137 //- Destructor
138 virtual ~qZeta() = default;
139
140
141 // Member Functions
142
143 //- Re-read model coefficients if they have changed
144 virtual bool read();
145
146 //- Return the lower allowable limit for q (default: SMALL)
147 const dimensionedScalar& qMin() const
148 {
149 return qMin_;
150 }
151
152 //- Return the lower allowable limit for zeta (default: SMALL)
153 const dimensionedScalar& zetaMin() const
154 {
155 return zetaMin_;
156 }
157
158 //- Allow qMin to be changed
160 {
161 return qMin_;
162 }
163
164 //- Allow zetaMin to be changed
166 {
167 return zetaMin_;
168 }
169
170 //- Return the effective diffusivity for q
172 {
174 (
175 new volScalarField("DqEff", nut_ + nu())
176 );
177 }
178
179 //- Return the effective diffusivity for epsilon
181 {
183 (
184 new volScalarField("DzetaEff", nut_/sigmaZeta_ + nu())
185 );
186 }
187
188 //- Return the turbulence kinetic energy
189 virtual tmp<volScalarField> k() const
190 {
191 return k_;
192 }
193
194 //- Return the turbulence kinetic energy dissipation rate
195 virtual tmp<volScalarField> epsilon() const
196 {
197 return epsilon_;
198 }
200 virtual const volScalarField& q() const
201 {
202 return q_;
203 }
205 virtual const volScalarField& zeta() const
206 {
207 return zeta_;
208 }
209
210 //- Solve the turbulence equations and correct the turbulence viscosity
211 virtual void correct();
212};
213
214
215// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216
217} // End namespace RASModels
218} // End namespace incompressible
219} // End namespace Foam
220
221// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222
223#endif
224
225// ************************************************************************* //
surfaceScalarField & phi
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:58
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Gibson and Dafa'Alla's q-zeta two-equation low-Re turbulence model for incompressible flows.
Definition: qZeta.H:79
const dimensionedScalar & qMin() const
Return the lower allowable limit for q (default: SMALL)
Definition: qZeta.H:146
tmp< volScalarField > f2() const
Definition: qZeta.C:66
virtual ~qZeta()=default
Destructor.
dimensionedScalar & zetaMin()
Allow zetaMin to be changed.
Definition: qZeta.H:164
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: qZeta.C:240
TypeName("qZeta")
Runtime type information.
dimensionedScalar sigmaZeta_
Definition: qZeta.H:90
dimensionedScalar & qMin()
Allow qMin to be changed.
Definition: qZeta.H:158
dimensionedScalar qMin_
Lower limit of q.
Definition: qZeta.H:94
tmp< volScalarField > DqEff() const
Return the effective diffusivity for q.
Definition: qZeta.H:170
const dimensionedScalar & zetaMin() const
Return the lower allowable limit for zeta (default: SMALL)
Definition: qZeta.H:152
tmp< volScalarField > DzetaEff() const
Return the effective diffusivity for epsilon.
Definition: qZeta.H:179
virtual const volScalarField & zeta() const
Definition: qZeta.H:204
dimensionedScalar zetaMin_
Lower limit of zeta.
Definition: qZeta.H:97
virtual const volScalarField & q() const
Definition: qZeta.H:199
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: qZeta.H:188
tmp< volScalarField > fMu() const
Definition: qZeta.C:49
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: qZeta.H:194
virtual bool read()
Re-read model coefficients if they have changed.
Definition: qZeta.C:220
BasicTurbulenceModel::transportModel transportModel
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:290
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
volScalarField & nu
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73