kOmegaSSTBase.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-2015 OpenFOAM Foundation
9 Copyright (C) 2016-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 "kOmegaSSTBase.H"
30#include "fvOptions.H"
31#include "bound.H"
32#include "wallDist.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
40
41template<class BasicEddyViscosityModel>
43(
44 const volScalarField& CDkOmega
45) const
46{
47 tmp<volScalarField> CDkOmegaPlus = max
48 (
49 CDkOmega,
50 dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
51 );
52
54 (
55 min
56 (
57 max
58 (
59 (scalar(1)/betaStar_)*sqrt(k_)/(omega_*y_),
60 scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
61 ),
62 (4*alphaOmega2_)*k_/(CDkOmegaPlus*sqr(y_))
63 ),
64 scalar(10)
65 );
66
67 return tanh(pow4(arg1));
68}
69
70
71template<class BasicEddyViscosityModel>
73{
75 (
76 max
77 (
78 (scalar(2)/betaStar_)*sqrt(k_)/(omega_*y_),
79 scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
80 ),
81 scalar(100)
82 );
83
84 return tanh(sqr(arg2));
85}
86
87
88template<class BasicEddyViscosityModel>
90{
92 (
93 150*(this->mu()/this->rho_)/(omega_*sqr(y_)),
94 scalar(10)
95 );
96
97 return 1 - tanh(pow4(arg3));
98}
99
100
101template<class BasicEddyViscosityModel>
103{
104 tmp<volScalarField> f23(F2());
105
106 if (F3_)
107 {
108 f23.ref() *= F3();
109 }
110
111 return f23;
112}
113
114
115template<class BasicEddyViscosityModel>
117(
118 const volScalarField& S2
119)
120{
121 // Correct the turbulence viscosity
122 this->nut_ = a1_*k_/max(a1_*omega_, b1_*F23()*sqrt(S2));
123 this->nut_.correctBoundaryConditions();
124 fv::options::New(this->mesh_).correct(this->nut_);
125}
126
127
128template<class BasicEddyViscosityModel>
130{
131 correctNut(2*magSqr(symm(fvc::grad(this->U_))));
132}
133
134
135template<class BasicEddyViscosityModel>
137(
139) const
140{
141 return min(G, (c1_*betaStar_)*this->k_()*this->omega_());
142}
143
144
145template<class BasicEddyViscosityModel>
148(
149 const volScalarField& F1,
150 const volTensorField& gradU
151) const
152{
153 return betaStar_*omega_();
154}
155
156
157template<class BasicEddyViscosityModel>
159(
160 const volScalarField::Internal& GbyNu0,
163) const
164{
165 return min
166 (
167 GbyNu0,
168 (c1_/a1_)*betaStar_*omega_()
169 *max(a1_*omega_(), b1_*F2*sqrt(S2))
170 );
171}
172
173
174template<class BasicEddyViscosityModel>
176{
178 (
180 (
181 k_,
182 dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
183 )
184 );
185}
186
187
188template<class BasicEddyViscosityModel>
195 omega_,
196 dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
197 )
198 );
199}
200
201
202template<class BasicEddyViscosityModel>
204(
205 const volScalarField::Internal& S2,
208) const
209{
211 (
213 (
214 omega_,
215 dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
216 )
217 );
218}
219
220
221// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
222
223template<class BasicEddyViscosityModel>
225(
226 const word& type,
227 const alphaField& alpha,
228 const rhoField& rho,
229 const volVectorField& U,
230 const surfaceScalarField& alphaRhoPhi,
231 const surfaceScalarField& phi,
232 const transportModel& transport,
233 const word& propertiesName
234)
235:
236 BasicEddyViscosityModel
237 (
238 type,
239 alpha,
240 rho,
241 U,
242 alphaRhoPhi,
243 phi,
244 transport,
245 propertiesName
246 ),
247
248 alphaK1_
249 (
250 dimensioned<scalar>::getOrAddToDict
252 "alphaK1",
253 this->coeffDict_,
254 0.85
255 )
256 ),
257 alphaK2_
258 (
259 dimensioned<scalar>::getOrAddToDict
260 (
261 "alphaK2",
262 this->coeffDict_,
263 1.0
264 )
265 ),
266 alphaOmega1_
267 (
268 dimensioned<scalar>::getOrAddToDict
270 "alphaOmega1",
271 this->coeffDict_,
272 0.5
273 )
274 ),
275 alphaOmega2_
277 dimensioned<scalar>::getOrAddToDict
279 "alphaOmega2",
280 this->coeffDict_,
281 0.856
282 )
283 ),
284 gamma1_
285 (
286 dimensioned<scalar>::getOrAddToDict
287 (
288 "gamma1",
289 this->coeffDict_,
290 5.0/9.0
291 )
292 ),
293 gamma2_
294 (
295 dimensioned<scalar>::getOrAddToDict
296 (
297 "gamma2",
298 this->coeffDict_,
299 0.44
300 )
301 ),
302 beta1_
303 (
304 dimensioned<scalar>::getOrAddToDict
305 (
306 "beta1",
307 this->coeffDict_,
308 0.075
309 )
310 ),
311 beta2_
312 (
313 dimensioned<scalar>::getOrAddToDict
314 (
315 "beta2",
316 this->coeffDict_,
317 0.0828
319 ),
320 betaStar_
321 (
322 dimensioned<scalar>::getOrAddToDict
323 (
324 "betaStar",
325 this->coeffDict_,
326 0.09
327 )
328 ),
329 a1_
330 (
331 dimensioned<scalar>::getOrAddToDict
332 (
333 "a1",
334 this->coeffDict_,
335 0.31
336 )
337 ),
338 b1_
339 (
340 dimensioned<scalar>::getOrAddToDict
341 (
342 "b1",
343 this->coeffDict_,
344 1.0
345 )
346 ),
347 c1_
348 (
349 dimensioned<scalar>::getOrAddToDict
350 (
351 "c1",
352 this->coeffDict_,
353 10.0
354 )
355 ),
356 F3_
357 (
358 Switch::getOrAddToDict
359 (
360 "F3",
361 this->coeffDict_,
362 false
363 )
364 ),
365
366 y_(wallDist::New(this->mesh_).y()),
367
368 k_
369 (
371 (
372 IOobject::groupName("k", alphaRhoPhi.group()),
373 this->runTime_.timeName(),
374 this->mesh_,
375 IOobject::MUST_READ,
376 IOobject::AUTO_WRITE
377 ),
378 this->mesh_
379 ),
380 omega_
381 (
383 (
384 IOobject::groupName("omega", alphaRhoPhi.group()),
385 this->runTime_.timeName(),
386 this->mesh_,
387 IOobject::MUST_READ,
388 IOobject::AUTO_WRITE
389 ),
390 this->mesh_
391 ),
392 decayControl_
393 (
394 Switch::getOrAddToDict
395 (
396 "decayControl",
397 this->coeffDict_,
398 false
399 )
400 ),
401 kInf_
402 (
403 dimensioned<scalar>::getOrAddToDict
404 (
405 "kInf",
406 this->coeffDict_,
407 k_.dimensions(),
408 0
409 )
410 ),
411 omegaInf_
412 (
413 dimensioned<scalar>::getOrAddToDict
414 (
415 "omegaInf",
416 this->coeffDict_,
417 omega_.dimensions(),
418 0
419 )
420 )
421{
422 bound(k_, this->kMin_);
423 bound(omega_, this->omegaMin_);
424
425 setDecayControl(this->coeffDict_);
426}
427
428
429// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
430
431template<class BasicEddyViscosityModel>
433(
434 const dictionary& dict
435)
436{
437 decayControl_.readIfPresent("decayControl", dict);
438
439 if (decayControl_)
440 {
441 kInf_.read(dict);
442 omegaInf_.read(dict);
443
444 Info<< " Employing decay control with kInf:" << kInf_
445 << " and omegaInf:" << omegaInf_ << endl;
446 }
447 else
448 {
449 kInf_.value() = 0;
450 omegaInf_.value() = 0;
451 }
452}
453
454
455template<class BasicEddyViscosityModel>
457{
458 if (BasicEddyViscosityModel::read())
459 {
460 alphaK1_.readIfPresent(this->coeffDict());
461 alphaK2_.readIfPresent(this->coeffDict());
462 alphaOmega1_.readIfPresent(this->coeffDict());
463 alphaOmega2_.readIfPresent(this->coeffDict());
464 gamma1_.readIfPresent(this->coeffDict());
465 gamma2_.readIfPresent(this->coeffDict());
466 beta1_.readIfPresent(this->coeffDict());
467 beta2_.readIfPresent(this->coeffDict());
468 betaStar_.readIfPresent(this->coeffDict());
469 a1_.readIfPresent(this->coeffDict());
470 b1_.readIfPresent(this->coeffDict());
471 c1_.readIfPresent(this->coeffDict());
472 F3_.readIfPresent("F3", this->coeffDict());
473
474 setDecayControl(this->coeffDict());
475
476 return true;
477 }
478
479 return false;
480}
481
482
483template<class BasicEddyViscosityModel>
485{
486 if (!this->turbulence_)
487 {
488 return;
489 }
490
491 // Local references
492 const alphaField& alpha = this->alpha_;
493 const rhoField& rho = this->rho_;
494 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
495 const volVectorField& U = this->U_;
496 volScalarField& nut = this->nut_;
498
499 BasicEddyViscosityModel::correct();
500
502
504 volScalarField S2(2*magSqr(symm(tgradU())));
506 (
507 this->type() + ":GbyNu",
508 (tgradU() && dev(twoSymm(tgradU())))
509 );
510 volScalarField::Internal G(this->GName(), nut*GbyNu0);
511
512 // Update omega and G at the wall
513 omega_.boundaryFieldRef().updateCoeffs();
514
515 volScalarField CDkOmega
516 (
517 (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_
518 );
519
520 volScalarField F1(this->F1(CDkOmega));
521 volScalarField F23(this->F23());
522
523 {
526
527 GbyNu0 = GbyNu(GbyNu0, F23(), S2());
528
529 // Turbulent frequency equation
530 tmp<fvScalarMatrix> omegaEqn
531 (
532 fvm::ddt(alpha, rho, omega_)
533 + fvm::div(alphaRhoPhi, omega_)
534 - fvm::laplacian(alpha*rho*DomegaEff(F1), omega_)
535 ==
536 alpha()*rho()*gamma*GbyNu0
537 - fvm::SuSp((2.0/3.0)*alpha()*rho()*gamma*divU, omega_)
538 - fvm::Sp(alpha()*rho()*beta*omega_(), omega_)
539 - fvm::SuSp
540 (
541 alpha()*rho()*(F1() - scalar(1))*CDkOmega()/omega_(),
542 omega_
543 )
544 + alpha()*rho()*beta*sqr(omegaInf_)
545 + Qsas(S2(), gamma, beta)
546 + omegaSource()
547 + fvOptions(alpha, rho, omega_)
548 );
549
550 omegaEqn.ref().relax();
551 fvOptions.constrain(omegaEqn.ref());
552 omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
553 solve(omegaEqn);
554 fvOptions.correct(omega_);
555 bound(omega_, this->omegaMin_);
556 }
557
558 // Turbulent kinetic energy equation
560 (
561 fvm::ddt(alpha, rho, k_)
562 + fvm::div(alphaRhoPhi, k_)
563 - fvm::laplacian(alpha*rho*DkEff(F1), k_)
564 ==
565 alpha()*rho()*Pk(G)
566 - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
567 - fvm::Sp(alpha()*rho()*epsilonByk(F1, tgradU()), k_)
568 + alpha()*rho()*betaStar_*omegaInf_*kInf_
569 + kSource()
570 + fvOptions(alpha, rho, k_)
571 );
572
573 tgradU.clear();
574
575 kEqn.ref().relax();
576 fvOptions.constrain(kEqn.ref());
577 solve(kEqn);
578 fvOptions.correct(k_);
579 bound(k_, this->kMin_);
580
581 correctNut(S2);
582}
583
584
585// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
586
587} // End namespace Foam
588
589// ************************************************************************* //
scalar y
#define F2(B, C, D)
Definition: SHA1.C:151
#define F1(B, C, D)
Definition: SHA1.C:150
#define F3(B, C, D)
Definition: SHA1.C:152
Bound the given scalar field if it has gone unbounded.
fv::options & fvOptions
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Generic dimensioned Type class.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Finite-volume options.
Definition: fvOptions.H:59
Base class implementation of the k-omega-SST turbulence model for incompressible and compressible flo...
void setDecayControl(const dictionary &dict)
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< volScalarField > F2() const
Definition: kOmegaSSTBase.C:72
virtual tmp< volScalarField > F23() const
BasicEddyViscosityModel::transportModel transportModel
virtual tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Return G/nu.
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
virtual tmp< fvScalarMatrix > omegaSource() const
virtual void correctNut()
virtual tmp< fvScalarMatrix > kSource() const
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Return epsilon/k which for standard RAS is betaStar*omega.
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
BasicEddyViscosityModel::rhoField rhoField
BasicEddyViscosityModel::alphaField alphaField
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< volScalarField > F3() const
Definition: kOmegaSSTBase.C:89
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
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
Definition: wallDist.H:78
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
const volScalarField & mu
const scalar gamma
Definition: EEqn.H:9
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)
const dimensionSet dimless
Dimensionless.
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)
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
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
dimensionedScalar pow4(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
volScalarField & alpha
dictionary dict
CEqn solve()
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)