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-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::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 
88  // Protected Member Functions
89 
90  //- Bound epsilon and return Cmu*sqr(k) for nut
92 
93  //- Correct the eddy-viscosity nut
94  virtual void correctNut();
95 
96 
97 public:
98 
99  //- Runtime type information
100  TypeName("LamBremhorstKE");
101 
102 
103  // Constructors
104 
105  //- Construct from components
107  (
108  const geometricOneField& alpha,
109  const geometricOneField& rho,
110  const volVectorField& U,
111  const surfaceScalarField& alphaRhoPhi,
112  const surfaceScalarField& phi,
113  const viscosity& viscosity,
114  const word& type = typeName
115  );
116 
117  //- Disallow default bitwise copy construction
118  LamBremhorstKE(const LamBremhorstKE&) = delete;
119 
120 
121  //- Destructor
122  virtual ~LamBremhorstKE()
123  {}
124 
125 
126  // Member Functions
127 
128  //- Read RASProperties dictionary
129  virtual bool read();
130 
131  //- Return the effective diffusivity for k
132  tmp<volScalarField> DkEff() const
133  {
134  return volScalarField::New
135  (
136  "DkEff",
137  nut_ + nu()
138  );
139  }
140 
141  //- Return the effective diffusivity for epsilon
143  {
144  return volScalarField::New
145  (
146  "DepsilonEff",
147  nut_/sigmaEps_ + nu()
148  );
149  }
150 
151  //- Return the turbulence kinetic energy
152  virtual tmp<volScalarField> k() const
153  {
154  return k_;
155  }
156 
157  //- Return the turbulence kinetic energy dissipation rate
158  virtual tmp<volScalarField> epsilon() const
159  {
160  return epsilon_;
161  }
162 
163  //- Return the turbulence specific dissipation rate
164  virtual tmp<volScalarField> omega() const
165  {
166  return volScalarField::New
167  (
168  "omega",
169  epsilon_/(Cmu_*k_)
170  );
171  }
172 
173  //- Solve the turbulence equations and correct the turbulence viscosity
174  virtual void correct();
175 
176 
177  // Member Operators
178 
179  //- Disallow default bitwise assignment
180  void operator=(const LamBremhorstKE&) = delete;
181 };
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 } // End namespace RASModels
187 } // End namespace incompressible
188 } // End namespace Foam
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 #endif
193 
194 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
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...
Lam and Bremhorst low-Reynolds number k-epsilon turbulence model for incompressible flows.
void operator=(const LamBremhorstKE &)=delete
Disallow default bitwise assignment.
virtual void correctNut()
Correct the eddy-viscosity nut.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
LamBremhorstKE(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 void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
TypeName("LamBremhorstKE")
Runtime type information.
virtual bool read()
Read RASProperties dictionary.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
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))
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