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) 2017-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 Class
28  Foam::kOmegaSSTBase
29 
30 Description
31  Base class implementation of the k-omega-SST turbulence model for
32  incompressible and compressible flows.
33 
34  Turbulence model described in:
35  \verbatim
36  Menter, F. R. & Esch, T. (2001).
37  Elements of Industrial Heat Transfer Prediction.
38  16th Brazilian Congress of Mechanical Engineering (COBEM).
39  \endverbatim
40 
41  with updated coefficients from
42  \verbatim
43  Menter, F. R., Kuntz, M., and Langtry, R. (2003).
44  Ten Years of Industrial Experience with the SST Turbulence Model.
45  Turbulence, Heat and Mass Transfer 4, ed: K. Hanjalic, Y. Nagano,
46  & M. Tummers, Begell House, Inc., 625 - 632.
47  \endverbatim
48 
49  but with the consistent production terms from the 2001 paper as form in the
50  2003 paper is a typo, see
51  \verbatim
52  http://turbmodels.larc.nasa.gov/sst.html
53  \endverbatim
54 
55  and the addition of the optional F3 term for rough walls from
56  \verbatim
57  Hellsten, A. (1998).
58  Some Improvements in Menter's k-omega-SST turbulence model
59  29th AIAA Fluid Dynamics Conference, AIAA-98-2554.
60  \endverbatim
61 
62  and the optional decay control from:
63  \verbatim
64  Spalart, P. R. and Rumsey, C. L. (2007).
65  Effective Inflow Conditions for Turbulence Models in Aerodynamic
66  Calculations
67  AIAA Journal, 45(10), 2544 - 2553.
68  \endverbatim
69 
70  Note that this implementation is written in terms of alpha diffusion
71  coefficients rather than the more traditional sigma (alpha = 1/sigma) so
72  that the blending can be applied to all coefficients in a consistent
73  manner. The paper suggests that sigma is blended but this would not be
74  consistent with the blending of the k-epsilon and k-omega models.
75 
76  Also note that the error in the last term of equation (2) relating to
77  sigma has been corrected.
78 
79  Wall-functions are applied in this implementation by using equations (14)
80  to specify the near-wall omega as appropriate.
81 
82  The blending functions (15) and (16) are not currently used because of the
83  uncertainty in their origin, range of applicability and that if y+ becomes
84  sufficiently small blending u_tau in this manner clearly becomes nonsense.
85 
86  The default model coefficients are
87  \verbatim
88  kOmegaSSTBaseCoeffs
89  {
90  alphaK1 0.85;
91  alphaK2 1.0;
92  alphaOmega1 0.5;
93  alphaOmega2 0.856;
94  beta1 0.075;
95  beta2 0.0828;
96  betaStar 0.09;
97  gamma1 5/9;
98  gamma2 0.44;
99  a1 0.31;
100  b1 1.0;
101  c1 10.0;
102  F3 no;
103 
104  // Optional decay control
105  decayControl yes;
106  kInf <far-field k value>;
107  omegaInf <far-field omega value>;
108  }
109  \endverbatim
110 
111 SourceFiles
112  kOmegaSSTBase.C
113 
114 \*---------------------------------------------------------------------------*/
115 
116 #ifndef kOmegaSSTBase_H
117 #define kOmegaSSTBase_H
118 
119 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120 
121 namespace Foam
122 {
123 
124 /*---------------------------------------------------------------------------*\
125  Class kOmegaSSTBase Declaration
126 \*---------------------------------------------------------------------------*/
127 
128 template<class BasicEddyViscosityModel>
129 class kOmegaSSTBase
130 :
131  public BasicEddyViscosityModel
132 {
133  // Private Member Functions
134 
135  //- No copy construct
136  kOmegaSSTBase(const kOmegaSSTBase&) = delete;
137 
138  //- No copy assignment
139  void operator=(const kOmegaSSTBase&) = delete;
140 
141 
142 protected:
143 
144  // Protected data
145 
146  // Model coefficients
147 
150 
153 
156 
159 
161 
165 
166  //- Flag to include the F3 term
167  Switch F3_;
168 
169 
170  // Fields
171 
172  //- Wall distance
173  // Note: different to wall distance in parent RASModel
174  // which is for near-wall cells only
175  const volScalarField& y_;
176 
179 
180 
181  // Decay control
182 
183  //- Flag to include the decay control
187 
188 
189  // Protected Member Functions
190 
191  void setDecayControl(const dictionary& dict);
192 
193  virtual tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
194  virtual tmp<volScalarField> F2() const;
195  virtual tmp<volScalarField> F3() const;
196  virtual tmp<volScalarField> F23() const;
197 
199  (
200  const volScalarField& F1,
201  const dimensionedScalar& psi1,
202  const dimensionedScalar& psi2
203  ) const
204  {
205  return F1*(psi1 - psi2) + psi2;
206  }
207 
209  (
211  const dimensionedScalar& psi1,
212  const dimensionedScalar& psi2
213  ) const
214  {
215  return F1*(psi1 - psi2) + psi2;
216  }
217 
219  {
220  return blend(F1, alphaK1_, alphaK2_);
221  }
222 
224  {
225  return blend(F1, alphaOmega1_, alphaOmega2_);
226  }
227 
229  (
231  ) const
232  {
234  (
235  this->type() + ":beta",
236  blend(F1, beta1_, beta2_)
237  );
238  }
239 
241  (
243  ) const
244  {
246  (
247  this->type() + ":gamma",
248  blend(F1, gamma1_, gamma2_)
249  );
250  }
251 
252  virtual void correctNut(const volScalarField& S2);
253 
254  virtual void correctNut();
255 
256  //- Return k production rate
258  (
260  ) const;
261 
262  //- Return epsilon/k which for standard RAS is betaStar*omega
264  (
265  const volScalarField& F1,
266  const volTensorField& gradU
267  ) const;
268 
269  //- Return G/nu
271  (
272  const volScalarField::Internal& GbyNu0,
274  const volScalarField::Internal& S2
275  ) const;
276 
277  virtual tmp<fvScalarMatrix> kSource() const;
278 
279  virtual tmp<fvScalarMatrix> omegaSource() const;
280 
281  virtual tmp<fvScalarMatrix> Qsas
282  (
283  const volScalarField::Internal& S2,
286  ) const;
287 
288 
289 public:
290 
291  typedef typename BasicEddyViscosityModel::alphaField alphaField;
292  typedef typename BasicEddyViscosityModel::rhoField rhoField;
293  typedef typename BasicEddyViscosityModel::transportModel transportModel;
294 
295 
296  // Constructors
297 
298  //- Construct from components
300  (
301  const word& type,
302  const alphaField& alpha,
303  const rhoField& rho,
304  const volVectorField& U,
305  const surfaceScalarField& alphaRhoPhi,
306  const surfaceScalarField& phi,
307  const transportModel& transport,
308  const word& propertiesName = turbulenceModel::propertiesName
309  );
310 
311 
312  //- Destructor
313  virtual ~kOmegaSSTBase() = default;
314 
315 
316  // Member Functions
317 
318  //- Re-read model coefficients if they have changed
319  virtual bool read();
320 
321  //- Return the effective diffusivity for k
323  {
324  return tmp<volScalarField>
325  (
326  new volScalarField("DkEff", alphaK(F1)*this->nut_ + this->nu())
327  );
328  }
329 
330  //- Return the effective diffusivity for omega
332  {
333  return tmp<volScalarField>
334  (
335  new volScalarField
336  (
337  "DomegaEff",
338  alphaOmega(F1)*this->nut_ + this->nu()
339  )
340  );
341  }
342 
343  //- Return the turbulence kinetic energy
344  virtual tmp<volScalarField> k() const
345  {
346  return k_;
347  }
348 
349  //- Return the turbulence kinetic energy dissipation rate
350  virtual tmp<volScalarField> omega() const
351  {
352  return omega_;
353  }
354 
355  //- Solve the turbulence equations and correct the turbulence viscosity
356  virtual void correct();
357 };
358 
359 
360 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361 
362 } // End namespace Foam
363 
364 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
365 
366 #ifdef NoRepository
367  #include "kOmegaSSTBase.C"
368 #endif
369 
370 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
371 
372 #endif
373 
374 // ************************************************************************* //
Foam::kOmegaSSTBase::Qsas
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
Definition: kOmegaSSTBase.C:204
Foam::kOmegaSSTBase::F2
virtual tmp< volScalarField > F2() const
Definition: kOmegaSSTBase.C:72
Foam::kOmegaSSTBase::omega
virtual tmp< volScalarField > omega() const
Return the turbulence kinetic energy dissipation rate.
Definition: kOmegaSSTBase.H:349
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::kOmegaSSTBase::rhoField
BasicEddyViscosityModel::rhoField rhoField
Definition: kOmegaSSTBase.H:291
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
kOmegaSSTBase.C
Foam::kOmegaSSTBase::decayControl_
Switch decayControl_
Flag to include the decay control.
Definition: kOmegaSSTBase.H:183
Foam::kOmegaSSTBase::alphaK2_
dimensionedScalar alphaK2_
Definition: kOmegaSSTBase.H:148
Foam::kOmegaSSTBase::kInf_
dimensionedScalar kInf_
Definition: kOmegaSSTBase.H:184
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::kOmegaSSTBase::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: kOmegaSSTBase.C:175
Foam::kOmegaSSTBase::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: kOmegaSSTBase.H:343
Foam::kOmegaSSTBase::F1
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Definition: kOmegaSSTBase.C:43
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::kOmegaSSTBase::omegaSource
virtual tmp< fvScalarMatrix > omegaSource() const
Definition: kOmegaSSTBase.C:189
Foam::kOmegaSSTBase::epsilonByk
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Return epsilon/k which for standard RAS is betaStar*omega.
Definition: kOmegaSSTBase.C:148
Foam::kOmegaSSTBase::omega_
volScalarField omega_
Definition: kOmegaSSTBase.H:177
Foam::kOmegaSSTBase::setDecayControl
void setDecayControl(const dictionary &dict)
Definition: kOmegaSSTBase.C:433
Foam::kOmegaSSTBase::alphaOmega1_
dimensionedScalar alphaOmega1_
Definition: kOmegaSSTBase.H:150
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::kOmegaSSTBase::beta1_
dimensionedScalar beta1_
Definition: kOmegaSSTBase.H:156
Foam::kOmegaSSTBase::F3_
Switch F3_
Flag to include the F3 term.
Definition: kOmegaSSTBase.H:166
Foam::kOmegaSSTBase::Pk
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
Definition: kOmegaSSTBase.C:137
Foam::kOmegaSSTBase::alphaOmega
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
Definition: kOmegaSSTBase.H:222
Foam::kOmegaSSTBase::betaStar_
dimensionedScalar betaStar_
Definition: kOmegaSSTBase.H:159
Foam::kOmegaSSTBase::gamma
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
Definition: kOmegaSSTBase.H:240
Foam::kOmegaSSTBase::alphaOmega2_
dimensionedScalar alphaOmega2_
Definition: kOmegaSSTBase.H:151
rho
rho
Definition: readInitialConditions.H:88
psi2
const volScalarField & psi2
Definition: setRegionFluidFields.H:31
Foam::kOmegaSSTBase::DomegaEff
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Return the effective diffusivity for omega.
Definition: kOmegaSSTBase.H:330
Foam::kOmegaSSTBase::blend
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
Definition: kOmegaSSTBase.H:198
Foam::kOmegaSSTBase::correctNut
virtual void correctNut()
Definition: kOmegaSSTBase.C:129
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::kOmegaSSTBase::GbyNu
virtual tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Return G/nu.
Definition: kOmegaSSTBase.C:159
Foam::kOmegaSSTBase::gamma2_
dimensionedScalar gamma2_
Definition: kOmegaSSTBase.H:154
Foam::kOmegaSSTBase::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTBase.C:456
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::kOmegaSSTBase::~kOmegaSSTBase
virtual ~kOmegaSSTBase()=default
Destructor.
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:98
Foam::kOmegaSSTBase::c1_
dimensionedScalar c1_
Definition: kOmegaSSTBase.H:163
Foam::kOmegaSSTBase::beta2_
dimensionedScalar beta2_
Definition: kOmegaSSTBase.H:157
Foam::dimensioned< scalar >
Foam::kOmegaSSTBase::b1_
dimensionedScalar b1_
Definition: kOmegaSSTBase.H:162
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::kOmegaSSTBase::y_
const volScalarField & y_
Wall distance.
Definition: kOmegaSSTBase.H:174
psi1
const volScalarField & psi1
Definition: setRegionFluidFields.H:28
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:97
Foam::kOmegaSSTBase::DkEff
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.
Definition: kOmegaSSTBase.H:321
Foam::kOmegaSSTBase::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTBase.C:484
Foam::kOmegaSSTBase::F23
virtual tmp< volScalarField > F23() const
Definition: kOmegaSSTBase.C:102
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:99
U
U
Definition: pEqn.H:72
Foam::kOmegaSSTBase::alphaK
tmp< volScalarField > alphaK(const volScalarField &F1) const
Definition: kOmegaSSTBase.H:217
Foam::kOmegaSSTBase::k_
volScalarField k_
Definition: kOmegaSSTBase.H:176
Foam::kOmegaSSTBase::alphaK1_
dimensionedScalar alphaK1_
Definition: kOmegaSSTBase.H:147
Foam::kOmegaSSTBase::transportModel
BasicEddyViscosityModel::transportModel transportModel
Definition: kOmegaSSTBase.H:292
Foam::kOmegaSSTBase::beta
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
Definition: kOmegaSSTBase.H:228
Foam::kOmegaSSTBase::gamma1_
dimensionedScalar gamma1_
Definition: kOmegaSSTBase.H:153
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::kOmegaSSTBase::omegaInf_
dimensionedScalar omegaInf_
Definition: kOmegaSSTBase.H:185
Foam::kOmegaSSTBase::alphaField
BasicEddyViscosityModel::alphaField alphaField
Definition: kOmegaSSTBase.H:290
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::kOmegaSSTBase::F3
virtual tmp< volScalarField > F3() const
Definition: kOmegaSSTBase.C:89
Foam::kOmegaSSTBase
Base class implementation of the k-omega-SST turbulence model for incompressible and compressible flo...
Definition: kOmegaSSTBase.H:128
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::kOmegaSSTBase::a1_
dimensionedScalar a1_
Definition: kOmegaSSTBase.H:161