LamBremhorstKE.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::LamBremhorstKE
26 
27 Description
28  Lam and Bremhorst low-Reynolds number k-epsilon turbulence model
29  for incompressible flows
30 
31  This turbulence model is described in:
32  \verbatim
33  Lam, C. K. G., & Bremhorst, K. (1981).
34  A modified form of the k-ε model for predicting wall turbulence.
35  Journal of Fluids Engineering, 103(3), 456-460.
36  \endverbatim
37 
38 SourceFiles
39  LamBremhorstKE.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef LamBremhorstKE_H
44 #define LamBremhorstKE_H
45 
47 #include "eddyViscosity.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace incompressible
54 {
55 namespace RASModels
56 {
57 
58 /*---------------------------------------------------------------------------*\
59  Class LamBremhorstKE Declaration
60 \*---------------------------------------------------------------------------*/
61 
62 class LamBremhorstKE
63 :
64  public eddyViscosity<incompressible::RASModel>
65 {
66  // Private Member Functions
67 
68  // Disallow default bitwise copy construct and assignment
70  void operator=(const LamBremhorstKE&);
71 
72  tmp<volScalarField> Rt() const;
73  tmp<volScalarField> fMu(const volScalarField& Rt) const;
74  tmp<volScalarField> f1(const volScalarField& fMu) const;
75  tmp<volScalarField> f2(const volScalarField& Rt) const;
76  void correctNut(const volScalarField& fMu);
77 
78 
79 protected:
80 
81  // Protected data
82 
87 
90 
91  //- Wall distance
92  // Note: different to wall distance in parent RASModel
93  // which is for near-wall cells only
94  const volScalarField& y_;
95 
96 
97  // Protected Member Functions
98 
99  virtual void correctNut();
100 
101 
102 public:
103 
104  //- Runtime type information
105  TypeName("LamBremhorstKE");
106 
107 
108  // Constructors
109 
110  //- Construct from components
112  (
113  const geometricOneField& alpha,
114  const geometricOneField& rho,
115  const volVectorField& U,
116  const surfaceScalarField& alphaRhoPhi,
117  const surfaceScalarField& phi,
118  const transportModel& transport,
119  const word& propertiesName = turbulenceModel::propertiesName,
120  const word& type = typeName
121  );
122 
123 
124  //- Destructor
125  virtual ~LamBremhorstKE()
126  {}
127 
128 
129  // Member Functions
130 
131  //- Read RASProperties dictionary
132  virtual bool read();
133 
134  //- Return the effective diffusivity for k
135  tmp<volScalarField> DkEff() const
136  {
137  return tmp<volScalarField>
138  (
139  new volScalarField("DkEff", nut_ + nu())
140  );
141  }
142 
143  //- Return the effective diffusivity for epsilon
145  {
146  return tmp<volScalarField>
147  (
148  new volScalarField("DepsilonEff", nut_/sigmaEps_ + nu())
149  );
150  }
151 
152  //- Return the turbulence kinetic energy
153  virtual tmp<volScalarField> k() const
154  {
155  return k_;
156  }
157 
158  //- Return the turbulence kinetic energy dissipation rate
159  virtual tmp<volScalarField> epsilon() const
160  {
161  return epsilon_;
162  }
163 
164  //- Solve the turbulence equations and correct the turbulence viscosity
165  virtual void correct();
166 };
167 
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
171 } // End namespace RASModels
172 } // End namespace incompressible
173 } // End namespace Foam
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 #endif
178 
179 // ************************************************************************* //
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual bool read()
Read RASProperties dictionary.
surfaceScalarField & phi
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Lam and Bremhorst low-Reynolds number k-epsilon turbulence model for incompressible flows...
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
TypeName("LamBremhorstKE")
Runtime type information.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const volScalarField & y_
Wall distance.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
U
Definition: pEqn.H:72
Base-class for all transport models used by the incompressible turbulence models. ...
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & nu
Namespace for OpenFOAM.