Loading...
Searching...
No Matches
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) 2022 Upstream CFD GmbH
10 Copyright (C) 2016-2023 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "kOmegaSSTBase.H"
31#include "fvOptions.H"
32#include "bound.H"
33#include "wallDist.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37namespace Foam
38{
39
40// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
41
42template<class BasicEddyViscosityModel>
44(
45 const volScalarField& CDkOmega
46) const
47{
48 tmp<volScalarField> CDkOmegaPlus = max
49 (
50 CDkOmega,
52 );
53
55 (
56 min
57 (
58 max
59 (
60 (scalar(1)/betaStar_)*sqrt(k_)/(omega_*y_),
61 scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
62 ),
63 (4*alphaOmega2_)*k_/(CDkOmegaPlus*sqr(y_))
64 ),
65 scalar(10)
66 );
67
68 return tanh(pow4(arg1));
69}
70
71
72template<class BasicEddyViscosityModel>
74{
76 (
77 max
78 (
79 (scalar(2)/betaStar_)*sqrt(k_)/(omega_*y_),
80 scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
81 ),
82 scalar(100)
83 );
84
85 return tanh(sqr(arg2));
86}
87
88
89template<class BasicEddyViscosityModel>
91{
93 (
94 150*(this->mu()/this->rho_)/(omega_*sqr(y_)),
95 scalar(10)
96 );
97
98 return 1 - tanh(pow4(arg3));
99}
100
101
102template<class BasicEddyViscosityModel>
104{
105 tmp<volScalarField> f23(F2());
106
107 if (F3_)
108 {
109 f23.ref() *= F3();
111
112 return f23;
113}
114
115
116template<class BasicEddyViscosityModel>
118(
119 const volScalarField& S2
120)
121{
122 // Correct the turbulence viscosity
123 this->nut_ = a1_*k_/max(a1_*omega_, b1_*F23()*sqrt(S2));
124 this->nut_.correctBoundaryConditions();
125 fv::options::New(this->mesh_).correct(this->nut_);
126}
127
128
129template<class BasicEddyViscosityModel>
131{
132 correctNut(2*magSqr(symm(fvc::grad(this->U_))));
133}
134
135
136template<class BasicEddyViscosityModel>
138(
139 const volTensorField& gradU
140) const
141{
142 return 2*magSqr(symm(gradU));
143}
144
145
146template<class BasicEddyViscosityModel>
148(
150) const
151{
152 return min(G, (c1_*betaStar_)*this->k_()*this->omega_());
153}
154
155
156template<class BasicEddyViscosityModel>
158(
159 const volScalarField& /* F1 not used */,
160 const volTensorField& /* gradU not used */
161) const
162{
163 return betaStar_*omega_();
164}
165
166
167template<class BasicEddyViscosityModel>
169(
170 const volTensorField& gradU,
171 const volScalarField& /* S2 not used */
172) const
173{
175 (
176 IOobject::scopedName(this->type(), "GbyNu"),
177 gradU() && devTwoSymm(gradU())
178 );
179}
180
181
182template<class BasicEddyViscosityModel>
184(
185 const volScalarField::Internal& GbyNu0,
188) const
189{
190 return min
191 (
194 );
196
198template<class BasicEddyViscosityModel>
202 (
204 dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
205 );
206}
207
208
209template<class BasicEddyViscosityModel>
211{
213 (
215 dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
216 );
217}
218
219
220template<class BasicEddyViscosityModel>
222(
223 const volScalarField::Internal& S2,
226) const
227{
229 (
230 omega_,
231 dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
232 );
233}
234
235
236// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
237
238template<class BasicEddyViscosityModel>
240(
241 const word& type,
242 const alphaField& alpha,
243 const rhoField& rho,
244 const volVectorField& U,
245 const surfaceScalarField& alphaRhoPhi,
246 const surfaceScalarField& phi,
247 const transportModel& transport,
248 const word& propertiesName
249)
250:
251 BasicEddyViscosityModel
252 (
253 type,
254 alpha,
255 rho,
256 U,
257 alphaRhoPhi,
258 phi,
259 transport,
260 propertiesName
261 ),
264 (
265 dimensioned<scalar>::getOrAddToDict
266 (
267 "alphaK1",
268 this->coeffDict_,
269 0.85
270 )
271 ),
273 (
274 dimensioned<scalar>::getOrAddToDict
275 (
276 "alphaK2",
277 this->coeffDict_,
278 1.0
279 )
280 ),
282 (
283 dimensioned<scalar>::getOrAddToDict
284 (
285 "alphaOmega1",
286 this->coeffDict_,
287 0.5
288 )
289 ),
291 (
292 dimensioned<scalar>::getOrAddToDict
293 (
294 "alphaOmega2",
295 this->coeffDict_,
296 0.856
297 )
298 ),
299 gamma1_
300 (
301 dimensioned<scalar>::getOrAddToDict
302 (
303 "gamma1",
304 this->coeffDict_,
305 5.0/9.0
306 )
307 ),
308 gamma2_
309 (
310 dimensioned<scalar>::getOrAddToDict
311 (
312 "gamma2",
313 this->coeffDict_,
314 0.44
315 )
316 ),
317 beta1_
318 (
319 dimensioned<scalar>::getOrAddToDict
320 (
321 "beta1",
322 this->coeffDict_,
323 0.075
324 )
325 ),
326 beta2_
327 (
328 dimensioned<scalar>::getOrAddToDict
329 (
330 "beta2",
331 this->coeffDict_,
332 0.0828
333 )
334 ),
336 (
337 dimensioned<scalar>::getOrAddToDict
338 (
339 "betaStar",
340 this->coeffDict_,
341 0.09
342 )
343 ),
344 a1_
345 (
346 dimensioned<scalar>::getOrAddToDict
347 (
348 "a1",
349 this->coeffDict_,
350 0.31
351 )
352 ),
353 b1_
354 (
355 dimensioned<scalar>::getOrAddToDict
356 (
357 "b1",
358 this->coeffDict_,
359 1.0
360 )
361 ),
362 c1_
363 (
364 dimensioned<scalar>::getOrAddToDict
365 (
366 "c1",
367 this->coeffDict_,
368 10.0
369 )
370 ),
372 (
373 Switch::getOrAddToDict
374 (
375 "F3",
376 this->coeffDict_,
377 false
378 )
379 ),
380
381 y_(wallDist::New(this->mesh_).y()),
382
383 k_
384 (
386 (
387 IOobject::groupName("k", alphaRhoPhi.group()),
388 this->runTime_.timeName(),
389 this->mesh_,
390 IOobject::MUST_READ,
391 IOobject::AUTO_WRITE
392 ),
393 this->mesh_
394 ),
395 omega_
396 (
398 (
399 IOobject::groupName("omega", alphaRhoPhi.group()),
400 this->runTime_.timeName(),
401 this->mesh_,
402 IOobject::MUST_READ,
403 IOobject::AUTO_WRITE
404 ),
405 this->mesh_
406 ),
407
409 (
410 Switch::getOrAddToDict
411 (
412 "decayControl",
413 this->coeffDict_,
414 false
415 )
416 ),
417 kInf_
418 (
419 dimensioned<scalar>::getOrAddToDict
420 (
421 "kInf",
422 this->coeffDict_,
423 k_.dimensions(),
424 0
425 )
426 ),
428 (
429 dimensioned<scalar>::getOrAddToDict
430 (
431 "omegaInf",
432 this->coeffDict_,
433 omega_.dimensions(),
434 0
435 )
436 )
437{
438 bound(k_, this->kMin_);
439 bound(omega_, this->omegaMin_);
440
441 setDecayControl(this->coeffDict_);
442}
443
444
445// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
446
447template<class BasicEddyViscosityModel>
449(
450 const dictionary& dict
451)
452{
453 decayControl_.readIfPresent("decayControl", dict);
454
455 if (decayControl_)
456 {
457 kInf_.read(dict);
458 omegaInf_.read(dict);
459
460 Info<< " Employing decay control with kInf:" << kInf_
461 << " and omegaInf:" << omegaInf_ << endl;
462 }
463 else
464 {
465 kInf_.value() = 0;
466 omegaInf_.value() = 0;
467 }
468}
469
470
471template<class BasicEddyViscosityModel>
473{
474 if (BasicEddyViscosityModel::read())
475 {
476 alphaK1_.readIfPresent(this->coeffDict());
477 alphaK2_.readIfPresent(this->coeffDict());
478 alphaOmega1_.readIfPresent(this->coeffDict());
479 alphaOmega2_.readIfPresent(this->coeffDict());
480 gamma1_.readIfPresent(this->coeffDict());
481 gamma2_.readIfPresent(this->coeffDict());
482 beta1_.readIfPresent(this->coeffDict());
483 beta2_.readIfPresent(this->coeffDict());
484 betaStar_.readIfPresent(this->coeffDict());
485 a1_.readIfPresent(this->coeffDict());
486 b1_.readIfPresent(this->coeffDict());
487 c1_.readIfPresent(this->coeffDict());
488 F3_.readIfPresent("F3", this->coeffDict());
489
490 setDecayControl(this->coeffDict());
491
492 return true;
494
495 return false;
496}
497
498
499template<class BasicEddyViscosityModel>
501{
502 if (!this->turbulence_)
503 {
504 return;
505 }
506
507 // Local references
508 const alphaField& alpha = this->alpha_;
509 const rhoField& rho = this->rho_;
510 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
511 const volVectorField& U = this->U_;
512 volScalarField& nut = this->nut_;
514
515 BasicEddyViscosityModel::correct();
516
518 (
519 fvc::div(fvc::absolute(this->phi(), U))
520 );
521
523 const volScalarField S2(this->S2(tgradU()));
524 volScalarField::Internal GbyNu0(this->GbyNu0(tgradU(), S2));
525 volScalarField::Internal G(this->GName(), nut*GbyNu0);
526
527
528 // - boundary condition changes a cell value
529 // - normally this would be triggered through correctBoundaryConditions
530 // - which would do
531 // - fvPatchField::evaluate() which calls
532 // - fvPatchField::updateCoeffs()
533 // - however any processor boundary conditions already start sending
534 // at initEvaluate so would send over the old value.
535 // - avoid this by explicitly calling updateCoeffs early and then
536 // only doing the boundary conditions that rely on initEvaluate
537 // (currently only coupled ones)
538
539 //- 1. Explicitly swap values on coupled boundary conditions
540 // Update omega and G at the wall
541 omega_.boundaryFieldRef().updateCoeffs();
542 // omegaWallFunctions change the cell value! Make sure to push these to
543 // coupled neighbours. Note that we want to avoid the re-updateCoeffs
544 // of the wallFunctions so make sure to bypass the evaluate on
545 // those patches and only do the coupled ones.
546 omega_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
547
552 //omega_.correctBoundaryConditions();
553
554
555 const volScalarField CDkOmega
556 (
557 (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_
558 );
559
560 const volScalarField F1(this->F1(CDkOmega));
561 const volScalarField F23(this->F23());
562
563 {
565 const volScalarField::Internal beta(this->beta(F1));
566
567 GbyNu0 = GbyNu(GbyNu0, F23(), S2());
568
569 // Turbulent frequency equation
570 tmp<fvScalarMatrix> omegaEqn
571 (
572 fvm::ddt(alpha, rho, omega_)
573 + fvm::div(alphaRhoPhi, omega_)
574 - fvm::laplacian(alpha*rho*DomegaEff(F1), omega_)
575 ==
576 alpha()*rho()*gamma*GbyNu0
577 - fvm::SuSp((2.0/3.0)*alpha()*rho()*gamma*divU, omega_)
578 - fvm::Sp(alpha()*rho()*beta*omega_(), omega_)
579 - fvm::SuSp
580 (
581 alpha()*rho()*(F1() - scalar(1))*CDkOmega()/omega_(),
582 omega_
583 )
584 + alpha()*rho()*beta*sqr(omegaInf_)
585 + Qsas(S2(), gamma, beta)
586 + omegaSource()
587 + fvOptions(alpha, rho, omega_)
588 );
589
590 omegaEqn.ref().relax();
591 fvOptions.constrain(omegaEqn.ref());
592 omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
593 solve(omegaEqn);
594 fvOptions.correct(omega_);
595 bound(omega_, this->omegaMin_);
596 }
597
598 {
599 // Turbulent kinetic energy equation
601 (
602 fvm::ddt(alpha, rho, k_)
603 + fvm::div(alphaRhoPhi, k_)
604 - fvm::laplacian(alpha*rho*DkEff(F1), k_)
605 ==
606 alpha()*rho()*Pk(G)
607 - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
608 - fvm::Sp(alpha()*rho()*epsilonByk(F1, tgradU()), k_)
609 + alpha()*rho()*betaStar_*omegaInf_*kInf_
610 + kSource()
611 + fvOptions(alpha, rho, k_)
612 );
613
614 tgradU.clear();
615
616 kEqn.ref().relax();
617 fvOptions.constrain(kEqn.ref());
618 solve(kEqn);
619 fvOptions.correct(k_);
620 bound(k_, this->kMin_);
621 }
622
623 correctNut(S2);
624}
625
626
627// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
628
629} // End namespace Foam
630
631// ************************************************************************* //
scalar y
#define F2(B, C, D)
Definition SHA1.C:150
#define F1(B, C, D)
Definition SHA1.C:149
#define F3(B, C, D)
Definition SHA1.C:151
Bound the given scalar field if it has gone unbounded.
fv::options & fvOptions
DimensionedField< scalar, volMesh > Internal
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
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Generic dimensioned Type class.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
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
volScalarField k_
Turbulent kinetic energy field [m^2/s^2].
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
Switch F3_
Flag to include the F3 term.
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
kOmegaSSTBase(const kOmegaSSTBase &)=delete
No copy construct.
virtual void correctNut(const volScalarField &S2)
dimensionedScalar alphaOmega2_
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< volScalarField > F2() const
virtual tmp< volScalarField > F23() const
BasicEddyViscosityModel::transportModel transportModel
dimensionedScalar b1_
virtual tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Return G/nu.
dimensionedScalar c1_
const volScalarField & y_
Wall distance.
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Return the effective diffusivity for omega.
dimensionedScalar a1_
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
virtual tmp< volScalarField::Internal > GbyNu0(const volTensorField &gradU, const volScalarField &S2) const
Return (G/nu)_0.
dimensionedScalar alphaK1_
volScalarField omega_
Specific dissipation rate field [1/s].
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.
virtual tmp< volScalarField > S2(const volTensorField &gradU) const
Return square of strain rate.
BasicEddyViscosityModel::rhoField rhoField
BasicEddyViscosityModel::alphaField alphaField
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.
dimensionedScalar betaStar_
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< volScalarField > F3() const
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.
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
Definition wallDist.H:74
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
const scalar gamma
Definition EEqn.H:9
scalar nut
zeroField divU
Definition alphaSuSp.H:3
word timeName
Definition getTimeIndex.H:3
const dimensionedScalar G
Newtonian constant of gravitation.
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.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
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: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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
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:26
dimensionedScalar pow4(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
volScalarField & alpha
dictionary dict
CEqn solve()
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)