RNGkEpsilon.C
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-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2020 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
27\*---------------------------------------------------------------------------*/
28
29#include "RNGkEpsilon.H"
30#include "fvOptions.H"
31#include "bound.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace RASModels
38{
39
40// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41
42template<class BasicTurbulenceModel>
44{
45 this->nut_ = Cmu_*sqr(k_)/epsilon_;
46 this->nut_.correctBoundaryConditions();
47 fv::options::New(this->mesh_).correct(this->nut_);
48
49 BasicTurbulenceModel::correctNut();
50}
51
52
53template<class BasicTurbulenceModel>
55{
57 (
59 (
60 k_,
61 dimVolume*this->rho_.dimensions()*k_.dimensions()
62 /dimTime
63 )
64 );
65}
66
67
68template<class BasicTurbulenceModel>
70{
72 (
74 (
75 epsilon_,
76 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
77 /dimTime
78 )
79 );
80}
81
82
83// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84
85template<class BasicTurbulenceModel>
87(
88 const alphaField& alpha,
89 const rhoField& rho,
90 const volVectorField& U,
91 const surfaceScalarField& alphaRhoPhi,
93 const transportModel& transport,
94 const word& propertiesName,
95 const word& type
96)
97:
98 eddyViscosity<RASModel<BasicTurbulenceModel>>
99 (
100 type,
101 alpha,
102 rho,
103 U,
104 alphaRhoPhi,
105 phi,
106 transport,
107 propertiesName
108 ),
109
110 Cmu_
111 (
112 dimensioned<scalar>::getOrAddToDict
113 (
114 "Cmu",
115 this->coeffDict_,
116 0.0845
117 )
118 ),
119 C1_
120 (
121 dimensioned<scalar>::getOrAddToDict
122 (
123 "C1",
124 this->coeffDict_,
125 1.42
126 )
127 ),
128 C2_
129 (
130 dimensioned<scalar>::getOrAddToDict
131 (
132 "C2",
133 this->coeffDict_,
134 1.68
135 )
136 ),
137 C3_
138 (
139 dimensioned<scalar>::getOrAddToDict
140 (
141 "C3",
142 this->coeffDict_,
143 0
144 )
145 ),
146 sigmak_
147 (
148 dimensioned<scalar>::getOrAddToDict
149 (
150 "sigmak",
151 this->coeffDict_,
152 0.71942
153 )
154 ),
155 sigmaEps_
156 (
157 dimensioned<scalar>::getOrAddToDict
158 (
159 "sigmaEps",
160 this->coeffDict_,
161 0.71942
162 )
163 ),
164 eta0_
165 (
166 dimensioned<scalar>::getOrAddToDict
167 (
168 "eta0",
169 this->coeffDict_,
170 4.38
171 )
172 ),
173 beta_
174 (
175 dimensioned<scalar>::getOrAddToDict
176 (
177 "beta",
178 this->coeffDict_,
179 0.012
180 )
181 ),
182
183 k_
184 (
186 (
187 IOobject::groupName("k", alphaRhoPhi.group()),
188 this->runTime_.timeName(),
189 this->mesh_,
190 IOobject::MUST_READ,
191 IOobject::AUTO_WRITE
192 ),
193 this->mesh_
194 ),
195 epsilon_
196 (
198 (
199 IOobject::groupName("epsilon", alphaRhoPhi.group()),
200 this->runTime_.timeName(),
201 this->mesh_,
202 IOobject::MUST_READ,
203 IOobject::AUTO_WRITE
204 ),
205 this->mesh_
206 )
207{
208 bound(k_, this->kMin_);
209 bound(epsilon_, this->epsilonMin_);
210
211 if (type == typeName)
212 {
213 this->printCoeffs(type);
214 }
215}
216
217
218// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
219
220template<class BasicTurbulenceModel>
222{
224 {
225 Cmu_.readIfPresent(this->coeffDict());
226 C1_.readIfPresent(this->coeffDict());
227 C2_.readIfPresent(this->coeffDict());
228 C3_.readIfPresent(this->coeffDict());
229 sigmak_.readIfPresent(this->coeffDict());
230 sigmaEps_.readIfPresent(this->coeffDict());
231 eta0_.readIfPresent(this->coeffDict());
232 beta_.readIfPresent(this->coeffDict());
233
234 return true;
235 }
236
237 return false;
238}
239
240
241template<class BasicTurbulenceModel>
243{
244 if (!this->turbulence_)
245 {
246 return;
247 }
248
249 // Local references
250 const alphaField& alpha = this->alpha_;
251 const rhoField& rho = this->rho_;
252 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
253 const volVectorField& U = this->U_;
254 const volScalarField& nut = this->nut_;
256
258
260 (
261 fvc::div(fvc::absolute(this->phi(), U))().v()
262 );
263
265 const volScalarField::Internal GbyNu
266 (
267 this->type() + ":GbyNu",
268 tgradU().v() && dev(twoSymm(tgradU().v()))
269 );
270 tgradU.clear();
271
272 const volScalarField::Internal G(this->GName(), nut()*GbyNu);
273
274 const volScalarField::Internal eta(sqrt(mag(GbyNu))*k_/epsilon_);
275 const volScalarField::Internal eta3(eta*sqr(eta));
276
278 (
279 ((eta*(-eta/eta0_ + scalar(1)))/(beta_*eta3 + scalar(1)))
280 );
281
282 // Update epsilon and G at the wall
283 epsilon_.boundaryFieldRef().updateCoeffs();
284
285 // Dissipation equation
287 (
288 fvm::ddt(alpha, rho, epsilon_)
289 + fvm::div(alphaRhoPhi, epsilon_)
290 - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
291 ==
292 (C1_ - R)*alpha()*rho()*GbyNu*Cmu_*k_()
293 - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
294 - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
295 + epsilonSource()
296 + fvOptions(alpha, rho, epsilon_)
297 );
298
299 epsEqn.ref().relax();
300 fvOptions.constrain(epsEqn.ref());
301 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
302 solve(epsEqn);
303 fvOptions.correct(epsilon_);
304 bound(epsilon_, this->epsilonMin_);
305
306
307 // Turbulent kinetic energy equation
308
310 (
311 fvm::ddt(alpha, rho, k_)
312 + fvm::div(alphaRhoPhi, k_)
313 - fvm::laplacian(alpha*rho*DkEff(), k_)
314 ==
315 alpha()*rho()*G
316 - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
317 - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
318 + kSource()
319 + fvOptions(alpha, rho, k_)
320 );
321
322 kEqn.ref().relax();
323 fvOptions.constrain(kEqn.ref());
324 solve(kEqn);
325 fvOptions.correct(k_);
326 bound(k_, this->kMin_);
327
328 correctNut();
329}
330
331
332// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333
334} // End namespace RASModels
335} // End namespace Foam
336
337// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
Bound the given scalar field if it has gone unbounded.
fv::options & fvOptions
surfaceScalarField & phi
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:55
Renormalization group k-epsilon turbulence model for incompressible and compressible flows.
Definition: RNGkEpsilon.H:91
BasicTurbulenceModel::alphaField alphaField
Definition: RNGkEpsilon.H:132
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: RNGkEpsilon.C:69
BasicTurbulenceModel::rhoField rhoField
Definition: RNGkEpsilon.H:133
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:242
virtual tmp< fvScalarMatrix > kSource() const
Definition: RNGkEpsilon.C:54
BasicTurbulenceModel::transportModel transportModel
Definition: RNGkEpsilon.H:134
virtual bool read()
Re-read model coefficients if they have changed.
Definition: RNGkEpsilon.C:221
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:58
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Finite-volume options.
Definition: fvOptions.H:59
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
thermo correct()
scalar nut
zeroField divU
Definition: alphaSuSp.H:3
word timeName
Definition: getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
volScalarField & alpha
CEqn solve()