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-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::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  corrected according to:
40 
41  \verbatim
42  Furst, J. (2013).
43  Numerical simulation of transitional flows with laminar kinetic energy.
44  Engineering Mechanics, 20(5), 379-388.
45  \endverbatim
46 
47  and includes the improvements proposed in:
48 
49  \verbatim
50  Lopez, M., and Keith Walters, D. (2016).
51  A Recommended Correction to the kT−kL−ω Transition-Sensitive
52  Eddy-Viscosity Model.
53  Journal of Fluids Engineering, 139(2), 024501.
54  \endverbatim
55 
56  The default model coefficients are
57  \verbatim
58  kkLOmegaCoeffs
59  {
60  A0 4.04
61  As 2.12
62  Av 6.75
63  Abp 0.6
64  Anat 200
65  Ats 200
66  CbpCrit 1.2
67  Cnc 0.1
68  CnatCrit 1250
69  Cint 0.75
70  CtsCrit 1000
71  CrNat 0.02
72  C11 3.4e-6
73  C12 1.0e-10
74  CR 0.12
75  CalphaTheta 0.035
76  Css 1.5
77  CtauL 4360
78  Cw1 0.44
79  Cw2 0.92
80  Cw3 0.3
81  CwR 1.5
82  Clambda 2.495
83  CmuStd 0.09
84  Prtheta 0.85
85  Sigmak 1
86  Sigmaw 1.17
87  }
88  \endverbatim
89 
90 SourceFiles
91  kkLOmega.C
92 
93 \*---------------------------------------------------------------------------*/
94 
95 #ifndef kkLOmega_H
96 #define kkLOmega_H
97 
99 #include "eddyViscosity.H"
100 
101 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
102 
103 namespace Foam
104 {
105 namespace incompressible
106 {
107 namespace RASModels
108 {
109 
110 /*---------------------------------------------------------------------------*\
111  Class kkLOmega Declaration
112 \*---------------------------------------------------------------------------*/
113 
114 class kkLOmega
115 :
116  public eddyViscosity<incompressible::RASModel>
117 {
118  // Private Member Functions
119 
120  tmp<volScalarField> fv(const volScalarField& Ret) const;
121 
122  tmp<volScalarField> fINT() const;
123 
124  tmp<volScalarField> fSS(const volScalarField& omega) const;
125 
126  tmp<volScalarField> Cmu(const volScalarField& S) const;
127 
128  tmp<volScalarField> BetaTS(const volScalarField& Rew) const;
129 
130  tmp<volScalarField> fTaul
131  (
132  const volScalarField& lambdaEff,
133  const volScalarField& ktL,
134  const volScalarField& omega
135  ) const;
136 
137  tmp<volScalarField> alphaT
138  (
139  const volScalarField& lambdaEff,
140  const volScalarField& fv,
141  const volScalarField& ktS
142  ) const;
143 
144  tmp<volScalarField> fOmega
145  (
146  const volScalarField& lambdaEff,
147  const volScalarField& lambdaT
148  ) const;
149 
150  tmp<volScalarField> phiBP(const volScalarField& omega) const;
151 
152  tmp<volScalarField> phiNAT
153  (
154  const volScalarField& ReOmega,
155  const volScalarField& fNatCrit
156  ) const;
157 
158  tmp<volScalarField> D(const volScalarField& k) const;
159 
160 
161 protected:
162 
163  // Protected data
164 
165  // Model coefficients
166 
194 
195 
196  // Fields
197 
202 
203  //- Wall distance
204  // Note: different to wall distance in parent RASModel
205  // which is for near-wall cells only
206  const volScalarField& y_;
207 
208 
209  // Protected Member Functions
210 
211  virtual void correctNut();
212 
213 
214 public:
215 
216  //- Runtime type information
217  TypeName("kkLOmega");
218 
219 
220  // Constructors
221 
222  //- Construct from components
223  kkLOmega
224  (
225  const geometricOneField& alpha,
226  const geometricOneField& rho,
227  const volVectorField& U,
228  const surfaceScalarField& alphaRhoPhi,
229  const surfaceScalarField& phi,
230  const viscosity& viscosity,
231  const word& type = typeName
232  );
233 
234 
235  //- Destructor
236  virtual ~kkLOmega()
237  {}
238 
239 
240  // Member Functions
241 
242  //- Read RASProperties dictionary
243  virtual bool read();
244 
245  //- Return the effective diffusivity for k
246  tmp<volScalarField> DkEff(const volScalarField& alphaT) const
247  {
248  return volScalarField::New
249  (
250  "DkEff",
251  alphaT/Sigmak_ + nu()
252  );
253  }
254 
255  //- Return the effective diffusivity for omega
256  tmp<volScalarField> DomegaEff(const volScalarField& alphaT) const
257  {
258  return volScalarField::New
259  (
260  "DomegaEff",
261  alphaT/Sigmaw_ + nu()
262  );
263  }
264 
265  //- Return the laminar kinetic energy
266  virtual tmp<volScalarField> kl() const
267  {
268  return kl_;
269  }
270 
271  //- Return the turbulence kinetic energy
272  virtual tmp<volScalarField> kt() const
273  {
274  return kt_;
275  }
276 
277  //- Return the turbulence specific dissipation rate
278  virtual tmp<volScalarField> omega() const
279  {
280  return omega_;
281  }
282 
283  //- Return the total fluctuation kinetic energy
284  virtual tmp<volScalarField> k() const
285  {
286  return volScalarField::New
287  (
288  "k",
289  kt_ + kl_
290  );
291  }
292 
293  //- Return the total fluctuation kinetic energy dissipation rate
294  virtual tmp<volScalarField> epsilon() const
295  {
296  return epsilon_;
297  }
298 
299  //- Validate the turbulence fields after construction
300  // Update turbulence viscosity and other derived fields as requires
301  virtual void validate();
302 
303  //- Solve the turbulence equations and correct the turbulence viscosity
304  virtual void correct();
305 };
306 
307 
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 
310 } // End namespace RASModels
311 } // End namespace incompressible
312 } // End namespace Foam
313 
314 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315 
316 #endif
317 
318 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
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...
Low Reynolds-number k-kl-omega turbulence model for incompressible flows.
Definition: kkLOmega.H:116
virtual void validate()
Validate the turbulence fields after construction.
TypeName("kkLOmega")
Runtime type information.
virtual tmp< volScalarField > k() const
Return the total fluctuation kinetic energy.
Definition: kkLOmega.H:283
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
Definition: kkLOmega.H:277
virtual tmp< volScalarField > kt() const
Return the turbulence kinetic energy.
Definition: kkLOmega.H:271
const volScalarField & y_
Wall distance.
Definition: kkLOmega.H:205
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
kkLOmega(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 bool read()
Read RASProperties dictionary.
virtual tmp< volScalarField > kl() const
Return the laminar kinetic energy.
Definition: kkLOmega.H:265
virtual tmp< volScalarField > epsilon() const
Return the total fluctuation kinetic energy dissipation rate.
Definition: kkLOmega.H:293
tmp< volScalarField > DomegaEff(const volScalarField &alphaT) const
Return the effective diffusivity for omega.
Definition: kkLOmega.H:255
tmp< volScalarField > DkEff(const volScalarField &alphaT) const
Return the effective diffusivity for k.
Definition: kkLOmega.H:245
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))
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
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
labelList fv(nPoints)