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-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::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 
88 protected:
89 
90  // Protected data
91 
92  // Model coefficients
93 
102 
103 
104  // Fields
105 
108 
109 
110  // Protected Member Functions
111 
112  virtual void correctNut();
113  virtual tmp<fvScalarMatrix> kSource() const;
114  virtual tmp<fvScalarMatrix> epsilonSource() const;
115 
116 
117 public:
118 
119  typedef typename BasicMomentumTransportModel::alphaField alphaField;
120  typedef typename BasicMomentumTransportModel::rhoField rhoField;
121 
122 
123  //- Runtime type information
124  TypeName("RNGkEpsilon");
125 
126 
127  // Constructors
128 
129  //- Construct from components
131  (
132  const alphaField& alpha,
133  const rhoField& rho,
134  const volVectorField& U,
135  const surfaceScalarField& alphaRhoPhi,
136  const surfaceScalarField& phi,
137  const viscosity& viscosity,
138  const word& type = typeName
139  );
140 
141  //- Disallow default bitwise copy construction
142  RNGkEpsilon(const RNGkEpsilon&) = delete;
143 
144 
145  //- Destructor
146  virtual ~RNGkEpsilon()
147  {}
148 
149 
150  // Member Functions
151 
152  //- Re-read model coefficients if they have changed
153  virtual bool read();
154 
155  //- Return the effective diffusivity for k
156  tmp<volScalarField> DkEff() const
157  {
158  return volScalarField::New
159  (
160  "DkEff",
161  (this->nut_/sigmak_ + this->nu())
162  );
163  }
164 
165  //- Return the effective diffusivity for epsilon
167  {
168  return volScalarField::New
169  (
170  "DepsilonEff",
171  (this->nut_/sigmaEps_ + this->nu())
172  );
173  }
174 
175  //- Return the turbulence kinetic energy
176  virtual tmp<volScalarField> k() const
177  {
178  return k_;
179  }
180 
181  //- Return the turbulence kinetic energy dissipation rate
182  virtual tmp<volScalarField> epsilon() const
183  {
184  return epsilon_;
185  }
186 
187  //- Return the turbulence specific dissipation rate
188  virtual tmp<volScalarField> omega() const
189  {
190  return volScalarField::New
191  (
192  "omega",
193  epsilon_/(Cmu_*k_)
194  );
195  }
196 
197  //- Solve the turbulence equations and correct the turbulence viscosity
198  virtual void correct();
199 
200 
201  // Member Operators
202 
203  //- Disallow default bitwise assignment
204  void operator=(const RNGkEpsilon&) = delete;
205 };
206 
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
210 } // End namespace RASModels
211 } // End namespace Foam
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 #ifdef NoRepository
216  #include "RNGkEpsilon.C"
217 #endif
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 #endif
222 
223 // ************************************************************************* //
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,.
Renormalisation group k-epsilon turbulence model for incompressible and compressible flows.
Definition: RNGkEpsilon.H:85
BasicMomentumTransportModel::alphaField alphaField
Definition: RNGkEpsilon.H:118
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: RNGkEpsilon.C:66
dimensionedScalar sigmak_
Definition: RNGkEpsilon.H:97
TypeName("RNGkEpsilon")
Runtime type information.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: RNGkEpsilon.H:155
dimensionedScalar C1_
Definition: RNGkEpsilon.H:94
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: RNGkEpsilon.H:175
dimensionedScalar sigmaEps_
Definition: RNGkEpsilon.H:98
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
Definition: RNGkEpsilon.H:187
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:239
virtual ~RNGkEpsilon()
Destructor.
Definition: RNGkEpsilon.H:145
void operator=(const RNGkEpsilon &)=delete
Disallow default bitwise assignment.
dimensionedScalar C2_
Definition: RNGkEpsilon.H:95
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: RNGkEpsilon.H:165
virtual tmp< fvScalarMatrix > kSource() const
Definition: RNGkEpsilon.C:50
dimensionedScalar eta0_
Definition: RNGkEpsilon.H:99
dimensionedScalar Cmu_
Definition: RNGkEpsilon.H:93
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: RNGkEpsilon.H:181
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 C3_
Definition: RNGkEpsilon.H:96
virtual bool read()
Re-read model coefficients if they have changed.
Definition: RNGkEpsilon.C:216
dimensionedScalar beta_
Definition: RNGkEpsilon.H:100
BasicMomentumTransportModel::rhoField rhoField
Definition: RNGkEpsilon.H:119
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