Loading...
Searching...
No Matches
RNGkEpsilon.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-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::RASModels::RNGkEpsilon
29
30Group
31 grpRASTurbulence
32
33Description
34 Renormalization group k-epsilon turbulence model for incompressible and
35 compressible flows.
36
37 Reference:
38 \verbatim
39 Yakhot, V., Orszag, S. A., Thangam, S.,
40 Gatski, T. B., & Speziale, C. G. (1992).
41 Development of turbulence models for shear flows
42 by a double expansion technique.
43 Physics of Fluids A: Fluid Dynamics (1989-1993), 4(7), 1510-1520.
44
45 For the RDT-based compression term:
46 El Tahry, S. H. (1983).
47 k-epsilon equation for compressible reciprocating engine flows.
48 Journal of Energy, 7(4), 345-353.
49 \endverbatim
50
51 The default model coefficients are
52 \verbatim
53 RNGkEpsilonCoeffs
54 {
55 Cmu 0.0845;
56 C1 1.42;
57 C2 1.68;
58 C3 -0.33;
59 sigmak 0.71942;
60 sigmaEps 0.71942;
61 eta0 4.38;
62 beta 0.012;
63 }
64 \endverbatim
65
66SourceFiles
67 RNGkEpsilon.C
68
69\*---------------------------------------------------------------------------*/
70
71#ifndef RNGkEpsilon_H
72#define RNGkEpsilon_H
73
74#include "RASModel.H"
75#include "eddyViscosity.H"
76
77// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78
79namespace Foam
80{
81namespace RASModels
82{
83
84/*---------------------------------------------------------------------------*\
85 Class RNGkEpsilon Declaration
86\*---------------------------------------------------------------------------*/
87
88template<class BasicTurbulenceModel>
89class RNGkEpsilon
90:
91 public eddyViscosity<RASModel<BasicTurbulenceModel>>
92{
93 // Private Member Functions
94
95 //- No copy construct
96 RNGkEpsilon(const RNGkEpsilon&) = delete;
97
98 //- No copy assignment
99 void operator=(const RNGkEpsilon&) = delete;
100
101
102protected:
103
104 // Protected data
105
106 // Model coefficients
116
117
118 // Fields
122
123
124 // Protected Member Functions
126 virtual void correctNut();
128 virtual tmp<fvScalarMatrix> epsilonSource() const;
129
130
131public:
133 typedef typename BasicTurbulenceModel::alphaField alphaField;
134 typedef typename BasicTurbulenceModel::rhoField rhoField;
135 typedef typename BasicTurbulenceModel::transportModel transportModel;
136
137
138 //- Runtime type information
139 TypeName("RNGkEpsilon");
141
142 // Constructors
143
144 //- Construct from components
145 RNGkEpsilon
146 (
147 const alphaField& alpha,
148 const rhoField& rho,
149 const volVectorField& U,
150 const surfaceScalarField& alphaRhoPhi,
151 const surfaceScalarField& phi,
152 const transportModel& transport,
153 const word& propertiesName = turbulenceModel::propertiesName,
154 const word& type = typeName
155 );
156
157
158 //- Destructor
159 virtual ~RNGkEpsilon() = default;
160
161
162 // Member Functions
163
164 //- Re-read model coefficients if they have changed
165 virtual bool read();
166
167 //- Return the effective diffusivity for k
169 {
171 (
173 (
174 "DkEff",
175 (this->nut_/sigmak_ + this->nu())
176 )
177 );
178 }
179
180 //- Return the effective diffusivity for epsilon
182 {
184 (
186 (
187 "DepsilonEff",
188 (this->nut_/sigmaEps_ + this->nu())
189 )
190 );
191 }
193 //- Return the turbulence kinetic energy
194 virtual tmp<volScalarField> k() const
195 {
196 return k_;
197 }
198
199 //- Return the turbulence kinetic energy dissipation rate
200 virtual tmp<volScalarField> epsilon() const
201 {
202 return epsilon_;
203 }
204
205 //- Solve the turbulence equations and correct the turbulence viscosity
206 virtual void correct();
208
209
210// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211
212} // End namespace RASModels
213} // End namespace Foam
214
215// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216
217#ifdef NoRepository
218 #include "RNGkEpsilon.C"
219#endif
220
221// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222
223#endif
224
225// ************************************************************************* //
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition RNGkEpsilon.C:58
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::rhoField rhoField
TypeName("RNGkEpsilon")
Runtime type information.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual ~RNGkEpsilon()=default
Destructor.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
virtual tmp< fvScalarMatrix > kSource() const
Definition RNGkEpsilon.C:47
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::transportModel transportModel
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual bool read()
Re-read model coefficients if they have changed.
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:72
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
volScalarField nut_
eddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
A class for managing temporary objects.
Definition tmp.H:75
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & nu
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68