kEpsilon.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::RASModels::kEpsilon
26 
27 Description
28  Standard k-epsilon turbulence model for incompressible and compressible
29  flows including rapid distortion theory (RDT) based compression term.
30 
31  Reference:
32  \verbatim
33  Standard model:
34  Launder, B. E., & Spalding, D. B. (1972).
35  Lectures in mathematical models of turbulence.
36 
37  Launder, B. E., & Spalding, D. B. (1974).
38  The numerical computation of turbulent flows.
39  Computer methods in applied mechanics and engineering,
40  3(2), 269-289.
41 
42  For the RDT-based compression term:
43  El Tahry, S. H. (1983).
44  k-epsilon equation for compressible reciprocating engine flows.
45  Journal of Energy, 7(4), 345-353.
46  \endverbatim
47 
48  The default model coefficients are
49  \verbatim
50  kEpsilonCoeffs
51  {
52  Cmu 0.09;
53  C1 1.44;
54  C2 1.92;
55  C3 0;
56  sigmak 1.0;
57  sigmaEps 1.3;
58  }
59  \endverbatim
60 
61 SourceFiles
62  kEpsilon.C
63 
64 \*---------------------------------------------------------------------------*/
65 
66 #ifndef kEpsilon_H
67 #define kEpsilon_H
68 
69 #include "RASModel.H"
70 #include "eddyViscosity.H"
71 
72 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74 namespace Foam
75 {
76 namespace RASModels
77 {
78 
79 /*---------------------------------------------------------------------------*\
80  Class kEpsilon Declaration
81 \*---------------------------------------------------------------------------*/
82 
83 template<class BasicMomentumTransportModel>
84 class kEpsilon
85 :
86  public eddyViscosity<RASModel<BasicMomentumTransportModel>>
87 {
88 protected:
89 
90  // Protected data
91 
92  // Model coefficients
93 
100 
101  // Fields
102 
105 
106 
107  // Protected Member Functions
108 
109  virtual void correctNut();
110  virtual tmp<fvScalarMatrix> kSource() const;
111  virtual tmp<fvScalarMatrix> epsilonSource() const;
112 
113 
114 public:
115 
116  typedef typename BasicMomentumTransportModel::alphaField alphaField;
117  typedef typename BasicMomentumTransportModel::rhoField rhoField;
118 
119 
120  //- Runtime type information
121  TypeName("kEpsilon");
122 
123 
124  // Constructors
125 
126  //- Construct from components
127  kEpsilon
128  (
129  const alphaField& alpha,
130  const rhoField& rho,
131  const volVectorField& U,
132  const surfaceScalarField& alphaRhoPhi,
133  const surfaceScalarField& phi,
134  const viscosity& viscosity,
135  const word& type = typeName
136  );
137 
138  //- Disallow default bitwise copy construction
139  kEpsilon(const kEpsilon&) = delete;
140 
141 
142  //- Destructor
143  virtual ~kEpsilon()
144  {}
145 
146 
147  // Member Functions
148 
149  //- Re-read model coefficients if they have changed
150  virtual bool read();
151 
152  //- Return the effective diffusivity for k
153  tmp<volScalarField> DkEff() const
154  {
155  return volScalarField::New
156  (
157  "DkEff",
158  (this->nut_/sigmak_ + this->nu())
159  );
160  }
161 
162  //- Return the effective diffusivity for epsilon
164  {
165  return volScalarField::New
166  (
167  "DepsilonEff",
168  (this->nut_/sigmaEps_ + this->nu())
169  );
170  }
171 
172  //- Return the turbulence kinetic energy
173  virtual tmp<volScalarField> k() const
174  {
175  return k_;
176  }
177 
178  //- Return the turbulence kinetic energy dissipation rate
179  virtual tmp<volScalarField> epsilon() const
180  {
181  return epsilon_;
182  }
183 
184  //- Return the turbulence specific dissipation rate
185  virtual tmp<volScalarField> omega() const
186  {
187  return volScalarField::New
188  (
189  "omega",
190  epsilon_/(Cmu_*k_)
191  );
192  }
193 
194  //- Solve the turbulence equations and correct the turbulence viscosity
195  virtual void correct();
196 
197 
198  // Member Operators
199 
200  //- Disallow default bitwise assignment
201  void operator=(const kEpsilon&) = delete;
202 };
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 } // End namespace RASModels
208 } // End namespace Foam
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 #ifdef NoRepository
213  #include "kEpsilon.C"
214 #endif
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 #endif
219 
220 // ************************************************************************* //
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,.
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:86
BasicMomentumTransportModel::alphaField alphaField
Definition: kEpsilon.H:115
volScalarField epsilon_
Definition: kEpsilon.H:103
virtual ~kEpsilon()
Destructor.
Definition: kEpsilon.H:142
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: kEpsilon.C:65
dimensionedScalar sigmak_
Definition: kEpsilon.H:97
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: kEpsilon.H:152
volScalarField k_
Definition: kEpsilon.H:102
void operator=(const kEpsilon &)=delete
Disallow default bitwise assignment.
dimensionedScalar C1_
Definition: kEpsilon.H:94
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: kEpsilon.H:172
dimensionedScalar sigmaEps_
Definition: kEpsilon.H:98
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
Definition: kEpsilon.H:184
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:218
dimensionedScalar C2_
Definition: kEpsilon.H:95
kEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: kEpsilon.C:83
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: kEpsilon.H:162
virtual void correctNut()
Definition: kEpsilon.C:41
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEpsilon.C:50
dimensionedScalar Cmu_
Definition: kEpsilon.H:93
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: kEpsilon.H:178
TypeName("kEpsilon")
Runtime type information.
dimensionedScalar C3_
Definition: kEpsilon.H:96
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kEpsilon.C:197
BasicMomentumTransportModel::rhoField rhoField
Definition: kEpsilon.H:116
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
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