Loading...
Searching...
No Matches
SpalartAllmarasBase.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-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2022 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 "SpalartAllmarasBase.H"
30#include "wallDist.H"
31#include "bound.H"
32#include "fvOptions.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36namespace Foam
37{
38
39// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40
41template<class BasicEddyViscosityModel>
50(
51 const volScalarField& chi
52) const
54 const volScalarField chi3("chi3", pow3(chi));
55 return chi3/(chi3 + pow3(Cv1_));
56}
57
58
59template<class BasicEddyViscosityModel>
61(
62 const volScalarField& chi,
63 const volScalarField& fv1
64) const
65{
66 return scalar(1) - chi/(scalar(1) + chi*fv1);
67}
68
69
70template<class BasicEddyViscosityModel>
72(
73 const volScalarField& chi
74) const
75{
76 if (ft2_)
77 {
78 return Ct3_*exp(-Ct4_*sqr(chi));
79 }
80
82 (
83 IOobject
84 (
85 "ft2",
86 this->runTime_.timeName(),
87 this->mesh_,
90 ),
91 this->mesh_,
93 );
94}
95
96
97template<class BasicEddyViscosityModel>
99(
100 const volTensorField& gradU
101) const
102{
103 return sqrt(2.0)*mag(skew(gradU));
104}
105
106
107template<class BasicEddyViscosityModel>
109(
110 const volScalarField& nur,
111 const volScalarField& Stilda,
112 const volScalarField& dTilda
113) const
114{
115 const dimensionedScalar eps(Stilda.dimensions(), SMALL);
116
118 min(nur/(max(Stilda, eps)*sqr(kappa_*dTilda)), scalar(10));
119
120 tr.ref().boundaryFieldRef() == 0;
121
122 return tr;
123}
124
125
126template<class BasicEddyViscosityModel>
128(
129 const volScalarField& Stilda,
130 const volScalarField& dTilda
131) const
132{
133 const volScalarField::Internal r(this->r(nuTilda_, Stilda, dTilda)()());
135
136 return g*pow((1 + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)), 1.0/6.0);
137}
139
140template<class BasicEddyViscosityModel>
142(
143 const volScalarField& chi,
145 const volTensorField& gradU,
147) const
149 const volScalarField Omega(this->Omega(gradU));
150
151 return
152 max
153 (
156 );
157}
158
159
160template<class BasicEddyViscosityModel>
162(
163 const volScalarField& fv1
164)
165{
166 this->nut_ = nuTilda_*fv1;
167 this->nut_.correctBoundaryConditions();
168 fv::options::New(this->mesh_).correct(this->nut_);
169}
170
171
172template<class BasicEddyViscosityModel>
174{
175 correctNut(fv1(this->chi()));
176}
177
178
179// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
180
181template<class BasicEddyViscosityModel>
182SpalartAllmarasBase<BasicEddyViscosityModel>::SpalartAllmarasBase
183(
184 const word& type,
185 const alphaField& alpha,
186 const rhoField& rho,
187 const volVectorField& U,
188 const surfaceScalarField& alphaRhoPhi,
189 const surfaceScalarField& phi,
190 const transportModel& transport,
191 const word& propertiesName
192)
193:
194 BasicEddyViscosityModel
195 (
196 type,
197 alpha,
198 rho,
199 U,
200 alphaRhoPhi,
201 phi,
202 transport,
203 propertiesName
204 ),
205
207 (
208 dimensioned<scalar>::getOrAddToDict
209 (
210 "sigmaNut",
211 this->coeffDict_,
212 0.66666
213 )
214 ),
215 kappa_
216 (
217 dimensioned<scalar>::getOrAddToDict
218 (
219 "kappa",
220 this->coeffDict_,
221 0.41
222 )
223 ),
225 (
226 dimensioned<scalar>::getOrAddToDict
227 (
228 "Cb1",
229 this->coeffDict_,
230 0.1355
231 )
232 ),
233 Cb2_
235 dimensioned<scalar>::getOrAddToDict
236 (
237 "Cb2",
238 this->coeffDict_,
239 0.622
240 )
241 ),
242 Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
243 Cw2_
244 (
245 dimensioned<scalar>::getOrAddToDict
246 (
247 "Cw2",
248 this->coeffDict_,
249 0.3
250 )
251 ),
253 (
254 dimensioned<scalar>::getOrAddToDict
255 (
256 "Cw3",
257 this->coeffDict_,
258 2.0
259 )
260 ),
261 Cv1_
262 (
263 dimensioned<scalar>::getOrAddToDict
264 (
265 "Cv1",
266 this->coeffDict_,
267 7.1
268 )
269 ),
270 Cs_
271 (
272 dimensioned<scalar>::getOrAddToDict
273 (
274 "Cs",
275 this->coeffDict_,
276 0.3
277 )
278 ),
279 ck_
280 (
281 dimensioned<scalar>::getOrAddToDict
282 (
283 "ck",
284 this->coeffDict_,
285 0.07
286 )
287 ),
288 ft2_
289 (
290 Switch::getOrAddToDict
291 (
292 "ft2",
293 this->coeffDict_,
294 false
295 )
296 ),
297 Ct3_
298 (
299 dimensioned<scalar>::getOrAddToDict
300 (
301 "Ct3",
302 this->coeffDict_,
303 1.2
304 )
305 ),
306 Ct4_
307 (
308 dimensioned<scalar>::getOrAddToDict
309 (
310 "Ct4",
311 this->coeffDict_,
312 0.5
313 )
314 ),
315
317 (
319 (
320 "nuTilda",
321 this->runTime_.timeName(),
322 this->mesh_,
323 IOobject::MUST_READ,
324 IOobject::AUTO_WRITE
325 ),
326 this->mesh_
327 ),
328
329 y_(wallDist::New(this->mesh_).y())
330{
331 if (ft2_)
332 {
333 Info<< "ft2 term: active" << nl;
334 }
335 else
336 {
337 Info<< "ft2 term: inactive" << nl;
339}
340
341
342// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
343
344template<class BasicEddyViscosityModel>
346{
347 if (BasicEddyViscosityModel::read())
348 {
349 sigmaNut_.readIfPresent(this->coeffDict());
350 kappa_.readIfPresent(this->coeffDict());
351
352 Cb1_.readIfPresent(this->coeffDict());
353 Cb2_.readIfPresent(this->coeffDict());
354 Cw1_ = Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
355 Cw2_.readIfPresent(this->coeffDict());
356 Cw3_.readIfPresent(this->coeffDict());
357 Cv1_.readIfPresent(this->coeffDict());
358 Cs_.readIfPresent(this->coeffDict());
359 ck_.readIfPresent(this->coeffDict());
360
361 ft2_.readIfPresent("ft2", this->coeffDict());
362 Ct3_.readIfPresent(this->coeffDict());
363 Ct4_.readIfPresent(this->coeffDict());
364
365 if (ft2_)
366 {
367 Info<< " ft2 term: active" << nl;
368 }
369 else
370 {
371 Info<< " ft2 term: inactive" << nl;
372 }
373
374 return true;
375 }
377 return false;
378}
379
380
381template<class BasicEddyViscosityModel>
384{
386 (
387 IOobject::groupName("DnuTildaEff", this->alphaRhoPhi_.group()),
388 (nuTilda_ + this->nu())/sigmaNut_
389 );
390}
391
392
393template<class BasicEddyViscosityModel>
395{
396 // (B:Eq. 4.50)
397 const scalar Cmu = 0.09;
398 const auto fv1 = this->fv1(chi());
399
401 (
402 IOobject::groupName("k", this->alphaRhoPhi_.group()),
403 cbrt(fv1)*nuTilda_*::sqrt(scalar(2)/Cmu)*mag(symm(fvc::grad(this->U_)))
404 );
405}
406
407template<class BasicEddyViscosityModel>
410{
411 // (B:Eq. 4.50)
412 const scalar Cmu = 0.09;
413 const auto fv1 = this->fv1(chi());
414 const dimensionedScalar nutSMALL(nuTilda_.dimensions(), SMALL);
415
417 (
418 IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
419 sqrt(fv1)*sqr(::sqrt(Cmu)*this->k())/(nuTilda_ + this->nut_ + nutSMALL)
420 );
421}
422
423
424template<class BasicEddyViscosityModel>
426{
427 // (P:p. 384)
428 const scalar betaStar = 0.09;
429 const dimensionedScalar k0(sqr(dimLength/dimTime), SMALL);
430
432 (
433 IOobject::groupName("omega", this->alphaRhoPhi_.group()),
434 this->epsilon()/(betaStar*(this->k() + k0))
435 );
436}
437
438
439template<class BasicEddyViscosityModel>
441{
442 if (!this->turbulence_)
443 {
444 return;
445 }
446
447 {
448 // Local references
449 const alphaField& alpha = this->alpha_;
450 const rhoField& rho = this->rho_;
451 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
452 const volVectorField& U = this->U_;
453 fv::options& fvOptions(fv::options::New(this->mesh_));
454
455 BasicEddyViscosityModel::correct();
456
457 const volScalarField chi(this->chi());
458 const volScalarField fv1(this->fv1(chi));
459 const volScalarField ft2(this->ft2(chi));
460
461 tmp<volTensorField> tgradU = fvc::grad(U);
462 volScalarField dTilda(this->dTilda(chi, fv1, tgradU()));
463 volScalarField Stilda(this->Stilda(chi, fv1, tgradU(), dTilda));
464 tgradU.clear();
465
466 tmp<fvScalarMatrix> nuTildaEqn
467 (
468 fvm::ddt(alpha, rho, nuTilda_)
469 + fvm::div(alphaRhoPhi, nuTilda_)
470 - fvm::laplacian(alpha*rho*DnuTildaEff(), nuTilda_)
471 - Cb2_/sigmaNut_*alpha()*rho()*magSqr(fvc::grad(nuTilda_)()())
472 ==
473 Cb1_*alpha()*rho()*Stilda()*nuTilda_()*(scalar(1) - ft2())
474 - fvm::Sp
475 (
476 (Cw1_*fw(Stilda, dTilda) - Cb1_/sqr(kappa_)*ft2())
477 *alpha()*rho()*nuTilda_()/sqr(dTilda()),
478 nuTilda_
479 )
480 + fvOptions(alpha, rho, nuTilda_)
481 );
482
483 nuTildaEqn.ref().relax();
484 fvOptions.constrain(nuTildaEqn.ref());
485 solve(nuTildaEqn);
486 fvOptions.correct(nuTilda_);
487 bound(nuTilda_, dimensionedScalar(nuTilda_.dimensions(), Zero));
488 nuTilda_.correctBoundaryConditions();
489 }
490
491 correctNut();
492}
493
494
495// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
496
497} // End namespace Foam
498
499// ************************************************************************* //
scalar y
label k
Bound the given scalar field if it has gone unbounded.
const uniformDimensionedVectorField & g
fv::options & fvOptions
const dimensionSet & dimensions() const noexcept
Return dimensions.
DimensionedField< scalar, volMesh > Internal
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
tmp< volScalarField > chi() const
void correctNut(const volScalarField &fv1)
virtual tmp< volScalarField > k() const
Return the (estimated) turbulent kinetic energy.
virtual void correct()
Correct nuTilda and related properties.
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
eddyViscosity< RASModel< BasicTurbulenceModel > >::transportModel transportModel
tmp< volScalarField > Omega(const volTensorField &gradU) const
volScalarField nuTilda_
Modified kinematic viscosity [m^2/s].
tmp< volScalarField > fv1(const volScalarField &chi) const
virtual tmp< volScalarField > epsilon() const
Return the (estimated) turbulent kinetic energy dissipation rate.
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &Stilda, const volScalarField &dTilda) const
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const =0
Length scale.
tmp< volScalarField::Internal > fw(const volScalarField &Stilda, const volScalarField &dTilda) const
tmp< volScalarField > ft2(const volScalarField &chi) const
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
eddyViscosity< RASModel< BasicTurbulenceModel > >::rhoField rhoField
BasicEddyViscosityModel::alphaField alphaField
virtual tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const
virtual tmp< volScalarField > omega() const
Return the (estimated) specific dissipation rate.
virtual bool read()
Re-read model coefficients if they have changed.
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.
bool readIfPresent(const dictionary &dict)
Update the value of dimensioned<Type> if found in the dictionary, lookup in dictionary with the name(...
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
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
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< 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 Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensionedScalar pow6(const dimensionedScalar &ds)
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.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
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
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedTensor skew(const dimensionedTensor &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & nu
volScalarField & alpha
CEqn solve()