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-2015 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::kOmegaSSTBase
26 
27 Description
28  Base class 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  kOmegaSSTBaseCoeffs
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  kOmegaSSTBase.C
97 
98 \*---------------------------------------------------------------------------*/
99 
100 #ifndef kOmegaSSTBase_H
101 #define kOmegaSSTBase_H
102 
103 #include "RASModel.H"
104 #include "eddyViscosity.H"
105 
106 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
107 
108 namespace Foam
109 {
110 
111 /*---------------------------------------------------------------------------*\
112  Class kOmegaSSTBase Declaration
113 \*---------------------------------------------------------------------------*/
114 
115 template<class BasicEddyViscosityModel>
116 class kOmegaSSTBase
117 :
118  public BasicEddyViscosityModel
119 {
120  // Private Member Functions
121 
122  // Disallow default bitwise copy construct and assignment
125 
126 
127 protected:
128 
129  // Protected data
130 
131  // Model coefficients
132 
135 
138 
141 
144 
146 
150 
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_;
160 
163 
164 
165  // Protected 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  }
181 
183  {
184  return blend(F1, alphaK1_, alphaK2_);
185  }
186 
188  {
189  return blend(F1, alphaOmega1_, alphaOmega2_);
190  }
191 
193  {
194  return blend(F1, beta1_, beta2_);
195  }
196 
198  {
199  return blend(F1, gamma1_, gamma2_);
200  }
201 
202  virtual void correctNut(const volScalarField& S2);
203 
204  virtual void correctNut();
205  virtual tmp<fvScalarMatrix> kSource() const;
206  virtual tmp<fvScalarMatrix> omegaSource() const;
207  virtual tmp<fvScalarMatrix> Qsas
208  (
209  const volScalarField& S2,
210  const volScalarField& gamma,
211  const volScalarField& beta
212  ) const;
213 
214 
215 public:
216 
217  typedef typename BasicEddyViscosityModel::alphaField alphaField;
218  typedef typename BasicEddyViscosityModel::rhoField rhoField;
219  typedef typename BasicEddyViscosityModel::transportModel transportModel;
220 
221 
222  // Constructors
223 
224  //- Construct from components
226  (
227  const word& type,
228  const alphaField& alpha,
229  const rhoField& rho,
230  const volVectorField& U,
231  const surfaceScalarField& alphaRhoPhi,
232  const surfaceScalarField& phi,
233  const transportModel& transport,
234  const word& propertiesName = turbulenceModel::propertiesName
235  );
236 
237 
238  //- Destructor
239  virtual ~kOmegaSSTBase()
240  {}
241 
242 
243  // Member Functions
244 
245  //- Re-read model coefficients if they have changed
246  virtual bool read();
247 
248  //- Return the effective diffusivity for k
250  {
251  return tmp<volScalarField>
252  (
253  new volScalarField("DkEff", alphaK(F1)*this->nut_ + this->nu())
254  );
255  }
256 
257  //- Return the effective diffusivity for omega
259  {
260  return tmp<volScalarField>
261  (
262  new volScalarField
263  (
264  "DomegaEff",
265  alphaOmega(F1)*this->nut_ + this->nu()
266  )
267  );
268  }
269 
270  //- Return the turbulence kinetic energy
271  virtual tmp<volScalarField> k() const
272  {
273  return k_;
274  }
275 
276  //- Return the turbulence kinetic energy dissipation rate
277  virtual tmp<volScalarField> epsilon() const
278  {
279  return tmp<volScalarField>
280  (
281  new volScalarField
282  (
283  IOobject
284  (
285  "epsilon",
286  this->mesh_.time().timeName(),
287  this->mesh_
288  ),
291  )
292  );
293  }
294 
295  //- Return the turbulence kinetic energy dissipation rate
296  virtual tmp<volScalarField> omega() const
297  {
298  return omega_;
299  }
300 
301  //- Solve the turbulence equations and correct the turbulence viscosity
302  virtual void correct();
303 };
304 
305 
306 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307 
308 } // End namespace Foam
309 
310 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 
312 #ifdef NoRepository
313 # include "kOmegaSSTBase.C"
314 #endif
315 
316 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
317 
318 #endif
319 
320 // ************************************************************************* //
Foam::kOmegaSSTBase::F2
tmp< volScalarField > F2() const
Definition: kOmegaSSTBase.C:68
Foam::kOmegaSSTBase::Qsas
virtual tmp< fvScalarMatrix > Qsas(const volScalarField &S2, const volScalarField &gamma, const volScalarField &beta) const
Definition: kOmegaSSTBase.C:160
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::kOmegaSSTBase::omega
virtual tmp< volScalarField > omega() const
Return the turbulence kinetic energy dissipation rate.
Definition: kOmegaSSTBase.H:295
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
Foam::kOmegaSSTBase::rhoField
BasicEddyViscosityModel::rhoField rhoField
Definition: kOmegaSSTBase.H:217
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
kOmegaSSTBase.C
Foam::kOmegaSSTBase::alphaK2_
dimensionedScalar alphaK2_
Definition: kOmegaSSTBase.H:133
Foam::kOmegaSSTBase::~kOmegaSSTBase
virtual ~kOmegaSSTBase()
Destructor.
Definition: kOmegaSSTBase.H:238
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::kOmegaSSTBase::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: kOmegaSSTBase.C:131
Foam::kOmegaSSTBase::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: kOmegaSSTBase.H:270
Foam::kOmegaSSTBase::F1
tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Definition: kOmegaSSTBase.C:39
Foam::kOmegaSSTBase::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: kOmegaSSTBase.H:276
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::kOmegaSSTBase::omegaSource
virtual tmp< fvScalarMatrix > omegaSource() const
Definition: kOmegaSSTBase.C:145
Foam::kOmegaSSTBase::operator=
kOmegaSSTBase & operator=(const kOmegaSSTBase &)
Foam::kOmegaSSTBase::omega_
volScalarField omega_
Definition: kOmegaSSTBase.H:161
Foam::kOmegaSSTBase::alphaOmega1_
dimensionedScalar alphaOmega1_
Definition: kOmegaSSTBase.H:135
psi1
psi1
Definition: TEqns.H:34
Foam::kOmegaSSTBase::beta
tmp< volScalarField > beta(const volScalarField &F1) const
Definition: kOmegaSSTBase.H:191
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::kOmegaSSTBase::beta1_
dimensionedScalar beta1_
Definition: kOmegaSSTBase.H:141
Foam::kOmegaSSTBase::F3_
Switch F3_
Definition: kOmegaSSTBase.H:150
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::kOmegaSSTBase::alphaOmega
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
Definition: kOmegaSSTBase.H:186
Foam::kOmegaSSTBase::betaStar_
dimensionedScalar betaStar_
Definition: kOmegaSSTBase.H:144
Foam::kOmegaSSTBase::alphaOmega2_
dimensionedScalar alphaOmega2_
Definition: kOmegaSSTBase.H:136
Foam::kOmegaSSTBase::DomegaEff
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Return the effective diffusivity for omega.
Definition: kOmegaSSTBase.H:257
Foam::kOmegaSSTBase::blend
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
Definition: kOmegaSSTBase.H:172
U
U
Definition: pEqn.H:46
Foam::kOmegaSSTBase::correctNut
virtual void correctNut()
Definition: kOmegaSSTBase.C:124
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
eddyViscosity.H
Foam::kOmegaSSTBase::gamma2_
dimensionedScalar gamma2_
Definition: kOmegaSSTBase.H:139
Foam::kOmegaSSTBase::kOmegaSSTBase
kOmegaSSTBase(const kOmegaSSTBase &)
Foam::kOmegaSSTBase::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTBase.C:357
psi2
psi2
Definition: TEqns.H:35
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
Foam::kOmegaSSTBase::c1_
dimensionedScalar c1_
Definition: kOmegaSSTBase.H:148
RASModel.H
Foam::kOmegaSSTBase::beta2_
dimensionedScalar beta2_
Definition: kOmegaSSTBase.H:142
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::kOmegaSSTBase::b1_
dimensionedScalar b1_
Definition: kOmegaSSTBase.H:147
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::kOmegaSSTBase::y_
const volScalarField & y_
Wall distance.
Definition: kOmegaSSTBase.H:158
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99
Foam::kOmegaSSTBase::DkEff
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.
Definition: kOmegaSSTBase.H:248
Foam::kOmegaSSTBase::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTBase.C:385
rho
rho
Definition: pEqn.H:3
Foam::kOmegaSSTBase::F23
tmp< volScalarField > F23() const
Definition: kOmegaSSTBase.C:98
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
Foam::GeometricField::GeometricBoundaryField::types
wordList types() const
Return a list of the patch types.
Definition: GeometricBoundaryField.C:534
Foam::kOmegaSSTBase::alphaK
tmp< volScalarField > alphaK(const volScalarField &F1) const
Definition: kOmegaSSTBase.H:181
Foam::kOmegaSSTBase::k_
volScalarField k_
Definition: kOmegaSSTBase.H:160
Foam::kOmegaSSTBase::alphaK1_
dimensionedScalar alphaK1_
Definition: kOmegaSSTBase.H:132
Foam::kOmegaSSTBase::transportModel
BasicEddyViscosityModel::transportModel transportModel
Definition: kOmegaSSTBase.H:218
Foam::kOmegaSSTBase::gamma1_
dimensionedScalar gamma1_
Definition: kOmegaSSTBase.H:138
Foam::kOmegaSSTBase::gamma
tmp< volScalarField > gamma(const volScalarField &F1) const
Definition: kOmegaSSTBase.H:196
Foam::kOmegaSSTBase::alphaField
BasicEddyViscosityModel::alphaField alphaField
Definition: kOmegaSSTBase.H:216
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::kOmegaSSTBase::F3
tmp< volScalarField > F3() const
Definition: kOmegaSSTBase.C:85
Foam::kOmegaSSTBase
Base class implementation of the k-omega-SST turbulence model for incompressible and compressible flo...
Definition: kOmegaSSTBase.H:115
Foam::kOmegaSSTBase::a1_
dimensionedScalar a1_
Definition: kOmegaSSTBase.H:146