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-2020 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& type = typeName
226  );
227 
228 
229  //- Destructor
230  virtual ~kkLOmega()
231  {}
232 
233 
234  // Member Functions
235 
236  //- Read RASProperties dictionary
237  virtual bool read();
238 
239  //- Return the effective diffusivity for k
240  tmp<volScalarField> DkEff(const volScalarField& alphaT) const
241  {
242  return volScalarField::New
243  (
244  "DkEff",
245  alphaT/Sigmak_ + nu()
246  );
247  }
248 
249  //- Return the effective diffusivity for omega
250  tmp<volScalarField> DomegaEff(const volScalarField& alphaT) const
251  {
252  return volScalarField::New
253  (
254  "DomegaEff",
255  alphaT/Sigmaw_ + nu()
256  );
257  }
258 
259  //- Return the laminar kinetic energy
260  virtual tmp<volScalarField> kl() const
261  {
262  return kl_;
263  }
264 
265  //- Return the turbulence kinetic energy
266  virtual tmp<volScalarField> kt() const
267  {
268  return kt_;
269  }
270 
271  //- Return the turbulence specific dissipation rate
272  virtual tmp<volScalarField> omega() const
273  {
274  return omega_;
275  }
276 
277  //- Return the total fluctuation kinetic energy
278  virtual tmp<volScalarField> k() const
279  {
280  return volScalarField::New
281  (
282  "k",
283  kt_ + kl_,
284  omega_.boundaryField().types()
285  );
286  }
287 
288  //- Return the total fluctuation kinetic energy dissipation rate
289  virtual tmp<volScalarField> epsilon() const
290  {
291  return epsilon_;
292  }
293 
294  //- Validate the turbulence fields after construction
295  // Update turbulence viscosity and other derived fields as requires
296  virtual void validate();
297 
298  //- Solve the turbulence equations and correct the turbulence viscosity
299  virtual void correct();
300 };
301 
302 
303 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
304 
305 } // End namespace RASModels
306 } // End namespace incompressible
307 } // End namespace Foam
308 
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 
311 #endif
312 
313 // ************************************************************************* //
virtual tmp< volScalarField > epsilon() const
Return the total fluctuation kinetic energy dissipation rate.
Definition: kkLOmega.H:288
const Boundary & boundaryField() const
Return const-reference to the boundary field.
kkLOmega(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: kkLOmega.C:219
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kkLOmega.C:594
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
Definition: kkLOmega.H:271
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
phi
Definition: pEqn.H:104
virtual tmp< volScalarField > kl() const
Return the laminar kinetic energy.
Definition: kkLOmega.H:259
virtual tmp< volScalarField > k() const
Return the total fluctuation kinetic energy.
Definition: kkLOmega.H:277
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:239
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read()
Read RASProperties dictionary.
Definition: kkLOmega.C:548
TypeName("kkLOmega")
Runtime type information.
virtual void validate()
Validate the turbulence fields after construction.
Definition: kkLOmega.C:590
virtual tmp< volScalarField > kt() const
Return the turbulence kinetic energy.
Definition: kkLOmega.H:265
U
Definition: pEqn.H:72
Base-class for all transport models used by the incompressible turbulence models. ...
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const volScalarField & y_
Wall distance.
Definition: kkLOmega.H:199
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
Namespace for OpenFOAM.