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