qZeta.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Class
25  Foam::incompressible::RASModels::qZeta
26 
27 Description
28  Gibson and Dafa'Alla's q-zeta two-equation low-Re turbulence model
29  for incompressible flows
30 
31  This turbulence model is described in:
32  \verbatim
33  Dafa'Alla, A.A., Juntasaro, E. & Gibson, M.M. (1996).
34  Calculation of oscillating boundary layers with the
35  q-zeta turbulence model.
36  Engineering Turbulence Modelling and Experiments 3:
37  Proceedings of the Third International Symposium,
38  Crete, Greece, May 27-29, 141.
39  \endverbatim
40  which is a development of the original q-zeta model described in:
41  \verbatim
42  Gibson, M. M., & Dafa'Alla, A. A. (1995).
43  Two-equation model for turbulent wall flow.
44  AIAA journal, 33(8), 1514-1518.
45  \endverbatim
46 
47 SourceFiles
48  qZeta.C
49 
50 \*---------------------------------------------------------------------------*/
51 
52 #ifndef qZeta_H
53 #define qZeta_H
54 
56 #include "eddyViscosity.H"
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 namespace Foam
61 {
62 namespace incompressible
63 {
64 namespace RASModels
65 {
66 
67 /*---------------------------------------------------------------------------*\
68  Class qZeta Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 class qZeta
72 :
73  public eddyViscosity<incompressible::RASModel>
74 {
75 
76 protected:
77 
78  // Protected data
79 
80  // Model coefficients
81 
87 
88  //- Lower limit of q
90 
91  // Fields
92 
95 
98 
99 
100  // Protected Member Functions
101 
103  tmp<volScalarField> f2() const;
104 
105  //- Bound epsilon and return Cmu*sqr(k) for nut
106  void boundZeta();
107 
108  //- Correct the eddy-viscosity nut
109  virtual void correctNut();
110 
111 
112 public:
113 
114  //- Runtime type information
115  TypeName("qZeta");
116 
117  // Constructors
118 
119  //- Construct from components
120  qZeta
121  (
122  const geometricOneField& alpha,
123  const geometricOneField& rho,
124  const volVectorField& U,
125  const surfaceScalarField& alphaRhoPhi,
126  const surfaceScalarField& phi,
127  const viscosity& viscosity,
128  const word& type = typeName
129  );
130 
131 
132  //- Destructor
133  virtual ~qZeta()
134  {}
135 
136 
137  // Member Functions
138 
139  //- Read RASProperties dictionary
140  virtual bool read();
141 
142  //- Return the effective diffusivity for q
143  tmp<volScalarField> DqEff() const
144  {
145  return volScalarField::New
146  (
147  "DqEff",
148  nut_ + nu()
149  );
150  }
151 
152  //- Return the effective diffusivity for epsilon
154  {
155  return volScalarField::New
156  (
157  "DzetaEff",
158  nut_/sigmaZeta_ + nu()
159  );
160  }
161 
162  //- Return the turbulence kinetic energy
163  virtual tmp<volScalarField> k() const
164  {
165  return k_;
166  }
167 
168  //- Return the turbulence kinetic energy dissipation rate
169  virtual tmp<volScalarField> epsilon() const
170  {
171  return epsilon_;
172  }
173 
174  //- Return the turbulence specific dissipation rate
175  virtual tmp<volScalarField> omega() const
176  {
177  return volScalarField::New
178  (
179  "omega",
180  epsilon_/(Cmu_*k_)
181  );
182  }
183 
184  virtual const volScalarField& q() const
185  {
186  return q_;
187  }
188 
189  virtual const volScalarField& zeta() const
190  {
191  return zeta_;
192  }
193 
194  //- Solve the turbulence equations and correct the turbulence viscosity
195  virtual void correct();
196 };
197 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 } // End namespace RASModels
202 } // End namespace incompressible
203 } // End namespace Foam
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 #endif
208 
209 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
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:73
virtual void correctNut()
Correct the eddy-viscosity nut.
tmp< volScalarField > f2() const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: qZeta.H:162
qZeta(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
Definition: qZeta.H:174
void boundZeta()
Bound epsilon and return Cmu*sqr(k) for nut.
TypeName("qZeta")
Runtime type information.
tmp< volScalarField > DqEff() const
Return the effective diffusivity for q.
Definition: qZeta.H:142
dimensionedScalar sigmaZeta_
Definition: qZeta.H:84
dimensionedScalar qMin_
Lower limit of q.
Definition: qZeta.H:88
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual const volScalarField & q() const
Definition: qZeta.H:183
virtual bool read()
Read RASProperties dictionary.
virtual const volScalarField & zeta() const
Definition: qZeta.H:188
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: qZeta.H:168
tmp< volScalarField > fMu() const
tmp< volScalarField > DzetaEff() const
Return the effective diffusivity for epsilon.
Definition: qZeta.H:152
virtual ~qZeta()
Destructor.
Definition: qZeta.H:132
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488