RNGkEpsilon.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-2022 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::RNGkEpsilon
26 
27 Description
28  Renormalisation group k-epsilon turbulence model for incompressible and
29  compressible flows.
30 
31  Reference:
32  \verbatim
33  Yakhot, V., Orszag, S. A., Thangam, S.,
34  Gatski, T. B., & Speziale, C. G. (1992).
35  Development of turbulence models for shear flows
36  by a double expansion technique.
37  Physics of Fluids A: Fluid Dynamics (1989-1993), 4(7), 1510-1520.
38 
39  For the RDT-based compression term:
40  El Tahry, S. H. (1983).
41  k-epsilon equation for compressible reciprocating engine flows.
42  Journal of Energy, 7(4), 345-353.
43  \endverbatim
44 
45  The default model coefficients are
46  \verbatim
47  RNGkEpsilonCoeffs
48  {
49  Cmu 0.0845;
50  C1 1.42;
51  C2 1.68;
52  C3 0;
53  sigmak 0.71942;
54  sigmaEps 0.71942;
55  eta0 4.38;
56  beta 0.012;
57  }
58  \endverbatim
59 
60 SourceFiles
61  RNGkEpsilon.C
62 
63 \*---------------------------------------------------------------------------*/
64 
65 #ifndef RNGkEpsilon_H
66 #define RNGkEpsilon_H
67 
68 #include "RASModel.H"
69 #include "eddyViscosity.H"
70 
71 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 
73 namespace Foam
74 {
75 namespace RASModels
76 {
77 
78 /*---------------------------------------------------------------------------*\
79  Class RNGkEpsilon Declaration
80 \*---------------------------------------------------------------------------*/
81 
82 template<class BasicMomentumTransportModel>
83 class RNGkEpsilon
84 :
85  public eddyViscosity<RASModel<BasicMomentumTransportModel>>
86 {
87  // Private Member Functions
88 
89  using IOobject::modelName;
90 
91 
92 protected:
93 
94  // Protected data
95 
96  // Model coefficients
97 
106 
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 BasicMomentumTransportModel::alphaField alphaField;
124  typedef typename BasicMomentumTransportModel::rhoField rhoField;
125 
126 
127  //- Runtime type information
128  TypeName("RNGkEpsilon");
129 
130 
131  // Constructors
132 
133  //- Construct from components
135  (
136  const alphaField& alpha,
137  const rhoField& rho,
138  const volVectorField& U,
139  const surfaceScalarField& alphaRhoPhi,
140  const surfaceScalarField& phi,
141  const viscosity& viscosity,
142  const word& type = typeName
143  );
144 
145  //- Disallow default bitwise copy construction
146  RNGkEpsilon(const RNGkEpsilon&) = delete;
147 
148 
149  //- Destructor
150  virtual ~RNGkEpsilon()
151  {}
152 
153 
154  // Member Functions
155 
156  //- Re-read model coefficients if they have changed
157  virtual bool read();
158 
159  //- Return the effective diffusivity for k
160  tmp<volScalarField> DkEff() const
161  {
162  return volScalarField::New
163  (
164  "DkEff",
165  (this->nut_/sigmak_ + this->nu())
166  );
167  }
168 
169  //- Return the effective diffusivity for epsilon
171  {
172  return volScalarField::New
173  (
174  "DepsilonEff",
175  (this->nut_/sigmaEps_ + this->nu())
176  );
177  }
178 
179  //- Return the turbulence kinetic energy
180  virtual tmp<volScalarField> k() const
181  {
182  return k_;
183  }
184 
185  //- Return the turbulence kinetic energy dissipation rate
186  virtual tmp<volScalarField> epsilon() const
187  {
188  return epsilon_;
189  }
190 
191  //- Return the turbulence specific dissipation rate
192  virtual tmp<volScalarField> omega() const
193  {
194  return volScalarField::New
195  (
196  "omega",
197  epsilon_/(Cmu_*k_),
198  epsilon_.boundaryField().types()
199  );
200  }
201 
202  //- Solve the turbulence equations and correct the turbulence viscosity
203  virtual void correct();
204 
205 
206  // Member Operators
207 
208  //- Disallow default bitwise assignment
209  void operator=(const RNGkEpsilon&) = delete;
210 };
211 
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 } // End namespace RASModels
216 } // End namespace Foam
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 #ifdef NoRepository
221  #include "RNGkEpsilon.C"
222 #endif
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 #endif
227 
228 // ************************************************************************* //
dimensionedScalar beta_
Definition: RNGkEpsilon.H:104
dimensionedScalar C1_
Definition: RNGkEpsilon.H:98
TypeName("RNGkEpsilon")
Runtime type information.
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
Definition: RNGkEpsilon.H:191
U
Definition: pEqn.H:72
const Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual tmp< fvScalarMatrix > kSource() const
Definition: RNGkEpsilon.C:50
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:239
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
wordList types() const
Return a list of the patch field types.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
BasicMomentumTransportModel::rhoField rhoField
Definition: RNGkEpsilon.H:123
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
BasicMomentumTransportModel::alphaField alphaField
Definition: RNGkEpsilon.H:122
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: RNGkEpsilon.C:66
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: RNGkEpsilon.H:159
dimensionedScalar sigmak_
Definition: RNGkEpsilon.H:101
A class for handling words, derived from string.
Definition: word.H:59
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: RNGkEpsilon.H:169
phi
Definition: correctPhi.H:3
virtual ~RNGkEpsilon()
Destructor.
Definition: RNGkEpsilon.H:149
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: RNGkEpsilon.H:179
dimensionedScalar Cmu_
Definition: RNGkEpsilon.H:97
void operator=(const RNGkEpsilon &)=delete
Disallow default bitwise assignment.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: RNGkEpsilon.H:185
dimensionedScalar eta0_
Definition: RNGkEpsilon.H:103
Abstract base class for all fluid physical properties.
Definition: viscosity.H:49
Renormalisation group k-epsilon turbulence model for incompressible and compressible flows...
Definition: RNGkEpsilon.H:82
virtual bool read()
Re-read model coefficients if they have changed.
Definition: RNGkEpsilon.C:216
RNGkEpsilon(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: RNGkEpsilon.C:84
dimensionedScalar sigmaEps_
Definition: RNGkEpsilon.H:102
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dimensionedScalar C2_
Definition: RNGkEpsilon.H:99
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Definition: RASModel.H:208
Namespace for OpenFOAM.