kOmegaSSTBase.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 Description
28  Implementation of the k-omega-SST turbulence model for
29  incompressible and compressible flows.
30 
31  Turbulence model described in:
32  \verbatim
33  Menter, F. R. & Esch, T. (2001).
34  Elements of Industrial Heat Transfer Prediction.
35  16th Brazilian Congress of Mechanical Engineering (COBEM).
36  \endverbatim
37 
38  with updated coefficients from
39  \verbatim
40  Menter, F. R., Kuntz, M., and Langtry, R. (2003).
41  Ten Years of Industrial Experience with the SST Turbulence Model.
42  Turbulence, Heat and Mass Transfer 4, ed: K. Hanjalic, Y. Nagano,
43  & M. Tummers, Begell House, Inc., 625 - 632.
44  \endverbatim
45 
46  but with the consistent production terms from the 2001 paper as form in the
47  2003 paper is a typo, see
48  \verbatim
49  http://turbmodels.larc.nasa.gov/sst.html
50  \endverbatim
51 
52  and the addition of the optional F3 term for rough walls from
53  \verbatim
54  Hellsten, A. (1998).
55  "Some Improvements in Menter’s k-omega-SST turbulence model"
56  29th AIAA Fluid Dynamics Conference, AIAA-98-2554.
57  \endverbatim
58 
59  Note that this implementation is written in terms of alpha diffusion
60  coefficients rather than the more traditional sigma (alpha = 1/sigma) so
61  that the blending can be applied to all coefficuients in a consistent
62  manner. The paper suggests that sigma is blended but this would not be
63  consistent with the blending of the k-epsilon and k-omega models.
64 
65  Also note that the error in the last term of equation (2) relating to
66  sigma has been corrected.
67 
68  Wall-functions are applied in this implementation by using equations (14)
69  to specify the near-wall omega as appropriate.
70 
71  The blending functions (15) and (16) are not currently used because of the
72  uncertainty in their origin, range of applicability and that if y+ becomes
73  sufficiently small blending u_tau in this manner clearly becomes nonsense.
74 
75  The default model coefficients are
76  \verbatim
77  kOmegaSSTCoeffs
78  {
79  alphaK1 0.85;
80  alphaK2 1.0;
81  alphaOmega1 0.5;
82  alphaOmega2 0.856;
83  beta1 0.075;
84  beta2 0.0828;
85  betaStar 0.09;
86  gamma1 5/9;
87  gamma2 0.44;
88  a1 0.31;
89  b1 1.0;
90  c1 10.0;
91  F3 no;
92  }
93  \endverbatim
94 
95 SourceFiles
96  kOmegaSST.C
97 
98 \*---------------------------------------------------------------------------*/
99 
100 #ifndef kOmegaSSTBase_H
101 #define kOmegaSSTBase_H
102 
103 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104 
105 namespace Foam
106 {
107 
108 /*---------------------------------------------------------------------------*\
109  Class kOmegaSST Declaration
110 \*---------------------------------------------------------------------------*/
111 
112 template<class TurbulenceModel, class BasicTurbulenceModel>
113 class kOmegaSST
114 :
115  public TurbulenceModel
116 {
117  // Private Member Functions
118 
119  // Disallow default bitwise copy construct and assignment
120  kOmegaSST(const kOmegaSST&);
121  void operator=(const kOmegaSST&);
122 
123 
124 protected:
125 
126  // Protected data
127 
128  // Model coefficients
148  Switch F3_;
149 
150 
151  // Fields
152 
153  //- Wall distance
154  // Note: different to wall distance in parent RASModel
155  // which is for near-wall cells only
156  const volScalarField& y_;
160 
161 
162  // Protected Member Functions
163 
164  virtual tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
165  virtual tmp<volScalarField> F2() const;
166  virtual tmp<volScalarField> F3() const;
167  virtual tmp<volScalarField> F23() const;
168 
170  (
171  const volScalarField& F1,
172  const dimensionedScalar& psi1,
173  const dimensionedScalar& psi2
174  ) const
175  {
176  return F1*(psi1 - psi2) + psi2;
177  }
178 
180  (
181  const volScalarField::Internal& F1,
182  const dimensionedScalar& psi1,
183  const dimensionedScalar& psi2
184  ) const
185  {
186  return F1*(psi1 - psi2) + psi2;
187  }
189  tmp<volScalarField> alphaK(const volScalarField& F1) const
190  {
191  return blend(F1, alphaK1_, alphaK2_);
192  }
195  {
196  return blend(F1, alphaOmega1_, alphaOmega2_);
197  }
198 
200  (
201  const volScalarField::Internal& F1
202  ) const
203  {
204  return blend(F1, beta1_, beta2_);
205  }
206 
208  (
209  const volScalarField::Internal& F1
210  ) const
211  {
212  return blend(F1, gamma1_, gamma2_);
213  }
214 
215  virtual void correctNut
216  (
217  const volScalarField& S2,
218  const volScalarField& F2
219  );
220 
221  virtual void correctNut();
222 
223  //- Return k production rate
225  (
227  ) const;
228 
229  //- Return epsilon/k which for standard RAS is betaStar*omega
231  (
232  const volScalarField::Internal& F1,
233  const volScalarField::Internal& F2
234  ) const;
235 
236  virtual tmp<fvScalarMatrix> kSource() const;
237 
238  virtual tmp<fvScalarMatrix> omegaSource() const;
239 
240  virtual tmp<fvScalarMatrix> Qsas
241  (
242  const volScalarField::Internal& S2,
245  ) const;
246 
247 
248 public:
250  typedef typename BasicTurbulenceModel::alphaField alphaField;
251  typedef typename BasicTurbulenceModel::rhoField rhoField;
252  typedef typename BasicTurbulenceModel::transportModel transportModel;
253 
254 
255  // Constructors
256 
257  //- Construct from components
258  kOmegaSST
259  (
260  const word& type,
261  const alphaField& alpha,
262  const rhoField& rho,
263  const volVectorField& U,
264  const surfaceScalarField& alphaRhoPhi,
265  const surfaceScalarField& phi,
266  const transportModel& transport,
267  const word& propertiesName = turbulenceModel::propertiesName
268  );
269 
270 
271  //- Destructor
272  virtual ~kOmegaSST()
273  {}
274 
275 
276  // Member Functions
277 
278  //- Re-read model coefficients if they have changed
279  virtual bool read();
280 
281  //- Return the effective diffusivity for k
282  tmp<volScalarField> DkEff(const volScalarField& F1) const
283  {
284  return tmp<volScalarField>
285  (
286  new volScalarField("DkEff", alphaK(F1)*this->nut_ + this->nu())
287  );
288  }
289 
290  //- Return the effective diffusivity for omega
292  {
293  return tmp<volScalarField>
294  (
295  new volScalarField
296  (
297  "DomegaEff",
298  alphaOmega(F1)*this->nut_ + this->nu()
299  )
300  );
301  }
302 
303  //- Return the turbulence kinetic energy
304  virtual tmp<volScalarField> k() const
305  {
306  return k_;
307  }
308 
309  //- Return the turbulence kinetic energy dissipation rate
310  virtual tmp<volScalarField> epsilon() const
311  {
312  return tmp<volScalarField>
313  (
314  new volScalarField
315  (
316  IOobject
317  (
318  "epsilon",
319  this->mesh_.time().timeName(),
320  this->mesh_
321  ),
322  betaStar_*k_*omega_,
323  omega_.boundaryField().types()
324  )
325  );
326  }
327 
328  //- Return the turbulence kinetic energy dissipation rate
329  virtual tmp<volScalarField> omega() const
330  {
331  return omega_;
332  }
333 
334  //- Solve the turbulence equations and correct the turbulence viscosity
335  virtual void correct();
336 };
337 
338 
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340 
341 } // End namespace Foam
342 
343 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
344 #ifdef NoRepository
345  #include "kOmegaSSTBase.C"
346 #endif
347 
348 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
349 #endif
350 
351 // ************************************************************************* //
virtual ~kOmegaSST()
Destructor.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
psi2
Definition: TEqns.H:35
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
tmp< volScalarField > alphaK(const volScalarField &F1) const
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.
dimensionedScalar a1_
dimensionedScalar gamma2_
surfaceScalarField & phi
dimensionedScalar c1_
dimensionedScalar betaStar_
dimensionedScalar beta1_
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const dimensionedScalar G
Newtonian constant of gravitation.
BasicTurbulenceModel::alphaField alphaField
Templated abstract base class for turbulence models.
virtual tmp< volScalarField > F2() const
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar gamma1_
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_
virtual tmp< fvScalarMatrix > omegaSource() const
volScalarField omega_
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
dimensionedScalar alphaOmega2_
const transportModel & transport() const
Access function to incompressible transport model.
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
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 > F23() const
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
const alphaField & alpha() const
Access function to phase fraction.
dimensionedScalar alphaOmega1_
dimensionedScalar alphaK2_
virtual tmp< volScalarField > omega() const
Return the turbulence kinetic energy dissipation rate.
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
dimensionedScalar alphaK1_
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
U
Definition: pEqn.H:72
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual tmp< volScalarField > F3() const
volScalarField k_
virtual tmp< fvScalarMatrix > kSource() const
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Return the effective diffusivity for omega.
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
dimensionedScalar beta2_
A class for managing temporary objects.
Definition: PtrList.H:53
const volScalarField & y_
Wall distance.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual bool read()
Re-read model coefficients if they have changed.
Namespace for OpenFOAM.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
virtual void correctNut()