LamBremhorstKE.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-2021 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::LamBremhorstKE
29 
30 Group
31  grpIcoRASTurbulence
32 
33 Description
34  Lam and Bremhorst low-Reynolds number k-epsilon turbulence model
35  for incompressible flows
36 
37  This turbulence model is described in:
38  \verbatim
39  Lam, C. K. G., & Bremhorst, K. (1981).
40  A modified form of the k-ε model for predicting wall turbulence.
41  Journal of Fluids Engineering, 103(3), 456-460.
42  \endverbatim
43 
44 SourceFiles
45  LamBremhorstKE.C
46 
47 \*---------------------------------------------------------------------------*/
48 
49 #ifndef LamBremhorstKE_H
50 #define LamBremhorstKE_H
51 
53 #include "eddyViscosity.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
59 namespace incompressible
60 {
61 namespace RASModels
62 {
63 
64 /*---------------------------------------------------------------------------*\
65  Class LamBremhorstKE Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 class LamBremhorstKE
69 :
70  public eddyViscosity<incompressible::RASModel>
71 {
72  // Private Member Functions
73 
74  //- No copy construct
75  LamBremhorstKE(const LamBremhorstKE&) = delete;
76 
77  //- No copy assignment
78  void operator=(const LamBremhorstKE&) = delete;
79 
80  tmp<volScalarField> Rt() const;
81  tmp<volScalarField> fMu(const volScalarField& Rt) const;
82  tmp<volScalarField> f1(const volScalarField& fMu) const;
83  tmp<volScalarField> f2(const volScalarField& Rt) const;
84  void correctNut(const volScalarField& fMu);
85 
86 
87 protected:
88 
89  // Protected data
90 
95 
98 
99  //- Wall distance
100  // Note: different to wall distance in parent RASModel
101  // which is for near-wall cells only
102  const volScalarField& y_;
103 
105  // Protected Member Functions
106 
107  virtual void correctNut();
108 
109 
110 public:
111 
112  //- Runtime type information
113  TypeName("LamBremhorstKE");
114 
115 
116  // Constructors
117 
118  //- Construct from components
120  (
121  const geometricOneField& alpha,
122  const geometricOneField& rho,
123  const volVectorField& U,
124  const surfaceScalarField& alphaRhoPhi,
125  const surfaceScalarField& phi,
126  const transportModel& transport,
127  const word& propertiesName = turbulenceModel::propertiesName,
128  const word& type = typeName
129  );
130 
131 
132  //- Destructor
133  virtual ~LamBremhorstKE() = default;
134 
135 
136  // Member Functions
137 
138  //- Re-read model coefficients if they have changed
139  virtual bool read();
140 
141  //- Return the effective diffusivity for k
142  tmp<volScalarField> DkEff() const
143  {
144  return tmp<volScalarField>
145  (
146  new volScalarField("DkEff", nut_ + nu())
147  );
148  }
149 
150  //- Return the effective diffusivity for epsilon
152  {
153  return tmp<volScalarField>
154  (
155  new volScalarField("DepsilonEff", nut_/sigmaEps_ + nu())
156  );
157  }
158 
159  //- Return the turbulence kinetic energy
160  virtual tmp<volScalarField> k() const
161  {
162  return k_;
163  }
164 
165  //- Return the turbulence kinetic energy dissipation rate
166  virtual tmp<volScalarField> epsilon() const
167  {
168  return epsilon_;
169  }
170 
171  //- Solve the turbulence equations and correct the turbulence viscosity
172  virtual void correct();
173 };
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 } // End namespace RASModels
179 } // End namespace incompressible
180 } // End namespace Foam
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 #endif
185 
186 // ************************************************************************* //
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual bool read()
Re-read model coefficients if they have changed.
virtual ~LamBremhorstKE()=default
Destructor.
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:51
TypeName("LamBremhorstKE")
Runtime type information.
incompressible::RASModel ::transportModel transportModel
Definition: eddyViscosity.H:74
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
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
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 Foam::string.
Definition: word.H:63
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const volScalarField & y_
Wall distance.
U
Definition: pEqn.H:72
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: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & nu
Namespace for OpenFOAM.