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-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::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  tmp<volScalarField> Rt() const;
69  tmp<volScalarField> fMu(const volScalarField& Rt) const;
70  tmp<volScalarField> f1(const volScalarField& fMu) const;
71  tmp<volScalarField> f2(const volScalarField& Rt) const;
72  void correctNut(const volScalarField& fMu);
73 
74 
75 protected:
76 
77  // Protected data
78 
83 
86 
87  //- Wall distance
88  // Note: different to wall distance in parent RASModel
89  // which is for near-wall cells only
90  const volScalarField& y_;
91 
92 
93  // Protected Member Functions
94 
95  virtual void correctNut();
96 
97 
98 public:
99 
100  //- Runtime type information
101  TypeName("LamBremhorstKE");
102 
103 
104  // Constructors
105 
106  //- Construct from components
108  (
109  const geometricOneField& alpha,
110  const geometricOneField& rho,
111  const volVectorField& U,
112  const surfaceScalarField& alphaRhoPhi,
113  const surfaceScalarField& phi,
114  const transportModel& transport,
115  const word& type = typeName
116  );
117 
118  //- Disallow default bitwise copy construction
119  LamBremhorstKE(const LamBremhorstKE&) = delete;
120 
121 
122  //- Destructor
123  virtual ~LamBremhorstKE()
124  {}
125 
126 
127  // Member Functions
128 
129  //- Read RASProperties dictionary
130  virtual bool read();
131 
132  //- Return the effective diffusivity for k
133  tmp<volScalarField> DkEff() const
134  {
135  return volScalarField::New
136  (
137  "DkEff",
138  nut_ + nu()
139  );
140  }
141 
142  //- Return the effective diffusivity for epsilon
144  {
145  return volScalarField::New
146  (
147  "DepsilonEff",
148  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  // Member Operators
169 
170  //- Disallow default bitwise assignment
171  void operator=(const LamBremhorstKE&) = delete;
172 };
173 
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 } // End namespace RASModels
178 } // End namespace incompressible
179 } // End namespace Foam
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 #endif
184 
185 // ************************************************************************* //
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual bool read()
Read RASProperties dictionary.
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...
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
TypeName("LamBremhorstKE")
Runtime type information.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
incompressible::RASModel ::transportModel transportModel
Definition: eddyViscosity.H:72
A class for handling words, derived from string.
Definition: word.H:59
phi
Definition: correctPhi.H:3
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const volScalarField & y_
Wall distance.
U
Definition: pEqn.H:72
void operator=(const LamBremhorstKE &)=delete
Disallow default bitwise assignment.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
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
Namespace for OpenFOAM.
LamBremhorstKE(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.