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-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::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.33;
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 BasicTurbulenceModel>
84 class kEpsilon
85 :
86  public eddyViscosity<RASModel<BasicTurbulenceModel>>
87 {
88  // Private Member Functions
89 
90  // Disallow default bitwise copy construct and assignment
91  kEpsilon(const kEpsilon&);
92  void operator=(const kEpsilon&);
93 
94 
95 protected:
96 
97  // Protected data
98 
99  // Model coefficients
107 
108  // Fields
112 
113 
114  // Protected Member Functions
115 
116  virtual void correctNut();
117  virtual tmp<fvScalarMatrix> kSource() const;
118  virtual tmp<fvScalarMatrix> epsilonSource() const;
119 
120 
121 public:
123  typedef typename BasicTurbulenceModel::alphaField alphaField;
124  typedef typename BasicTurbulenceModel::rhoField rhoField;
125  typedef typename BasicTurbulenceModel::transportModel transportModel;
126 
127 
128  //- Runtime type information
129  TypeName("kEpsilon");
130 
131 
132  // Constructors
133 
134  //- Construct from components
135  kEpsilon
136  (
137  const alphaField& alpha,
138  const rhoField& rho,
139  const volVectorField& U,
140  const surfaceScalarField& alphaRhoPhi,
141  const surfaceScalarField& phi,
142  const transportModel& transport,
143  const word& propertiesName = turbulenceModel::propertiesName,
144  const word& type = typeName
145  );
146 
147 
148  //- Destructor
149  virtual ~kEpsilon()
150  {}
151 
152 
153  // Member Functions
154 
155  //- Re-read model coefficients if they have changed
156  virtual bool read();
157 
158  //- Return the effective diffusivity for k
159  tmp<volScalarField> DkEff() const
160  {
161  return tmp<volScalarField>
162  (
163  new volScalarField
164  (
165  "DkEff",
166  (this->nut_/sigmak_ + this->nu())
167  )
168  );
169  }
170 
171  //- Return the effective diffusivity for epsilon
173  {
174  return tmp<volScalarField>
175  (
176  new volScalarField
177  (
178  "DepsilonEff",
179  (this->nut_/sigmaEps_ + this->nu())
180  )
181  );
182  }
183 
184  //- Return the turbulence kinetic energy
185  virtual tmp<volScalarField> k() const
186  {
187  return k_;
188  }
189 
190  //- Return the turbulence kinetic energy dissipation rate
191  virtual tmp<volScalarField> epsilon() const
192  {
193  return epsilon_;
194  }
195 
196  //- Solve the turbulence equations and correct the turbulence viscosity
197  virtual void correct();
198 };
199 
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 } // End namespace RASModels
204 } // End namespace Foam
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 #ifdef NoRepository
209  #include "kEpsilon.C"
210 #endif
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 #endif
215 
216 // ************************************************************************* //
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: kEpsilon.C:66
volScalarField epsilon_
Definition: kEpsilon.H:110
surfaceScalarField & phi
TypeName("kEpsilon")
Runtime type information.
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEpsilon.C:51
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: kEpsilon.H:184
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
dimensionedScalar Cmu_
Definition: kEpsilon.H:100
dimensionedScalar sigmaEps_
Definition: kEpsilon.H:105
BasicTurbulenceModel::alphaField alphaField
Definition: kEpsilon.H:122
dimensionedScalar sigmak_
Definition: kEpsilon.H:104
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar C2_
Definition: kEpsilon.H:102
virtual ~kEpsilon()
Destructor.
Definition: kEpsilon.H:148
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:83
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: kEpsilon.H:171
dimensionedScalar C1_
Definition: kEpsilon.H:101
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kEpsilon.C:200
U
Definition: pEqn.H:72
virtual void correctNut()
Definition: kEpsilon.C:40
dimensionedScalar C3_
Definition: kEpsilon.H:103
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:221
BasicTurbulenceModel::rhoField rhoField
Definition: kEpsilon.H:123
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
BasicTurbulenceModel::transportModel transportModel
Definition: kEpsilon.H:124
volScalarField k_
Definition: kEpsilon.H:109
volScalarField & nu
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: kEpsilon.H:158
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: kEpsilon.H:190
Namespace for OpenFOAM.