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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::kOmegaSST
26 
27 Group
28  grpTurbulence
29 
30 Description
31  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  Note that this implementation is written in terms of alpha diffusion
63  coefficients rather than the more traditional sigma (alpha = 1/sigma) so
64  that the blending can be applied to all coefficuients in a consistent
65  manner. The paper suggests that sigma is blended but this would not be
66  consistent with the blending of the k-epsilon and k-omega models.
67 
68  Also note that the error in the last term of equation (2) relating to
69  sigma has been corrected.
70 
71  Wall-functions are applied in this implementation by using equations (14)
72  to specify the near-wall omega as appropriate.
73 
74  The blending functions (15) and (16) are not currently used because of the
75  uncertainty in their origin, range of applicability and that if y+ becomes
76  sufficiently small blending u_tau in this manner clearly becomes nonsense.
77 
78  The default model coefficients are
79  \verbatim
80  kOmegaSSTCoeffs
81  {
82  alphaK1 0.85;
83  alphaK2 1.0;
84  alphaOmega1 0.5;
85  alphaOmega2 0.856;
86  beta1 0.075;
87  beta2 0.0828;
88  betaStar 0.09;
89  gamma1 5/9;
90  gamma2 0.44;
91  a1 0.31;
92  b1 1.0;
93  c1 10.0;
94  F3 no;
95  }
96  \endverbatim
97 
98 SourceFiles
99  kOmegaSST.C
100 
101 \*---------------------------------------------------------------------------*/
102 
103 #ifndef kOmegaSSTBase_H
104 #define kOmegaSSTBase_H
105 
106 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
107 
108 namespace Foam
109 {
110 
111 /*---------------------------------------------------------------------------*\
112  Class kOmegaSST Declaration
113 \*---------------------------------------------------------------------------*/
114 
115 template<class TurbulenceModel, class BasicTurbulenceModel>
116 class kOmegaSST
117 :
118  public TurbulenceModel
119 {
120  // Private Member Functions
121 
122  // Disallow default bitwise copy construct and assignment
123  kOmegaSST(const kOmegaSST&);
124  void operator=(const kOmegaSST&);
125 
126 
127 protected:
128 
129  // Protected data
130 
131  // Model coefficients
151  Switch F3_;
152 
153 
154  // Fields
155 
156  //- Wall distance
157  // Note: different to wall distance in parent RASModel
158  // which is for near-wall cells only
159  const volScalarField& y_;
163 
164 
165  // Private Member Functions
166 
167  tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
168  tmp<volScalarField> F2() const;
169  tmp<volScalarField> F3() const;
170  tmp<volScalarField> F23() const;
171 
173  (
174  const volScalarField& F1,
175  const dimensionedScalar& psi1,
176  const dimensionedScalar& psi2
177  ) const
178  {
179  return F1*(psi1 - psi2) + psi2;
180  }
183  {
184  return blend(F1, alphaK1_, alphaK2_);
185  }
188  {
189  return blend(F1, alphaOmega1_, alphaOmega2_);
190  }
193  {
194  return blend(F1, beta1_, beta2_);
195  }
198  {
199  return blend(F1, gamma1_, gamma2_);
200  }
201 
202  void correctNut(const volScalarField& S2, const volScalarField& F2);
203 
204 
205  // Protected Member Functions
206 
207  virtual void correctNut();
208 
209  //- Return epsilon/k which for standard RAS is betaStar*omega
211  (
212  const volScalarField& F1,
213  const volScalarField& F2
214  ) const;
215 
216  virtual tmp<fvScalarMatrix> kSource() const;
217 
218  virtual tmp<fvScalarMatrix> omegaSource() const;
219 
220  virtual tmp<fvScalarMatrix> Qsas
221  (
222  const volScalarField& S2,
223  const volScalarField& gamma,
224  const volScalarField& beta
225  ) const;
226 
227 
228 public:
230  typedef typename BasicTurbulenceModel::alphaField alphaField;
231  typedef typename BasicTurbulenceModel::rhoField rhoField;
232  typedef typename BasicTurbulenceModel::transportModel transportModel;
233 
234 
235  // Constructors
236 
237  //- Construct from components
238  kOmegaSST
239  (
240  const word& type,
241  const alphaField& alpha,
242  const rhoField& rho,
243  const volVectorField& U,
244  const surfaceScalarField& alphaRhoPhi,
245  const surfaceScalarField& phi,
246  const transportModel& transport,
247  const word& propertiesName = turbulenceModel::propertiesName
248  );
249 
250 
251  //- Destructor
252  virtual ~kOmegaSST()
253  {}
254 
255 
256  // Member Functions
257 
258  //- Re-read model coefficients if they have changed
259  virtual bool read();
260 
261  //- Return the effective diffusivity for k
262  tmp<volScalarField> DkEff(const volScalarField& F1) const
263  {
264  return tmp<volScalarField>
265  (
266  new volScalarField("DkEff", alphaK(F1)*this->nut_ + this->nu())
267  );
268  }
269 
270  //- Return the effective diffusivity for omega
272  {
273  return tmp<volScalarField>
274  (
275  new volScalarField
276  (
277  "DomegaEff",
278  alphaOmega(F1)*this->nut_ + this->nu()
279  )
280  );
281  }
282 
283  //- Return the turbulence kinetic energy
284  virtual tmp<volScalarField> k() const
285  {
286  return k_;
287  }
288 
289  //- Return the turbulence kinetic energy dissipation rate
290  virtual tmp<volScalarField> epsilon() const
291  {
292  return tmp<volScalarField>
293  (
294  new volScalarField
295  (
296  IOobject
297  (
298  "epsilon",
299  this->mesh_.time().timeName(),
300  this->mesh_
301  ),
302  betaStar_*k_*omega_,
303  omega_.boundaryField().types()
304  )
305  );
306  }
307 
308  //- Return the turbulence kinetic energy dissipation rate
309  virtual tmp<volScalarField> omega() const
310  {
311  return omega_;
312  }
313 
314  //- Solve the turbulence equations and correct the turbulence viscosity
315  virtual void correct();
316 };
317 
318 
319 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
320 
321 } // End namespace Foam
322 
323 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
324 #ifdef NoRepository
325  #include "kOmegaSSTBase.C"
326 #endif
327 
328 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329 #endif
330 
331 // ************************************************************************* //
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
surfaceScalarField & phi
tmp< volScalarField > beta(const volScalarField &F1) const
virtual ~kOmegaSST()
Destructor.
psi2
Definition: TEqns.H:35
const transportModel & transport() const
Access function to incompressible transport model.
U
Definition: pEqn.H:83
dimensionedScalar a1_
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
dimensionedScalar gamma2_
tmp< volScalarField > nu() const
Return the laminar viscosity.
dimensionedScalar c1_
dimensionedScalar betaStar_
dimensionedScalar beta1_
tmp< volScalarField > alphaK(const volScalarField &F1) const
tmp< volScalarField > gamma(const volScalarField &F1) const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::alphaField alphaField
Templated abstract base class for turbulence models.
BasicTurbulenceModel::rhoField rhoField
virtual tmp< fvScalarMatrix > omegaSource() const
dimensionedScalar gamma1_
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Return the effective diffusivity for omega.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
dimensionedScalar b1_
volScalarField omega_
dimensionedScalar alphaOmega2_
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows...
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
tmp< volScalarField > F1(const volScalarField &CDkOmega) const
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
virtual tmp< volScalarField > epsilonByk(const volScalarField &F1, const volScalarField &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
tmp< volScalarField > F2() const
virtual tmp< volScalarField > omega() const
Return the turbulence kinetic energy dissipation rate.
dimensionedScalar alphaOmega1_
dimensionedScalar alphaK2_
dimensionedScalar alphaK1_
const alphaField & alpha() const
Access function to phase fraction.
psi1
Definition: TEqns.H:34
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
virtual tmp< fvScalarMatrix > Qsas(const volScalarField &S2, const volScalarField &gamma, const volScalarField &beta) const
tmp< volScalarField > F23() const
volScalarField k_
dimensionedScalar beta2_
A class for managing temporary objects.
Definition: PtrList.H:54
const volScalarField & y_
Wall distance.
tmp< volScalarField > F3() const
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
virtual bool read()
Re-read model coefficients if they have changed.
Namespace for OpenFOAM.
BasicTurbulenceModel::transportModel transportModel
virtual void correctNut()
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.