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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::incompressible::RASModels::kkLOmega
29 
30 Group
31  grpIcoRASTurbulence
32 
33 Description
34  Low Reynolds-number k-kl-omega turbulence model for
35  incompressible flows.
36 
37  Reference:
38  \verbatim
39  Standard model:
40  Walters, D. K., & Cokljat, D. (2008).
41  A three-equation eddy-viscosity model for Reynolds-averaged
42  Navier–Stokes simulations of transitional flow.
43  Journal of Fluids Engineering, 130(12), 121401.
44  DOI:10.1115/1.2979230
45 
46  Typo corrections:
47  Furst, J. (2013).
48  Numerical simulation of transitional flows
49  with laminar kinetic energy.
50  Engineering MECHANICS, 20(5), 379-388.
51  \endverbatim
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 {
107 /*---------------------------------------------------------------------------*\
108  Class kkLOmega Declaration
109 \*---------------------------------------------------------------------------*/
110 
111 class kkLOmega
112 :
113  public eddyViscosity<incompressible::RASModel>
114 {
115  // Private member 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:
160  // Protected data
162  // Model coefficients
193  // Fields
194 
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() = default;
235 
236 
237  // Member Functions
238 
239  //- Re-read model coefficients if they have changed
240  virtual bool read();
241 
242  //- Return the effective diffusivity for k
243  tmp<volScalarField> DkEff(const volScalarField& alphaT) const
244  {
245  return tmp<volScalarField>
246  (
247  new volScalarField("DkEff", alphaT/Sigmak_ + nu())
248  );
249  }
250 
251  //- Return the effective diffusivity for omega
252  tmp<volScalarField> DomegaEff(const volScalarField& alphaT) const
253  {
254  return tmp<volScalarField>
255  (
256  new volScalarField("DomegaEff", alphaT/Sigmaw_ + nu())
257  );
258  }
259 
260  //- Return the laminar kinetic energy
261  virtual tmp<volScalarField> kl() const
262  {
263  return kl_;
264  }
265 
266  //- Return the turbulence kinetic energy
267  virtual tmp<volScalarField> kt() const
268  {
269  return kt_;
270  }
271 
272  //- Return the turbulence specific dissipation rate
273  virtual tmp<volScalarField> omega() const
274  {
275  return omega_;
276  }
277 
278  //- Return the total fluctuation kinetic energy
279  virtual tmp<volScalarField> k() const
280  {
282  (
283  new volScalarField
284  (
285  IOobject
286  (
287  "k",
288  mesh_.time().timeName(),
289  mesh_
290  ),
291  kt_ + kl_,
293  )
294  );
295  }
296 
297  //- Return the total fluctuation kinetic energy dissipation rate
298  virtual tmp<volScalarField> epsilon() const
299  {
300  return epsilon_;
301  }
302 
303  //- Validate the turbulence fields after construction
304  // Update turbulence viscosity and other derived fields as requires
305  virtual void validate();
306 
307  //- Solve the turbulence equations and correct the turbulence viscosity
308  virtual void correct();
309 };
310 
311 
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313 
314 } // End namespace RASModels
315 } // End namespace incompressible
316 } // End namespace Foam
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 #endif
321 
322 // ************************************************************************* //
virtual tmp< volScalarField > epsilon() const
Return the total fluctuation kinetic energy dissipation rate.
Definition: kkLOmega.H:318
wordList types() const
Return a list of the patch types.
virtual ~kkLOmega()=default
Destructor.
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:289
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:51
incompressible::RASModel ::transportModel transportModel
Definition: eddyViscosity.H:74
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
virtual tmp< volScalarField > kl() const
Return the laminar kinetic energy.
Definition: kkLOmega.H:273
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual tmp< volScalarField > k() const
Return the total fluctuation kinetic energy.
Definition: kkLOmega.H:297
Low Reynolds-number k-kl-omega turbulence model for incompressible flows.
Definition: kkLOmega.H:106
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:251
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kkLOmega.C:550
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:219
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:281
U
Definition: pEqn.H:72
const volScalarField & y_
Wall distance.
Definition: kkLOmega.H:201
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< volScalarField > DomegaEff(const volScalarField &alphaT) const
Return the effective diffusivity for omega.
Definition: kkLOmega.H:262
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
volScalarField & nu
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Namespace for OpenFOAM.