realizableKE.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 "realizableKE.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 const volTensorField& gradU,
46 const volScalarField& S2,
47 const volScalarField& magS
48)
49{
50 tmp<volSymmTensorField> tS = dev(symm(gradU));
51 const volSymmTensorField& S = tS();
52
54 (
55 (2*sqrt(2.0))*((S&S)&&S)
56 /(
57 magS*S2
58 + dimensionedScalar("small", dimensionSet(0, 0, -3, 0, 0), SMALL)
59 )
60 );
61
62 tS.clear();
63
65 (
66 (1.0/3.0)*acos(min(max(sqrt(6.0)*W, -scalar(1)), scalar(1)))
67 );
68 volScalarField As(sqrt(6.0)*cos(phis));
69 volScalarField Us(sqrt(S2/2.0 + magSqr(skew(gradU))));
70
71 return 1.0/(A0_ + As*Us*k_/epsilon_);
72}
73
74
75template<class BasicTurbulenceModel>
77(
78 const volTensorField& gradU,
79 const volScalarField& S2,
80 const volScalarField& magS
81)
82{
83 this->nut_ = rCmu(gradU, S2, magS)*sqr(k_)/epsilon_;
84 this->nut_.correctBoundaryConditions();
85 fv::options::New(this->mesh_).correct(this->nut_);
86
87 BasicTurbulenceModel::correctNut();
88}
89
90
91template<class BasicTurbulenceModel>
93{
94 tmp<volTensorField> tgradU = fvc::grad(this->U_);
95 volScalarField S2(2*magSqr(dev(symm(tgradU()))));
96 volScalarField magS(sqrt(S2));
97 correctNut(tgradU(), S2, magS);
98}
99
100
101template<class BasicTurbulenceModel>
103{
105 (
107 (
108 k_,
109 dimVolume*this->rho_.dimensions()*k_.dimensions()
110 /dimTime
111 )
112 );
113}
114
115
116template<class BasicTurbulenceModel>
118{
120 (
122 (
123 epsilon_,
124 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
125 /dimTime
126 )
127 );
128}
129
130
131// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
132
133template<class BasicTurbulenceModel>
135(
136 const alphaField& alpha,
137 const rhoField& rho,
138 const volVectorField& U,
139 const surfaceScalarField& alphaRhoPhi,
140 const surfaceScalarField& phi,
141 const transportModel& transport,
142 const word& propertiesName,
143 const word& type
144)
145:
146 eddyViscosity<RASModel<BasicTurbulenceModel>>
147 (
148 type,
149 alpha,
150 rho,
151 U,
152 alphaRhoPhi,
153 phi,
154 transport,
155 propertiesName
156 ),
157 A0_
158 (
159 dimensioned<scalar>::getOrAddToDict
160 (
161 "A0",
162 this->coeffDict_,
163 4.0
164 )
165 ),
166 C2_
167 (
168 dimensioned<scalar>::getOrAddToDict
169 (
170 "C2",
171 this->coeffDict_,
172 1.9
173 )
174 ),
175 sigmak_
176 (
177 dimensioned<scalar>::getOrAddToDict
178 (
179 "sigmak",
180 this->coeffDict_,
181 1.0
182 )
183 ),
184 sigmaEps_
185 (
186 dimensioned<scalar>::getOrAddToDict
187 (
188 "sigmaEps",
189 this->coeffDict_,
190 1.2
191 )
192 ),
193
194 k_
195 (
197 (
198 IOobject::groupName("k", U.group()),
199 this->runTime_.timeName(),
200 this->mesh_,
201 IOobject::MUST_READ,
202 IOobject::AUTO_WRITE
203 ),
204 this->mesh_
205 ),
206 epsilon_
207 (
209 (
210 IOobject::groupName("epsilon", U.group()),
211 this->runTime_.timeName(),
212 this->mesh_,
213 IOobject::MUST_READ,
214 IOobject::AUTO_WRITE
215 ),
216 this->mesh_
217 )
218{
219 bound(k_, this->kMin_);
220 bound(epsilon_, this->epsilonMin_);
221
222 if (type == typeName)
223 {
224 this->printCoeffs(type);
225 }
226}
227
228
229// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
230
231template<class BasicTurbulenceModel>
233{
235 {
236 A0_.readIfPresent(this->coeffDict());
237 C2_.readIfPresent(this->coeffDict());
238 sigmak_.readIfPresent(this->coeffDict());
239 sigmaEps_.readIfPresent(this->coeffDict());
240
241 return true;
242 }
243
244 return false;
245}
246
247
248template<class BasicTurbulenceModel>
250{
251 if (!this->turbulence_)
252 {
253 return;
254 }
255
256 // Local references
257 const alphaField& alpha = this->alpha_;
258 const rhoField& rho = this->rho_;
259 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
260 const volVectorField& U = this->U_;
261 volScalarField& nut = this->nut_;
263
265
267
269 volScalarField S2(2*magSqr(dev(symm(tgradU()))));
270 volScalarField magS(sqrt(S2));
271
272 volScalarField eta(magS*k_/epsilon_);
273 volScalarField C1(max(eta/(scalar(5) + eta), scalar(0.43)));
274
275 volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
276
277 // Update epsilon and G at the wall
278 epsilon_.boundaryFieldRef().updateCoeffs();
279
280 // SAF: limiting thermo->nu(). If psiThermo is used rho might be < 0
281 // temporarily when p < 0 then nu < 0 which needs limiting
282 volScalarField nuLimited
283 (
284 max
285 (
286 this->nu(),
287 dimensionedScalar(this->nu()().dimensions(), Zero)
288 )
289 );
290
291 // Dissipation equation
293 (
294 fvm::ddt(alpha, rho, epsilon_)
295 + fvm::div(alphaRhoPhi, epsilon_)
296 - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
297 ==
298 C1*alpha*rho*magS*epsilon_
299 - fvm::Sp
300 (
301 C2_*alpha*rho*epsilon_/(k_ + sqrt(nuLimited*epsilon_)),
302 epsilon_
303 )
304 + epsilonSource()
305 + fvOptions(alpha, rho, epsilon_)
306 );
307
308 epsEqn.ref().relax();
309 fvOptions.constrain(epsEqn.ref());
310 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
311 solve(epsEqn);
312 fvOptions.correct(epsilon_);
313 bound(epsilon_, this->epsilonMin_);
314
315
316 // Turbulent kinetic energy equation
317
319 (
320 fvm::ddt(alpha, rho, k_)
321 + fvm::div(alphaRhoPhi, k_)
322 - fvm::laplacian(alpha*rho*DkEff(), k_)
323 ==
324 alpha*rho*G
325 - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
326 - fvm::Sp(alpha*rho*epsilon_/k_, k_)
327 + kSource()
328 + fvOptions(alpha, rho, k_)
329 );
330
331 kEqn.ref().relax();
332 fvOptions.constrain(kEqn.ref());
333 solve(kEqn);
334 fvOptions.correct(k_);
335 bound(k_, this->kMin_);
336
337 correctNut(tgradU(), S2, magS);
338}
339
340
341// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
342
343} // End namespace RASModels
344} // End namespace Foam
345
346// ************************************************************************* //
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
Realizable k-epsilon turbulence model for incompressible and compressible flows.
Definition: realizableKE.H:86
BasicTurbulenceModel::alphaField alphaField
Definition: realizableKE.H:129
tmp< volScalarField > rCmu(const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS)
Definition: realizableKE.C:44
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: realizableKE.C:117
BasicTurbulenceModel::rhoField rhoField
Definition: realizableKE.H:130
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: realizableKE.C:249
virtual tmp< fvScalarMatrix > kSource() const
Definition: realizableKE.C:102
BasicTurbulenceModel::transportModel transportModel
Definition: realizableKE.H:131
virtual bool read()
Re-read model coefficients if they have changed.
Definition: realizableKE.C:232
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
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)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor symm(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)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedTensor skew(const dimensionedTensor &dt)
dimensionedScalar acos(const dimensionedScalar &ds)
volScalarField & nu
volScalarField & alpha
CEqn solve()
edgeScalarField phis(IOobject("phis", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), linearEdgeInterpolate(Us) &aMesh.Le())
Us