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-2023 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 MomentumTransportModel, class BasicMomentumTransportModel>
113 class kOmegaSST
114 :
115  public MomentumTransportModel
116 {
117 protected:
118 
119  // Protected data
120 
121  // Model coefficients
122 
125 
128 
131 
134 
136 
140 
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_;
150 
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  (
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  (
195  ) const
196  {
197  return blend(F1, beta1_, beta2_);
198  }
199 
201  (
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  (
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:
242 
243  typedef typename BasicMomentumTransportModel::alphaField alphaField;
244  typedef typename BasicMomentumTransportModel::rhoField rhoField;
245 
246 
247  // Constructors
248 
249  //- Construct from components
250  kOmegaSST
251  (
252  const word& type,
253  const alphaField& alpha,
254  const rhoField& rho,
255  const volVectorField& U,
256  const surfaceScalarField& alphaRhoPhi,
257  const surfaceScalarField& phi,
258  const viscosity& viscosity
259  );
260 
261  //- Disallow default bitwise copy construction
262  kOmegaSST(const kOmegaSST&) = delete;
263 
264 
265  //- Destructor
266  virtual ~kOmegaSST()
267  {}
268 
269 
270  // Member Functions
271 
272  //- Re-read model coefficients if they have changed
273  virtual bool read();
274 
275  //- Return the effective diffusivity for k
277  {
278  return volScalarField::New
279  (
280  "DkEff",
281  alphaK(F1)*this->nut_ + this->nu()
282  );
283  }
284 
285  //- Return the effective diffusivity for omega
287  {
288  return volScalarField::New
289  (
290  "DomegaEff",
291  alphaOmega(F1)*this->nut_ + this->nu()
292  );
293  }
294 
295  //- Return the turbulence kinetic energy
296  virtual tmp<volScalarField> k() const
297  {
298  return k_;
299  }
300 
301  //- Return the turbulence kinetic energy dissipation rate
302  virtual tmp<volScalarField> epsilon() const
303  {
304  return volScalarField::New
305  (
306  "epsilon",
308  );
309  }
310 
311  //- Return the turbulence kinetic energy dissipation rate
312  virtual tmp<volScalarField> omega() const
313  {
314  return omega_;
315  }
316 
317  //- Solve the turbulence equations and correct the turbulence viscosity
318  virtual void correct();
319 
320 
321  // Member Operators
322 
323  //- Disallow default bitwise assignment
324  void operator=(const kOmegaSST&) = delete;
325 };
326 
327 
328 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329 
330 } // End namespace Foam
331 
332 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333 #ifdef NoRepository
334  #include "kOmegaSSTBase.C"
335 #endif
336 
337 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
338 #endif
339 
340 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows.
BasicMomentumTransportModel::alphaField alphaField
virtual tmp< volScalarField > F23() const
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
dimensionedScalar beta2_
volScalarField k_
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Return the effective diffusivity for omega.
virtual tmp< volScalarField > F2() const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual tmp< volScalarField > omega() const
Return the turbulence kinetic energy dissipation rate.
dimensionedScalar alphaOmega2_
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedScalar gamma2_
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
kOmegaSST(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
dimensionedScalar b1_
dimensionedScalar beta1_
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
dimensionedScalar c1_
dimensionedScalar alphaK2_
const volScalarField & y_
Wall distance.
dimensionedScalar a1_
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
virtual ~kOmegaSST()
Destructor.
dimensionedScalar alphaK1_
volScalarField omega_
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
virtual tmp< volScalarField > F3() const
virtual tmp< fvScalarMatrix > omegaSource() const
virtual void correctNut()
void operator=(const kOmegaSST &)=delete
Disallow default bitwise assignment.
virtual tmp< fvScalarMatrix > kSource() const
dimensionedScalar alphaOmega1_
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.
dimensionedScalar gamma1_
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
tmp< volScalarField > alphaK(const volScalarField &F1) const
dimensionedScalar betaStar_
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar G
Newtonian constant of gravitation.
Namespace for OpenFOAM.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488