Loading...
Searching...
No Matches
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-2023 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#include "wallDist.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
38namespace RASModels
39{
40
41// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42
43template<class BasicTurbulenceModel>
45{
46 this->nut_ = Cmu_*sqr(k_)/epsilon_;
47
49 {
50 const volScalarField& y = wallDist::New(this->mesh_).y();
51
52 // Using the TKE of previous time-step to freeze the turbulent Reynolds
53 // number during the current time-step and hence freeze the separation
54 // of inner and outer layer.
55 const volScalarField Rey(mag(y*sqrt(k_.oldTime())/this->nu()));
56
57 const dimensionedScalar ClStar(0.41*pow(Cmu_, -0.75));
58 auto& lEps = *lEpsPtr_;
59 lEps = y*ClStar*(1.0 - exp(-1.0*Rey/(2*ClStar)));
60
62 auto& lambdaEps = *lambdaEpsPtr_;
63 lambdaEps = 0.5*(1.0 + tanh((Rey - ReyStar_)/A));
64
65 // High-Re part
66 this->nut_ *= lambdaEps;
67
68 // Low-Re part
69 const scalar Amu = 70.0;
70 const volScalarField lMu(y*ClStar*(1.0 - exp(-1.0*Rey/Amu)));
71 this->nut_ += (1.0-lambdaEps)*Cmu_*lMu*sqrt(k_);
72 }
73
74 this->nut_.correctBoundaryConditions();
75 fv::options::New(this->mesh_).correct(this->nut_);
76
77 BasicTurbulenceModel::correctNut();
78}
79
80
81template<class BasicTurbulenceModel>
83{
84 // Source term for k equation (added to RHS)
85
86 if (twoLayerTreatment_)
87 {
88 // Local references
89 const alphaField& alpha = this->alpha_;
90 const rhoField& rho = this->rho_;
91 const volScalarField::Internal& lEps = *lEpsPtr_;
92 const volScalarField& lambdaEps = *lambdaEpsPtr_;
93
94 return
95 -fvm::Sp
96 (
97 alpha()*rho()*(1.0-lambdaEps())
98 *(sqrt(k_())/lEps - epsilon_()/k_()),
99 k_
100 );
101 }
102
104 (
106 dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
107 );
108}
109
110
111template<class BasicTurbulenceModel>
113{
115 (
116 epsilon_,
117 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()/dimTime
118 );
119}
120
121
122// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123
124template<class BasicTurbulenceModel>
125kEpsilon<BasicTurbulenceModel>::kEpsilon
126(
127 const alphaField& alpha,
128 const rhoField& rho,
129 const volVectorField& U,
130 const surfaceScalarField& alphaRhoPhi,
131 const surfaceScalarField& phi,
132 const transportModel& transport,
133 const word& propertiesName,
134 const word& type
135)
136:
137 eddyViscosity<RASModel<BasicTurbulenceModel>>
138 (
139 type,
140 alpha,
141 rho,
142 U,
143 alphaRhoPhi,
144 phi,
145 transport,
146 propertiesName
147 ),
148
149 Cmu_
150 (
151 dimensioned<scalar>::getOrAddToDict
153 "Cmu",
154 this->coeffDict_,
155 0.09
156 )
157 ),
158 C1_
159 (
160 dimensioned<scalar>::getOrAddToDict
161 (
162 "C1",
163 this->coeffDict_,
164 1.44
165 )
166 ),
167 C2_
168 (
169 dimensioned<scalar>::getOrAddToDict
170 (
171 "C2",
172 this->coeffDict_,
173 1.92
174 )
175 ),
176 C3_
177 (
178 dimensioned<scalar>::getOrAddToDict
179 (
180 "C3",
181 this->coeffDict_,
182 0
183 )
184 ),
185 sigmak_
186 (
187 dimensioned<scalar>::getOrAddToDict
188 (
189 "sigmak",
190 this->coeffDict_,
191 1.0
192 )
193 ),
195 (
196 dimensioned<scalar>::getOrAddToDict
197 (
198 "sigmaEps",
199 this->coeffDict_,
200 1.3
201 )
202 ),
203 k_
204 (
206 (
208 this->runTime_.timeName(),
209 this->mesh_,
212 ),
213 this->mesh_
214 ),
216 (
218 (
219 IOobject::groupName("epsilon", alphaRhoPhi.group()),
220 this->runTime_.timeName(),
221 this->mesh_,
224 ),
225 this->mesh_
226 ),
228 (
229 Switch::getOrAddToDict
230 (
231 "twoLayerTreatment",
232 this->coeffDict_,
233 false
234 )
235 ),
237 (
238 dimensioned<scalar>::getOrAddToDict
239 (
240 "ReyStar",
241 this->coeffDict_,
242 200.0
244 ),
246 (
247 dimensioned<scalar>::getOrAddToDict
248 (
249 "ReyFactor",
250 this->coeffDict_,
251 0.05
252 )
253 )
254{
256 {
257 Info<< "Two-layer wall treatment activated" << endl;
258
259 lEpsPtr_ = std::make_unique<volScalarField::Internal>
260 (
262 (
263 "lEps",
264 this->runTime_.name(),
265 this->mesh_,
268 ),
269 this->mesh_,
271 );
272
273 lambdaEpsPtr_ = std::make_unique<volScalarField>
274 (
276 (
277 "lambdaEps",
278 this->runTime_.name(),
279 this->mesh_,
282 ),
283 this->mesh_,
284 dimensionedScalar("lambdaEps", dimless, 1.0)
285 );
286 }
287
288 bound(k_, this->kMin_);
289 bound(epsilon_, this->epsilonMin_);
290
291 if (type == typeName)
292 {
293 this->printCoeffs(type);
295}
296
297
298// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
299
300template<class BasicTurbulenceModel>
302{
304 {
305 Cmu_.readIfPresent(this->coeffDict());
306 C1_.readIfPresent(this->coeffDict());
307 C2_.readIfPresent(this->coeffDict());
308 C3_.readIfPresent(this->coeffDict());
309 sigmak_.readIfPresent(this->coeffDict());
310 sigmaEps_.readIfPresent(this->coeffDict());
311
312 // Two-layer wall treatment
313 // Note: these are the only parameters that can be changed at runtime
314 ReyStar_.readIfPresent(this->coeffDict());
315 ReyFactor_.readIfPresent(this->coeffDict());
316
317 return true;
319
320 return false;
321}
322
323
324template<class BasicTurbulenceModel>
326{
327 if (!this->turbulence_)
328 {
329 return;
330 }
331
332 // Local references
333 const alphaField& alpha = this->alpha_;
334 const rhoField& rho = this->rho_;
335 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
336 const volVectorField& U = this->U_;
337 const volScalarField& nut = this->nut_;
338
339 fv::options& fvOptions(fv::options::New(this->mesh_));
340
341 eddyViscosity<RASModel<BasicTurbulenceModel>>::correct();
342
344 (
345 fvc::div(fvc::absolute(this->phi(), U))().v()
346 );
347
348 tmp<volTensorField> tgradU = fvc::grad(U);
349 const volScalarField::Internal GbyNu
350 (
351 IOobject::scopedName(this->type(), "GbyNu"),
352 tgradU().v() && devTwoSymm(tgradU().v())
353 );
354 const volScalarField::Internal G(this->GName(), nut()*GbyNu);
355 tgradU.clear();
356
357 // Update epsilon and G at the wall
358 epsilon_.boundaryFieldRef().updateCoeffs();
359 // Push any changed cell values to coupled neighbours
360 epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
361
362 // Dissipation equation
363 tmp<fvScalarMatrix> epsEqn
364 (
365 fvm::ddt(alpha, rho, epsilon_)
366 + fvm::div(alphaRhoPhi, epsilon_)
367 - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
368 ==
369 C1_*alpha()*rho()*GbyNu*Cmu_*k_()
370 - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
371 - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
372 + epsilonSource()
373 + fvOptions(alpha, rho, epsilon_)
374 );
375
376 epsEqn.ref().relax();
377 fvOptions.constrain(epsEqn.ref());
378 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
379 solve(epsEqn);
380 fvOptions.correct(epsilon_);
381 bound(epsilon_, this->epsilonMin_);
382
383 // Turbulent kinetic energy equation
384 tmp<fvScalarMatrix> kEqn
385 (
386 fvm::ddt(alpha, rho, k_)
387 + fvm::div(alphaRhoPhi, k_)
388 - fvm::laplacian(alpha*rho*DkEff(), k_)
389 ==
390 alpha()*rho()*G
391 - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
392 - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
393 + kSource()
394 + fvOptions(alpha, rho, k_)
395 );
396
397 kEqn.ref().relax();
398 fvOptions.constrain(kEqn.ref());
399 solve(kEqn);
400 fvOptions.correct(k_);
401 bound(k_, this->kMin_);
402
403 correctNut();
404}
405
406
407// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
408
409} // End namespace RASModels
410} // End namespace Foam
411
412// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
scalar y
Bound the given scalar field if it has gone unbounded.
fv::options & fvOptions
DimensionedField< scalar, volMesh > Internal
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
static word group(const word &name)
Return group (extension part of name).
Definition IOobject.C:300
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
static FOAM_NO_DANGLING_REFERENCE const wallDist & New(const fvMesh &mesh, Args &&... args)
Templated abstract base class for RAS turbulence models.
Definition RASModel.H:81
volScalarField epsilon_
Definition kEpsilon.H:137
BasicTurbulenceModel::alphaField alphaField
Definition kEpsilon.H:158
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition kEpsilon.C:105
BasicTurbulenceModel::rhoField rhoField
Definition kEpsilon.H:159
std::unique_ptr< volScalarField > lambdaEpsPtr_
Definition kEpsilon.H:146
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition kEpsilon.C:318
dimensionedScalar ReyStar_
Definition kEpsilon.H:143
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition kEpsilon.H:215
std::unique_ptr< volScalarField::Internal > lEpsPtr_
Definition kEpsilon.H:145
dimensionedScalar ReyFactor_
Definition kEpsilon.H:144
virtual void correctNut()
Definition kEpsilon.C:37
virtual tmp< fvScalarMatrix > kSource() const
Definition kEpsilon.C:75
BasicTurbulenceModel::transportModel transportModel
Definition kEpsilon.H:160
dimensionedScalar Cmu_
Definition kEpsilon.H:126
virtual bool read()
Re-read model coefficients if they have changed.
Definition kEpsilon.C:294
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition kEpsilon.H:203
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
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)
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Base-class for all transport models used by the incompressible turbulence models.
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
const volScalarField & y() const noexcept
Return reference to cached distance-to-wall field.
Definition wallDist.H:220
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
thermo correct()
scalar Rey
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:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition fvcMeshPhi.C:183
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
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.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition bound.C:29
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar atanh(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & alpha
CEqn solve()