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-2019 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 protected:
118 
119  // Protected data
120 
121  // Model coefficients
141  Switch F3_;
142 
143 
144  // Fields
145 
146  //- Wall distance
147  // Note: different to wall distance in parent RASModel
148  // which is for near-wall cells only
149  const volScalarField& y_;
153 
154 
155  // Protected Member Functions
156 
157  virtual tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
158  virtual tmp<volScalarField> F2() const;
159  virtual tmp<volScalarField> F3() const;
160  virtual tmp<volScalarField> F23() const;
161 
163  (
164  const volScalarField& F1,
165  const dimensionedScalar& psi1,
166  const dimensionedScalar& psi2
167  ) const
168  {
169  return F1*(psi1 - psi2) + psi2;
170  }
171 
173  (
174  const volScalarField::Internal& F1,
175  const dimensionedScalar& psi1,
176  const dimensionedScalar& psi2
177  ) const
178  {
179  return F1*(psi1 - psi2) + psi2;
180  }
182  tmp<volScalarField> alphaK(const volScalarField& F1) const
183  {
184  return blend(F1, alphaK1_, alphaK2_);
185  }
188  {
189  return blend(F1, alphaOmega1_, alphaOmega2_);
190  }
191 
193  (
194  const volScalarField::Internal& F1
195  ) const
196  {
197  return blend(F1, beta1_, beta2_);
198  }
199 
201  (
202  const volScalarField::Internal& F1
203  ) const
204  {
205  return blend(F1, gamma1_, gamma2_);
206  }
207 
208  virtual void correctNut
209  (
210  const volScalarField& S2,
211  const volScalarField& F2
212  );
213 
214  virtual void correctNut();
215 
216  //- Return k production rate
218  (
220  ) const;
221 
222  //- Return epsilon/k which for standard RAS is betaStar*omega
224  (
225  const volScalarField::Internal& F1,
226  const volScalarField::Internal& F2
227  ) const;
228 
229  virtual tmp<fvScalarMatrix> kSource() const;
230 
231  virtual tmp<fvScalarMatrix> omegaSource() const;
232 
233  virtual tmp<fvScalarMatrix> Qsas
234  (
235  const volScalarField::Internal& S2,
238  ) const;
239 
240 
241 public:
243  typedef typename BasicTurbulenceModel::alphaField alphaField;
244  typedef typename BasicTurbulenceModel::rhoField rhoField;
245  typedef typename BasicTurbulenceModel::transportModel transportModel;
246 
247 
248  // Constructors
249 
250  //- Construct from components
251  kOmegaSST
252  (
253  const word& type,
254  const alphaField& alpha,
255  const rhoField& rho,
256  const volVectorField& U,
257  const surfaceScalarField& alphaRhoPhi,
258  const surfaceScalarField& phi,
259  const transportModel& transport,
260  const word& propertiesName = turbulenceModel::propertiesName
261  );
262 
263  //- Disallow default bitwise copy construction
264  kOmegaSST(const kOmegaSST&) = delete;
265 
266 
267  //- Destructor
268  virtual ~kOmegaSST()
269  {}
270 
271 
272  // Member Functions
273 
274  //- Re-read model coefficients if they have changed
275  virtual bool read();
276 
277  //- Return the effective diffusivity for k
278  tmp<volScalarField> DkEff(const volScalarField& F1) const
279  {
280  return volScalarField::New
281  (
282  "DkEff",
283  alphaK(F1)*this->nut_ + this->nu()
284  );
285  }
286 
287  //- Return the effective diffusivity for omega
289  {
290  return volScalarField::New
291  (
292  "DomegaEff",
293  alphaOmega(F1)*this->nut_ + this->nu()
294  );
295  }
296 
297  //- Return the turbulence kinetic energy
298  virtual tmp<volScalarField> k() const
299  {
300  return k_;
301  }
302 
303  //- Return the turbulence kinetic energy dissipation rate
304  virtual tmp<volScalarField> epsilon() const
305  {
306  return volScalarField::New
307  (
308  "epsilon",
309  betaStar_*k_*omega_,
310  omega_.boundaryField().types()
311  );
312  }
313 
314  //- Return the turbulence kinetic energy dissipation rate
315  virtual tmp<volScalarField> omega() const
316  {
317  return omega_;
318  }
319 
320  //- Solve the turbulence equations and correct the turbulence viscosity
321  virtual void correct();
322 
323 
324  // Member Operators
325 
326  //- Disallow default bitwise assignment
327  void operator=(const kOmegaSST&) = delete;
328 };
329 
330 
331 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
332 
333 } // End namespace Foam
334 
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336 #ifdef NoRepository
337  #include "kOmegaSSTBase.C"
338 #endif
339 
340 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
341 #endif
342 
343 // ************************************************************************* //
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
kOmegaSST(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
Construct from components.
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/any.
Definition: Switch.H:60
dimensionedScalar b1_
virtual tmp< fvScalarMatrix > omegaSource() const
volScalarField omega_
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
dimensionedScalar alphaOmega2_
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
const transportModel & transport() const
Access function to incompressible transport model.
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.
U
Definition: pEqn.H:72
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
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.
void operator=(const kOmegaSST &)=delete
Disallow default bitwise assignment.
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.
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()