SpalartAllmaras.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) 2019-2021 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 "SpalartAllmaras.H"
30#include "fvOptions.H"
31#include "bound.H"
32#include "wallDist.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace RASModels
39{
40
41// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42
43template<class BasicTurbulenceModel>
45{
46 return nuTilda_/this->nu();
47}
48
49
50template<class BasicTurbulenceModel>
52(
53 const volScalarField& chi
54) const
55{
56 const volScalarField chi3(pow3(chi));
57
58 return chi3/(chi3 + pow3(Cv1_));
59}
60
61
62template<class BasicTurbulenceModel>
64(
65 const volScalarField::Internal& chi,
67) const
68{
69 return scalar(1) - chi/(scalar(1) + chi*fv1);
70}
71
72
73template<class BasicTurbulenceModel>
75const
76{
77 const volScalarField chi(this->chi());
78
79 const volScalarField fv1(this->fv1(chi));
80
81 const volScalarField::Internal Omega
82 (
83 ::sqrt(scalar(2))*mag(skew(fvc::grad(this->U_)().v()))
84 );
85
86 return
87 (
88 max
89 (
90 Omega + fv2(chi(), fv1())*nuTilda_()/sqr(kappa_*y_),
91 Cs_*Omega
92 )
93 );
94}
95
96
97template<class BasicTurbulenceModel>
99(
100 const volScalarField::Internal& Stilda
101) const
102{
104 (
105 min
106 (
107 nuTilda_
108 /(
109 max
110 (
111 Stilda,
112 dimensionedScalar(Stilda.dimensions(), SMALL)
113 )
114 *sqr(kappa_*y_)
115 ),
116 scalar(10)
117 )
118 );
119
120 const volScalarField::Internal g(r + Cw2_*(pow6(r) - r));
121
122 return
123 g*pow
124 (
125 (scalar(1) + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)),
126 scalar(1)/scalar(6)
127 );
128}
129
130
131template<class BasicTurbulenceModel>
133{
134 this->nut_ = nuTilda_*this->fv1(this->chi());
135 this->nut_.correctBoundaryConditions();
136 fv::options::New(this->mesh_).correct(this->nut_);
137
138 BasicTurbulenceModel::correctNut();
139}
140
141
142// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
143
144template<class BasicTurbulenceModel>
146(
147 const alphaField& alpha,
148 const rhoField& rho,
149 const volVectorField& U,
150 const surfaceScalarField& alphaRhoPhi,
151 const surfaceScalarField& phi,
152 const transportModel& transport,
153 const word& propertiesName,
154 const word& type
155)
156:
157 eddyViscosity<RASModel<BasicTurbulenceModel>>
158 (
159 type,
160 alpha,
161 rho,
162 U,
163 alphaRhoPhi,
164 phi,
165 transport,
166 propertiesName
167 ),
168
169 sigmaNut_
170 (
171 dimensioned<scalar>::getOrAddToDict
172 (
173 "sigmaNut",
174 this->coeffDict_,
175 scalar(2)/scalar(3)
176 )
177 ),
178 kappa_
179 (
180 dimensioned<scalar>::getOrAddToDict
181 (
182 "kappa",
183 this->coeffDict_,
184 0.41
185 )
186 ),
187
188 Cb1_
189 (
190 dimensioned<scalar>::getOrAddToDict
191 (
192 "Cb1",
193 this->coeffDict_,
194 0.1355
195 )
196 ),
197 Cb2_
198 (
199 dimensioned<scalar>::getOrAddToDict
200 (
201 "Cb2",
202 this->coeffDict_,
203 0.622
204 )
205 ),
206 Cw1_(Cb1_/sqr(kappa_) + (scalar(1) + Cb2_)/sigmaNut_),
207 Cw2_
208 (
209 dimensioned<scalar>::getOrAddToDict
210 (
211 "Cw2",
212 this->coeffDict_,
213 0.3
214 )
215 ),
216 Cw3_
217 (
218 dimensioned<scalar>::getOrAddToDict
219 (
220 "Cw3",
221 this->coeffDict_,
222 2.0
223 )
224 ),
225 Cv1_
226 (
227 dimensioned<scalar>::getOrAddToDict
228 (
229 "Cv1",
230 this->coeffDict_,
231 7.1
232 )
233 ),
234 Cs_
235 (
236 dimensioned<scalar>::getOrAddToDict
237 (
238 "Cs",
239 this->coeffDict_,
240 0.3
241 )
242 ),
243
244 nuTilda_
245 (
247 (
248 "nuTilda",
249 this->runTime_.timeName(),
250 this->mesh_,
251 IOobject::MUST_READ,
252 IOobject::AUTO_WRITE
253 ),
254 this->mesh_
255 ),
256
257 y_(wallDist::New(this->mesh_).y())
258{
259 if (type == typeName)
260 {
261 this->printCoeffs(type);
262 }
263}
264
265
266// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
267
268template<class BasicTurbulenceModel>
270{
272 {
273 sigmaNut_.readIfPresent(this->coeffDict());
274 kappa_.readIfPresent(this->coeffDict());
275
276 Cb1_.readIfPresent(this->coeffDict());
277 Cb2_.readIfPresent(this->coeffDict());
278 Cw1_ = Cb1_/sqr(kappa_) + (scalar(1) + Cb2_)/sigmaNut_;
279 Cw2_.readIfPresent(this->coeffDict());
280 Cw3_.readIfPresent(this->coeffDict());
281 Cv1_.readIfPresent(this->coeffDict());
282 Cs_.readIfPresent(this->coeffDict());
283
284 return true;
285 }
286
287 return false;
288}
289
290
291template<class BasicTurbulenceModel>
293{
295 (
296 "DnuTildaEff",
297 (nuTilda_ + this->nu())/sigmaNut_
298 );
299}
300
301
302template<class BasicTurbulenceModel>
304{
305 // (B:Eq. 4.50)
306 const scalar Cmu = 0.09;
307
309 (
311 (
312 IOobject::groupName("k", this->alphaRhoPhi_.group()),
313 this->runTime_.timeName(),
314 this->mesh_
315 ),
316 cbrt(this->fv1(this->chi()))
317 *nuTilda_
318 *::sqrt(scalar(2)/Cmu)
319 *mag(symm(fvc::grad(this->U_))),
320 this->nut_.boundaryField().types()
321 );
322}
323
324
325template<class BasicTurbulenceModel>
327{
328 // (B:Eq. 4.50)
329 const scalar Cmu = 0.09;
330 const dimensionedScalar nutSMALL(sqr(dimLength)/dimTime, SMALL);
331
333 (
335 (
336 IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
337 this->runTime_.timeName(),
338 this->mesh_
339 ),
340 pow(this->fv1(this->chi()), 0.5)
341 *pow(::sqrt(Cmu)*this->k(), 2)
342 /(nuTilda_ + this->nut_ + nutSMALL),
343 this->nut_.boundaryField().types()
344 );
345}
346
347
348template<class BasicTurbulenceModel>
350{
351 // (P:p. 384)
352 const scalar betaStar = 0.09;
353 const dimensionedScalar k0(sqr(dimLength/dimTime), SMALL);
354
356 (
358 (
359 IOobject::groupName("omega", this->alphaRhoPhi_.group()),
360 this->runTime_.timeName(),
361 this->mesh_
362 ),
363 this->epsilon()/(betaStar*(this->k() + k0)),
364 this->nut_.boundaryField().types()
365 );
366}
367
368
369template<class BasicTurbulenceModel>
371{
372 if (!this->turbulence_)
373 {
374 return;
375 }
376
377 {
378 // Construct local convenience references
379 const alphaField& alpha = this->alpha_;
380 const rhoField& rho = this->rho_;
381 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
383
385
386 const volScalarField::Internal Stilda(this->Stilda());
387
388 tmp<fvScalarMatrix> nuTildaEqn
389 (
390 fvm::ddt(alpha, rho, nuTilda_)
391 + fvm::div(alphaRhoPhi, nuTilda_)
392 - fvm::laplacian(alpha*rho*DnuTildaEff(), nuTilda_)
393 - Cb2_/sigmaNut_*alpha*rho*magSqr(fvc::grad(nuTilda_))
394 ==
395 Cb1_*alpha()*rho()*Stilda*nuTilda_()
396 - fvm::Sp(Cw1_*alpha()*rho()*fw(Stilda)*nuTilda_()/sqr(y_), nuTilda_)
397 + fvOptions(alpha, rho, nuTilda_)
398 );
399
400 nuTildaEqn.ref().relax();
401 fvOptions.constrain(nuTildaEqn.ref());
402 solve(nuTildaEqn);
403 fvOptions.correct(nuTilda_);
404 bound(nuTilda_, dimensionedScalar(nuTilda_.dimensions(), Zero));
405 nuTilda_.correctBoundaryConditions();
406 }
407
408 // Update nut with latest available nuTilda
409 correctNut();
410}
411
412
413// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
414
415} // End namespace RASModels
416} // End namespace Foam
417
418// ************************************************************************* //
scalar y
label k
Bound the given scalar field if it has gone unbounded.
const uniformDimensionedVectorField & g
fv::options & fvOptions
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:55
BasicTurbulenceModel::alphaField alphaField
tmp< volScalarField > chi() const
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > k() const
Return the (estimated) turbulent kinetic energy.
virtual void correct()
Solve the turbulence equations and correct the turbulent viscosity.
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
tmp< volScalarField > fv1(const volScalarField &chi) const
tmp< volScalarField::Internal > fv2(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
virtual tmp< volScalarField > epsilon() const
Return the (estimated) turbulent kinetic energy dissipation rate.
virtual void correctNut()
Update nut with the latest available nuTilda.
BasicTurbulenceModel::transportModel transportModel
tmp< volScalarField::Internal > fw(const volScalarField::Internal &Stilda) const
virtual tmp< volScalarField > omega() const
Return the (estimated) specific dissipation rate.
tmp< volScalarField::Internal > Stilda() const
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Generic dimensioned Type class.
bool readIfPresent(const dictionary &dict)
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:58
Finite-volume options.
Definition: fvOptions.H:59
A class for managing temporary objects.
Definition: tmp.H:65
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
thermo correct()
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< 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 Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensionedScalar pow6(const dimensionedScalar &ds)
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)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
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
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 pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedTensor skew(const dimensionedTensor &dt)
volScalarField & nu
volScalarField & alpha
CEqn solve()