Loading...
Searching...
No Matches
kOmegaSSTBase.H
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) 2022 Upstream CFD GmbH
10 Copyright (C) 2017-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
28Class
29 Foam::kOmegaSSTBase
30
31Description
32 Base class implementation of the k-omega-SST turbulence model for
33 incompressible and compressible flows.
34
35 Turbulence model described in:
36 \verbatim
37 Menter, F. R. & Esch, T. (2001).
38 Elements of Industrial Heat Transfer Prediction.
39 16th Brazilian Congress of Mechanical Engineering (COBEM).
40 \endverbatim
41
42 with updated coefficients from
43 \verbatim
44 Menter, F. R., Kuntz, M., and Langtry, R. (2003).
45 Ten Years of Industrial Experience with the SST Turbulence Model.
46 Turbulence, Heat and Mass Transfer 4, ed: K. Hanjalic, Y. Nagano,
47 & M. Tummers, Begell House, Inc., 625 - 632.
48 \endverbatim
49
50 but with the consistent production terms from the 2001 paper as form in the
51 2003 paper is a typo, see
52 \verbatim
53 http://turbmodels.larc.nasa.gov/sst.html
54 \endverbatim
55
56 and the addition of the optional F3 term for rough walls from
57 \verbatim
58 Hellsten, A. (1998).
59 Some Improvements in Menter's k-omega-SST turbulence model
60 29th AIAA Fluid Dynamics Conference, AIAA-98-2554.
61 \endverbatim
62
63 and the optional decay control from:
64 \verbatim
65 Spalart, P. R. and Rumsey, C. L. (2007).
66 Effective Inflow Conditions for Turbulence Models in Aerodynamic
67 Calculations
68 AIAA Journal, 45(10), 2544 - 2553.
69 \endverbatim
70
71 Note that this implementation is written in terms of alpha diffusion
72 coefficients rather than the more traditional sigma (alpha = 1/sigma) so
73 that the blending can be applied to all coefficients in a consistent
74 manner. The paper suggests that sigma is blended but this would not be
75 consistent with the blending of the k-epsilon and k-omega models.
76
77 Also note that the error in the last term of equation (2) relating to
78 sigma has been corrected.
79
80 Wall-functions are applied in this implementation by using equations (14)
81 to specify the near-wall omega as appropriate.
82
83 The blending functions (15) and (16) are not currently used because of the
84 uncertainty in their origin, range of applicability and that if y+ becomes
85 sufficiently small blending u_tau in this manner clearly becomes nonsense.
86
87 The default model coefficients are
88 \verbatim
89 kOmegaSSTBaseCoeffs
90 {
91 alphaK1 0.85;
92 alphaK2 1.0;
93 alphaOmega1 0.5;
94 alphaOmega2 0.856;
95 beta1 0.075;
96 beta2 0.0828;
97 betaStar 0.09;
98 gamma1 5/9;
99 gamma2 0.44;
100 a1 0.31;
101 b1 1.0;
102 c1 10.0;
103 F3 no;
104
105 // Optional decay control
106 decayControl yes;
107 kInf <far-field k value>;
108 omegaInf <far-field omega value>;
109 }
110 \endverbatim
111
112SourceFiles
113 kOmegaSSTBase.C
114
115\*---------------------------------------------------------------------------*/
116
117#ifndef kOmegaSSTBase_H
118#define kOmegaSSTBase_H
119
120// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121
122namespace Foam
123{
124
125/*---------------------------------------------------------------------------*\
126 Class kOmegaSSTBase Declaration
127\*---------------------------------------------------------------------------*/
128
129template<class BasicEddyViscosityModel>
130class kOmegaSSTBase
131:
132 public BasicEddyViscosityModel
133{
134protected:
136 // Protected Data
137
138 // Model coefficients
148
153
157
158 //- Flag to include the F3 term
159 Switch F3_;
160
161
162 // Fields
163
164 //- Wall distance
165 // Note: different to wall distance in parent RASModel
166 // which is for near-wall cells only
168
169 //- Turbulent kinetic energy field [m^2/s^2]
171
172 //- Specific dissipation rate field [1/s]
174
175
176 // Decay control
178 //- Flag to include the decay control
182
183
184 // Protected Member Functions
186 //- Set decay control with kInf and omegaInf
188
189 virtual tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
190 virtual tmp<volScalarField> F2() const;
191 virtual tmp<volScalarField> F3() const;
192 virtual tmp<volScalarField> F23() const;
193
194 //- Return the blended field
196 (
197 const volScalarField& F1,
198 const dimensionedScalar& psi1,
200 ) const
201 {
202 return F1*(psi1 - psi2) + psi2;
203 }
204
205 //- Return the internal blended field
207 (
209 const dimensionedScalar& psi1,
211 ) const
212 {
213 return F1*(psi1 - psi2) + psi2;
214 }
215
217 {
219 }
220
222 {
224 }
225
226 tmp<volScalarField::Internal> beta
227 (
232 (
233 IOobject::scopedName(this->type(), "beta"),
234 blend(F1, beta1_, beta2_)
235 );
236 }
237
239 (
241 ) const
242 {
244 (
245 IOobject::scopedName(this->type(), "gamma"),
246 blend(F1, gamma1_, gamma2_)
247 );
248 }
249
250 virtual void correctNut(const volScalarField& S2);
251
252 virtual void correctNut();
253
254 //- Return square of strain rate
255 virtual tmp<volScalarField> S2
256 (
257 const volTensorField& gradU
258 ) const;
259
260 //- Return k production rate
262 (
264 ) const;
265
266 //- Return epsilon/k which for standard RAS is betaStar*omega
268 (
269 const volScalarField& F1,
270 const volTensorField& gradU
271 ) const;
272
273 //- Return (G/nu)_0
275 (
276 const volTensorField& gradU,
277 const volScalarField& S2
278 ) const;
279
280 //- Return G/nu
282 (
286 ) const;
287
288 virtual tmp<fvScalarMatrix> kSource() const;
289
290 virtual tmp<fvScalarMatrix> omegaSource() const;
291
293 (
297 ) const;
298
299
300public:
301
302 typedef typename BasicEddyViscosityModel::alphaField alphaField;
303 typedef typename BasicEddyViscosityModel::rhoField rhoField;
304 typedef typename BasicEddyViscosityModel::transportModel transportModel;
305
306
307 // Generated Methods
308
309 //- No copy construct
310 kOmegaSSTBase(const kOmegaSSTBase&) = delete;
311
312 //- No copy assignment
313 void operator=(const kOmegaSSTBase&) = delete;
314
315
316 // Constructors
317
318 //- Construct from components
320 (
321 const word& type,
322 const alphaField& alpha,
323 const rhoField& rho,
325 const surfaceScalarField& alphaRhoPhi,
327 const transportModel& transport,
328 const word& propertiesName = turbulenceModel::propertiesName
329 );
330
331
332 //- Destructor
333 virtual ~kOmegaSSTBase() = default;
335
336 // Member Functions
337
338 //- Re-read model coefficients if they have changed
339 virtual bool read();
340
341 //- Return the effective diffusivity for k
343 {
345 (
346 "DkEff",
347 alphaK(F1)*this->nut_ + this->nu()
348 );
349 }
350
351 //- Return the effective diffusivity for omega
353 {
355 (
356 "DomegaEff",
357 alphaOmega(F1)*this->nut_ + this->nu()
358 );
359 }
360
361 //- Return the turbulence kinetic energy
362 virtual tmp<volScalarField> k() const
364 return k_;
365 }
366
367 //- Return the turbulence kinetic energy dissipation rate
368 virtual tmp<volScalarField> omega() const
369 {
370 return omega_;
371 }
372
373 //- Solve the turbulence equations and correct the turbulence viscosity
374 virtual void correct();
375};
377
378// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
379
380} // End namespace Foam
381
382// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
383
384#ifdef NoRepository
385 #include "kOmegaSSTBase.C"
386#endif
387
388// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
389
390#endif
391
392// ************************************************************************* //
#define F2(B, C, D)
Definition SHA1.C:150
#define F1(B, C, D)
Definition SHA1.C:149
const volScalarField & psi2
const volScalarField & psi1
DimensionedField< scalar, volMesh > Internal
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
tmp< volScalarField::Internal > blend(const volScalarField::Internal &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
Return the internal blended field.
void setDecayControl(const dictionary &dict)
Set decay control with kInf and omegaInf.
dimensionedScalar beta2_
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.
dimensionedScalar alphaOmega2_
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< volScalarField > F2() const
dimensionedScalar gamma2_
virtual tmp< volScalarField > F23() const
BasicEddyViscosityModel::transportModel transportModel
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
Return the blended field.
dimensionedScalar b1_
dimensionedScalar beta1_
dimensionedScalar omegaInf_
virtual tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Return G/nu.
tmp< volScalarField > alphaK(const volScalarField &F1) const
dimensionedScalar c1_
dimensionedScalar alphaK2_
const volScalarField & y_
Wall distance.
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Return the effective diffusivity for omega.
dimensionedScalar a1_
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
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
Switch decayControl_
Flag to include the decay control.
virtual tmp< fvScalarMatrix > omegaSource() const
virtual void correctNut()
virtual tmp< fvScalarMatrix > kSource() const
virtual ~kOmegaSSTBase()=default
Destructor.
dimensionedScalar alphaOmega1_
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 > omega() const
Return the turbulence kinetic energy dissipation rate.
void operator=(const kOmegaSSTBase &)=delete
No copy assignment.
dimensionedScalar gamma1_
virtual tmp< volScalarField > S2(const volTensorField &gradU) const
Return square of strain rate.
BasicEddyViscosityModel::rhoField rhoField
BasicEddyViscosityModel::alphaField alphaField
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.
dimensionedScalar betaStar_
dimensionedScalar kInf_
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
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
static const word propertiesName
Default name of the turbulence properties dictionary.
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
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & nu
volScalarField & alpha
dictionary dict
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)