kkLOmega.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 | Copyright (C) 2011-2015 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::kkLOmega
26 
27 Group
28  grpIcoRASTurbulence
29 
30 Description
31  Low Reynolds-number k-kl-omega turbulence model for
32  incompressible flows.
33 
34  This turbulence model is described in:
35  \verbatim
36  Walters, D. K., & Cokljat, D. (2008).
37  A three-equation eddy-viscosity model for Reynolds-averaged
38  Navier–Stokes simulations of transitional flow.
39  Journal of Fluids Engineering, 130(12), 121401.
40  \endverbatim
41 
42  however the paper contains several errors which must be corrected for the
43  model to operation correctly as explained in
44 
45  \verbatim
46  Furst, J. (2013).
47  Numerical simulation of transitional flows with laminar kinetic energy.
48  Engineering MECHANICS, 20(5), 379-388.
49  \endverbatim
50 
51  All these corrections and updates are included in this implementation.
52 
53  The default model coefficients are
54  \verbatim
55  kkLOmegaCoeffs
56  {
57  A0 4.04
58  As 2.12
59  Av 6.75
60  Abp 0.6
61  Anat 200
62  Ats 200
63  CbpCrit 1.2
64  Cnc 0.1
65  CnatCrit 1250
66  Cint 0.75
67  CtsCrit 1000
68  CrNat 0.02
69  C11 3.4e-6
70  C12 1.0e-10
71  CR 0.12
72  CalphaTheta 0.035
73  Css 1.5
74  CtauL 4360
75  Cw1 0.44
76  Cw2 0.92
77  Cw3 0.3
78  CwR 1.5
79  Clambda 2.495
80  CmuStd 0.09
81  Prtheta 0.85
82  Sigmak 1
83  Sigmaw 1.17
84  }
85  \endverbatim
86 
87 SourceFiles
88  kkLOmega.C
89 
90 \*---------------------------------------------------------------------------*/
91 
92 #ifndef kkLOmega_H
93 #define kkLOmega_H
94 
96 #include "eddyViscosity.H"
97 
98 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
99 
100 namespace Foam
101 {
102 namespace incompressible
103 {
104 namespace RASModels
105 {
106 
107 /*---------------------------------------------------------------------------*\
108  Class kkLOmega Declaration
109 \*---------------------------------------------------------------------------*/
111 class kkLOmega
112 :
113  public eddyViscosity<incompressible::RASModel>
114 {
115  // Private memmber functions
116 
117  tmp<volScalarField> fv(const volScalarField& Ret) const;
118 
119  tmp<volScalarField> fINT() const;
120 
121  tmp<volScalarField> fSS(const volScalarField& omega) const;
122 
123  tmp<volScalarField> Cmu(const volScalarField& S) const;
124 
125  tmp<volScalarField> BetaTS(const volScalarField& Rew) const;
126 
127  tmp<volScalarField> fTaul
128  (
129  const volScalarField& lambdaEff,
130  const volScalarField& ktL,
131  const volScalarField& omega
132  ) const;
133 
134  tmp<volScalarField> alphaT
135  (
136  const volScalarField& lambdaEff,
137  const volScalarField& fv,
138  const volScalarField& ktS
139  ) const;
140 
141  tmp<volScalarField> fOmega
142  (
143  const volScalarField& lambdaEff,
144  const volScalarField& lambdaT
145  ) const;
146 
147  tmp<volScalarField> phiBP(const volScalarField& omega) const;
148 
149  tmp<volScalarField> phiNAT
150  (
151  const volScalarField& ReOmega,
152  const volScalarField& fNatCrit
153  ) const;
154 
155  tmp<volScalarField> D(const volScalarField& k) const;
156 
157 
158 protected:
159 
160  // Protected data
161 
162  // Model coefficients
191 
192 
193  // Fields
199 
200  //- Wall distance
201  // Note: different to wall distance in parent RASModel
202  // which is for near-wall cells only
203  const volScalarField& y_;
204 
205 
206  // Protected Member Functions
207 
208  virtual void correctNut();
209 
210 
211 public:
212 
213  //- Runtime type information
214  TypeName("kkLOmega");
215 
216 
217  // Constructors
218 
219  //- Construct from components
220  kkLOmega
221  (
222  const geometricOneField& alpha,
223  const geometricOneField& rho,
224  const volVectorField& U,
225  const surfaceScalarField& alphaRhoPhi,
226  const surfaceScalarField& phi,
227  const transportModel& transport,
228  const word& propertiesName = turbulenceModel::propertiesName,
229  const word& type = typeName
230  );
231 
232 
233  //- Destructor
234  virtual ~kkLOmega()
235  {}
236 
237 
238  // Member Functions
239 
240  //- Read RASProperties dictionary
241  virtual bool read();
242 
243  //- Return the effective diffusivity for k
244  tmp<volScalarField> DkEff(const volScalarField& alphaT) const
245  {
246  return tmp<volScalarField>
247  (
248  new volScalarField("DkEff", alphaT/Sigmak_ + nu())
249  );
250  }
251 
252  //- Return the effective diffusivity for omega
253  tmp<volScalarField> DomegaEff(const volScalarField& alphaT) const
254  {
255  return tmp<volScalarField>
256  (
257  new volScalarField("DomegaEff", alphaT/Sigmaw_ + nu())
258  );
259  }
260 
261  //- Return the laminar kinetic energy
262  virtual tmp<volScalarField> kl() const
263  {
264  return kl_;
265  }
266 
267  //- Return the turbulence kinetic energy
268  virtual tmp<volScalarField> kt() const
269  {
270  return kt_;
271  }
272 
273  //- Return the turbulence specific dissipation rate
274  virtual tmp<volScalarField> omega() const
275  {
276  return omega_;
277  }
278 
279  //- Return the total fluctuation kinetic energy
280  virtual tmp<volScalarField> k() const
281  {
282  return tmp<volScalarField>
283  (
284  new volScalarField
285  (
286  IOobject
287  (
288  "k",
289  mesh_.time().timeName(),
290  mesh_
291  ),
292  kt_ + kl_,
293  omega_.boundaryField().types()
294  )
295  );
296  }
297 
298  //- Return the total fluctuation kinetic energy dissipation rate
299  virtual tmp<volScalarField> epsilon() const
300  {
301  return epsilon_;
302  }
303 
304  //- Validate the turbulence fields after construction
305  // Update turbulence viscosity and other derived fields as requires
306  virtual void validate();
307 
308  //- Solve the turbulence equations and correct the turbulence viscosity
309  virtual void correct();
310 };
311 
312 
313 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314 
315 } // End namespace RASModels
316 } // End namespace incompressible
317 } // End namespace Foam
318 
319 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
320 
321 #endif
322 
323 // ************************************************************************* //
surfaceScalarField & phi
U
Definition: pEqn.H:83
virtual tmp< volScalarField > kt() const
Return the turbulence kinetic energy.
Definition: kkLOmega.H:267
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kkLOmega.C:600
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
incompressible::RASModel::transportModel transportModel
Definition: eddyViscosity.H:75
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Low Reynolds-number k-kl-omega turbulence model for incompressible flows.
Definition: kkLOmega.H:110
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read()
Read RASProperties dictionary.
Definition: kkLOmega.C:554
kkLOmega(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: kkLOmega.C:223
TypeName("kkLOmega")
Runtime type information.
virtual void validate()
Validate the turbulence fields after construction.
Definition: kkLOmega.C:596
virtual tmp< volScalarField > kl() const
Return the laminar kinetic energy.
Definition: kkLOmega.H:261
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
virtual tmp< volScalarField > k() const
Return the total fluctuation kinetic energy.
Definition: kkLOmega.H:279
virtual tmp< volScalarField > epsilon() const
Return the total fluctuation kinetic energy dissipation rate.
Definition: kkLOmega.H:298
tmp< volScalarField > DkEff(const volScalarField &alphaT) const
Return the effective diffusivity for k.
Definition: kkLOmega.H:243
const volScalarField & y_
Wall distance.
Definition: kkLOmega.H:202
A class for managing temporary objects.
Definition: PtrList.H:54
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
tmp< volScalarField > DomegaEff(const volScalarField &alphaT) const
Return the effective diffusivity for omega.
Definition: kkLOmega.H:252
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
volScalarField & nu
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
Definition: kkLOmega.H:273
Namespace for OpenFOAM.