kEpsilon.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::RASModels::kEpsilon
29 
30 Group
31  grpRASTurbulence
32 
33 Description
34  Standard k-epsilon turbulence model for incompressible and compressible
35  flows including rapid distortion theory (RDT) based compression term.
36 
37  Reference:
38  \verbatim
39  Standard model:
40  Launder, B. E., & Spalding, D. B. (1972).
41  Lectures in mathematical models of turbulence.
42 
43  Launder, B. E., & Spalding, D. B. (1974).
44  The numerical computation of turbulent flows.
45  Computer methods in applied mechanics and engineering,
46  3(2), 269-289.
47 
48  For the RDT-based compression term:
49  El Tahry, S. H. (1983).
50  k-epsilon equation for compressible reciprocating engine flows.
51  Journal of Energy, 7(4), 345-353.
52 
53  Two-layer wall treatment to allow operation on low-Re grids:
54  Jongen, T. and Marx, Y.P. (1997) “Design of an Unconditionally
55  Stable, Positive Scheme for the k−epsilon and Two-Layer Turbulence
56  Models”, Computers & Fluids, 26(5), 469-487.
57  \endverbatim
58 
59  The default model coefficients are
60  \verbatim
61  kEpsilonCoeffs
62  {
63  Cmu 0.09;
64  C1 1.44;
65  C2 1.92;
66  C3 -0.33;
67  sigmak 1.0;
68  sigmaEps 1.3;
69  twoLayerTreatment false;
70 
71  // If two twoLayerTreatment = true,
72 
73  // Critical turbulence Reynolds number marking the separation
74  // of inner-layer and outer-layer. Can be adjusted from 50.0
75  // (low-speed flows) to 200.0 (default value).
76  ReyStar 200.0;
77 
78  // Smoothness of transition from inner-layer to outer-layer.
79  // Typically from 0.05 (sharp transition) to 0.2 (smooth transition).
80  // For unstable cases, prefer larger values.
81  ReyFactor 0.05;
82  }
83  \endverbatim
84 
85 SourceFiles
86  kEpsilon.C
87 
88 \*---------------------------------------------------------------------------*/
89 
90 #ifndef kEpsilon_H
91 #define kEpsilon_H
92 
93 #include "RASModel.H"
94 #include "eddyViscosity.H"
95 
96 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97 
98 namespace Foam
99 {
100 namespace RASModels
101 {
102 
103 /*---------------------------------------------------------------------------*\
104  Class kEpsilon Declaration
105 \*---------------------------------------------------------------------------*/
106 
107 template<class BasicTurbulenceModel>
108 class kEpsilon
109 :
110  public eddyViscosity<RASModel<BasicTurbulenceModel>>
111 {
112  // Private Member Functions
113 
114  //- No copy construct
115  kEpsilon(const kEpsilon&) = delete;
116 
117  //- No copy assignment
118  void operator=(const kEpsilon&) = delete;
119 
120 
121 protected:
122 
123  // Protected data
124 
125  // Model coefficients
133 
134 
135  // Fields
139 
140 
141  // Two-layer wall treatment
146  std::unique_ptr<volScalarField::Internal> lEpsPtr_;
147  std::unique_ptr<volScalarField> lambdaEpsPtr_;
148 
149 
150  // Protected Member Functions
151 
152  virtual void correctNut();
153  virtual tmp<fvScalarMatrix> kSource() const;
154  virtual tmp<fvScalarMatrix> epsilonSource() const;
155 
156 
157 public:
159  typedef typename BasicTurbulenceModel::alphaField alphaField;
160  typedef typename BasicTurbulenceModel::rhoField rhoField;
161  typedef typename BasicTurbulenceModel::transportModel transportModel;
162 
163 
164  //- Runtime type information
165  TypeName("kEpsilon");
166 
167 
168  // Constructors
169 
170  //- Construct from components
171  kEpsilon
172  (
173  const alphaField& alpha,
174  const rhoField& rho,
175  const volVectorField& U,
176  const surfaceScalarField& alphaRhoPhi,
177  const surfaceScalarField& phi,
178  const transportModel& transport,
179  const word& propertiesName = turbulenceModel::propertiesName,
180  const word& type = typeName
181  );
182 
183 
184  //- Destructor
185  virtual ~kEpsilon() = default;
186 
187 
188  // Member Functions
189 
190  //- Re-read model coefficients if they have changed
191  virtual bool read();
192 
193  //- Return the effective diffusivity for k
194  tmp<volScalarField> DkEff() const
195  {
197  (
198  "DkEff",
199  (this->nut_/sigmak_ + this->nu())
200  );
201  }
202 
203  //- Return the effective diffusivity for epsilon
205  {
207  (
208  "DepsilonEff",
209  (this->nut_/sigmaEps_ + this->nu())
210  );
211  }
212 
213  //- Return the turbulence kinetic energy
214  virtual tmp<volScalarField> k() const
215  {
216  return k_;
217  }
218 
219  //- Return the turbulence kinetic energy dissipation rate
220  virtual tmp<volScalarField> epsilon() const
221  {
222  return epsilon_;
223  }
224 
225  //- Solve the turbulence equations and correct the turbulence viscosity
226  virtual void correct();
227 };
228 
229 
230 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 
232 } // End namespace RASModels
233 } // End namespace Foam
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 #ifdef NoRepository
238  #include "kEpsilon.C"
239 #endif
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 #endif
244 
245 // ************************************************************************* //
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: kEpsilon.C:105
volScalarField epsilon_
Definition: kEpsilon.H:137
TypeName("kEpsilon")
Runtime type information.
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEpsilon.C:75
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: kEpsilon.H:227
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:51
dimensionedScalar Cmu_
Definition: kEpsilon.H:126
dimensionedScalar sigmaEps_
Definition: kEpsilon.H:131
BasicTurbulenceModel::alphaField alphaField
Definition: kEpsilon.H:158
BasicTurbulenceModel::alphaField alphaField
dimensionedScalar sigmak_
Definition: kEpsilon.H:130
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:801
dimensionedScalar C2_
Definition: kEpsilon.H:128
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:103
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: kEpsilon.H:215
dimensionedScalar C1_
Definition: kEpsilon.H:127
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:215
std::unique_ptr< volScalarField::Internal > lEpsPtr_
Definition: kEpsilon.H:145
std::unique_ptr< volScalarField > lambdaEpsPtr_
Definition: kEpsilon.H:146
dimensionedScalar ReyFactor_
Definition: kEpsilon.H:144
virtual ~kEpsilon()=default
Destructor.
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kEpsilon.C:294
U
Definition: pEqn.H:72
BasicTurbulenceModel::transportModel transportModel
virtual void correctNut()
Definition: kEpsilon.C:37
dimensionedScalar C3_
Definition: kEpsilon.H:129
dimensionedScalar ReyStar_
Definition: kEpsilon.H:143
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:318
BasicTurbulenceModel::rhoField rhoField
Definition: kEpsilon.H:159
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
BasicTurbulenceModel::transportModel transportModel
Definition: kEpsilon.H:160
volScalarField k_
Definition: kEpsilon.H:136
volScalarField & nu
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: kEpsilon.H:203
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: kEpsilon.H:235
Namespace for OpenFOAM.