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-2019 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  Renormalization 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.33;
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 BasicTurbulenceModel>
83 class RNGkEpsilon
84 :
85  public eddyViscosity<RASModel<BasicTurbulenceModel>>
86 {
87 protected:
88 
89  // Protected data
90 
91  // Model coefficients
92 
101 
102 
103  // Fields
107 
108 
109  // Protected Member Functions
110 
111  virtual void correctNut();
112  virtual tmp<fvScalarMatrix> kSource() const;
113  virtual tmp<fvScalarMatrix> epsilonSource() const;
114 
115 
116 public:
118  typedef typename BasicTurbulenceModel::alphaField alphaField;
119  typedef typename BasicTurbulenceModel::rhoField rhoField;
120  typedef typename BasicTurbulenceModel::transportModel transportModel;
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 transportModel& transport,
138  const word& propertiesName = turbulenceModel::propertiesName,
139  const word& type = typeName
140  );
141 
142  //- Disallow default bitwise copy construction
143  RNGkEpsilon(const RNGkEpsilon&) = delete;
144 
145 
146  //- Destructor
147  virtual ~RNGkEpsilon()
148  {}
149 
150 
151  // Member Functions
152 
153  //- Re-read model coefficients if they have changed
154  virtual bool read();
155 
156  //- Return the effective diffusivity for k
157  tmp<volScalarField> DkEff() const
158  {
159  return volScalarField::New
160  (
161  "DkEff",
162  (this->nut_/sigmak_ + this->nu())
163  );
164  }
165 
166  //- Return the effective diffusivity for epsilon
168  {
169  return volScalarField::New
170  (
171  "DepsilonEff",
172  (this->nut_/sigmaEps_ + this->nu())
173  );
174  }
175 
176  //- Return the turbulence kinetic energy
177  virtual tmp<volScalarField> k() const
178  {
179  return k_;
180  }
181 
182  //- Return the turbulence kinetic energy dissipation rate
183  virtual tmp<volScalarField> epsilon() const
184  {
185  return epsilon_;
186  }
187 
188  //- Solve the turbulence equations and correct the turbulence viscosity
189  virtual void correct();
190 
191 
192  // Member Operators
193 
194  //- Disallow default bitwise assignment
195  void operator=(const RNGkEpsilon&) = delete;
196 };
197 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 } // End namespace RASModels
202 } // End namespace Foam
203 
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 
206 #ifdef NoRepository
207  #include "RNGkEpsilon.C"
208 #endif
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 #endif
213 
214 // ************************************************************************* //
dimensionedScalar beta_
Definition: RNGkEpsilon.H:99
dimensionedScalar C1_
Definition: RNGkEpsilon.H:93
surfaceScalarField & phi
TypeName("RNGkEpsilon")
Runtime type information.
virtual tmp< fvScalarMatrix > kSource() const
Definition: RNGkEpsilon.C:51
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:241
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
dimensionedScalar C3_
Definition: RNGkEpsilon.H:95
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: RNGkEpsilon.C:66
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: RNGkEpsilon.H:156
dimensionedScalar sigmak_
Definition: RNGkEpsilon.H:96
static const word propertiesName
Default name of the turbulence properties dictionary.
BasicTurbulenceModel::alphaField alphaField
Definition: RNGkEpsilon.H:117
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:166
RNGkEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: RNGkEpsilon.C:84
virtual ~RNGkEpsilon()
Destructor.
Definition: RNGkEpsilon.H:146
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: RNGkEpsilon.H:176
dimensionedScalar Cmu_
Definition: RNGkEpsilon.H:92
void operator=(const RNGkEpsilon &)=delete
Disallow default bitwise assignment.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: RNGkEpsilon.H:182
dimensionedScalar eta0_
Definition: RNGkEpsilon.H:98
Renormalization 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:218
U
Definition: pEqn.H:72
BasicTurbulenceModel::rhoField rhoField
Definition: RNGkEpsilon.H:118
BasicTurbulenceModel::transportModel transportModel
Definition: RNGkEpsilon.H:119
dimensionedScalar sigmaEps_
Definition: RNGkEpsilon.H:97
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:94
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
volScalarField & nu
Namespace for OpenFOAM.