kEpsilon.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 "kEpsilon.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.09
117 )
118 ),
119 C1_
120 (
121 dimensioned<scalar>::getOrAddToDict
122 (
123 "C1",
124 this->coeffDict_,
125 1.44
126 )
127 ),
128 C2_
129 (
130 dimensioned<scalar>::getOrAddToDict
131 (
132 "C2",
133 this->coeffDict_,
134 1.92
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 1.0
153 )
154 ),
155 sigmaEps_
156 (
157 dimensioned<scalar>::getOrAddToDict
158 (
159 "sigmaEps",
160 this->coeffDict_,
161 1.3
163 ),
164
165 k_
166 (
168 (
169 IOobject::groupName("k", alphaRhoPhi.group()),
170 this->runTime_.timeName(),
171 this->mesh_,
172 IOobject::MUST_READ,
173 IOobject::AUTO_WRITE
174 ),
175 this->mesh_
176 ),
177 epsilon_
178 (
180 (
181 IOobject::groupName("epsilon", alphaRhoPhi.group()),
182 this->runTime_.timeName(),
183 this->mesh_,
184 IOobject::MUST_READ,
185 IOobject::AUTO_WRITE
186 ),
187 this->mesh_
188 )
189{
190 bound(k_, this->kMin_);
191 bound(epsilon_, this->epsilonMin_);
192
193 if (type == typeName)
194 {
195 this->printCoeffs(type);
196 }
197}
198
199
200// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
201
202template<class BasicTurbulenceModel>
204{
206 {
207 Cmu_.readIfPresent(this->coeffDict());
208 C1_.readIfPresent(this->coeffDict());
209 C2_.readIfPresent(this->coeffDict());
210 C3_.readIfPresent(this->coeffDict());
211 sigmak_.readIfPresent(this->coeffDict());
212 sigmaEps_.readIfPresent(this->coeffDict());
213
214 return true;
215 }
216
217 return false;
218}
219
220
221template<class BasicTurbulenceModel>
223{
224 if (!this->turbulence_)
225 {
226 return;
227 }
228
229 // Local references
230 const alphaField& alpha = this->alpha_;
231 const rhoField& rho = this->rho_;
232 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
233 const volVectorField& U = this->U_;
234 const volScalarField& nut = this->nut_;
235
236 fv::options& fvOptions(fv::options::New(this->mesh_));
237
239
241 (
242 fvc::div(fvc::absolute(this->phi(), U))().v()
243 );
244
245 tmp<volTensorField> tgradU = fvc::grad(U);
246 const volScalarField::Internal GbyNu
247 (
248 this->type() + ":GbyNu",
249 tgradU().v() && dev(twoSymm(tgradU().v()))
250 );
251 const volScalarField::Internal G(this->GName(), nut()*GbyNu);
252 tgradU.clear();
253
254 // Update epsilon and G at the wall
255 epsilon_.boundaryFieldRef().updateCoeffs();
256
257 // Dissipation equation
259 (
260 fvm::ddt(alpha, rho, epsilon_)
261 + fvm::div(alphaRhoPhi, epsilon_)
262 - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
263 ==
264 C1_*alpha()*rho()*GbyNu*Cmu_*k_()
265 - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
266 - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
267 + epsilonSource()
268 + fvOptions(alpha, rho, epsilon_)
269 );
270
271 epsEqn.ref().relax();
272 fvOptions.constrain(epsEqn.ref());
273 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
274 solve(epsEqn);
275 fvOptions.correct(epsilon_);
276 bound(epsilon_, this->epsilonMin_);
277
278 // Turbulent kinetic energy equation
280 (
281 fvm::ddt(alpha, rho, k_)
282 + fvm::div(alphaRhoPhi, k_)
283 - fvm::laplacian(alpha*rho*DkEff(), k_)
284 ==
285 alpha()*rho()*G
286 - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
287 - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
288 + kSource()
289 + fvOptions(alpha, rho, k_)
290 );
291
292 kEqn.ref().relax();
293 fvOptions.constrain(kEqn.ref());
294 solve(kEqn);
295 fvOptions.correct(k_);
296 bound(k_, this->kMin_);
297
298 correctNut();
299}
300
301
302// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303
304} // End namespace RASModels
305} // End namespace Foam
306
307// ************************************************************************* //
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
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:92
BasicTurbulenceModel::alphaField alphaField
Definition: kEpsilon.H:130
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: kEpsilon.C:69
BasicTurbulenceModel::rhoField rhoField
Definition: kEpsilon.H:131
virtual void correctNut()
Definition: kEpsilon.C:43
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEpsilon.C:54
BasicTurbulenceModel::transportModel transportModel
Definition: kEpsilon.H:132
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
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
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
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
volScalarField & alpha
CEqn solve()