kOmegaSSTLM.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) 2016-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::RASModels::kOmegaSSTLM
26 
27 Description
28  Langtry-Menter 4-equation transitional SST model
29  based on the k-omega-SST RAS model.
30 
31  References:
32  \verbatim
33  Langtry, R. B., & Menter, F. R. (2009).
34  Correlation-based transition modeling for unstructured parallelized
35  computational fluid dynamics codes.
36  AIAA journal, 47(12), 2894-2906.
37 
38  Menter, F. R., Langtry, R., & Volker, S. (2006).
39  Transition modelling for general purpose CFD codes.
40  Flow, turbulence and combustion, 77(1-4), 277-303.
41 
42  Langtry, R. B. (2006).
43  A correlation-based transition model using local variables for
44  unstructured parallelized CFD codes.
45  Phd. Thesis, Universität Stuttgart.
46  \endverbatim
47 
48  The model coefficients are
49  \verbatim
50  kOmegaSSTCoeffs
51  {
52  // Default SST coefficients
53  alphaK1 0.85;
54  alphaK2 1;
55  alphaOmega1 0.5;
56  alphaOmega2 0.856;
57  beta1 0.075;
58  beta2 0.0828;
59  betaStar 0.09;
60  gamma1 5/9;
61  gamma2 0.44;
62  a1 0.31;
63  b1 1;
64  c1 10;
65  F3 no;
66 
67  // Default LM coefficients
68  ca1 2;
69  ca2 0.06;
70  ce1 1;
71  ce2 50;
72  cThetat 0.03;
73  sigmaThetat 2;
74 
75  lambdaErr 1e-6;
76  maxLambdaIter 10;
77  }
78  \endverbatim
79 
80 SourceFiles
81  kOmegaSSTLM.C
82 
83 \*---------------------------------------------------------------------------*/
84 
85 #ifndef kOmegaSSTLM_H
86 #define kOmegaSSTLM_H
87 
88 #include "kOmegaSST.H"
89 
90 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91 
92 namespace Foam
93 {
94 namespace RASModels
95 {
96 
97 /*---------------------------------------------------------------------------*\
98  Class kOmegaSSTLM Declaration
99 \*---------------------------------------------------------------------------*/
100 
101 template<class BasicTurbulenceModel>
102 class kOmegaSSTLM
103 :
104  public kOmegaSST<BasicTurbulenceModel>
105 {
106  // Private Member Functions
107 
108  // Disallow default bitwise copy construct and assignment
109  kOmegaSSTLM(const kOmegaSSTLM&);
110  void operator=(const kOmegaSSTLM&);
111 
112 
113 protected:
114 
115  // Protected data
116 
117  // Model constants
127 
128  //- Convergence criterion for the lambda/thetat loop
129  scalar lambdaErr_;
130 
131  //- Maximum number of iterations to converge the lambda/thetat loop
133 
134  //- Stabilization for division by the magnitude of the velocity
136 
137 
138  // Fields
139 
140  //- Transition onset momentum-thickness Reynolds number
142 
143  //- Intermittency
145 
146  //- Effective intermittency
148 
149 
150  // Protected Member Functions
151 
152  //- Modified form of the k-omega SST F1 function
153  virtual tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
154 
155  //- Modified form of the k-omega SST k production rate
157  (
159  ) const;
160 
161  //- Modified form of the k-omega SST epsilon/k
163  (
166  ) const;
167 
168  //- Freestream blending-function
170  (
171  const volScalarField::Internal& Us,
172  const volScalarField::Internal& Omega,
174  ) const;
175 
176  //- Empirical correlation for critical Reynolds number where the
177  // intermittency first starts to increase in the boundary layer
179 
180  //- Empirical correlation that controls the length of the
181  // transition region
183  (
185  ) const;
186 
187  //- Transition onset location control function
189  (
190  const volScalarField::Internal& Rev,
192  const volScalarField::Internal& RT
193  ) const;
194 
195  //- Return the transition onset momentum-thickness Reynolds number
196  // (based on freestream conditions)
198  (
199  const volScalarField::Internal& Us,
200  const volScalarField::Internal& dUsds,
202  ) const;
203 
204  //- Solve the turbulence equations and correct the turbulence viscosity
206 
207 
208 public:
210  typedef typename BasicTurbulenceModel::alphaField alphaField;
211  typedef typename BasicTurbulenceModel::rhoField rhoField;
212  typedef typename BasicTurbulenceModel::transportModel transportModel;
213 
214 
215  //- Runtime type information
216  TypeName("kOmegaSSTLM");
217 
218 
219  // Constructors
220 
221  //- Construct from components
223  (
224  const alphaField& alpha,
225  const rhoField& rho,
226  const volVectorField& U,
227  const surfaceScalarField& alphaRhoPhi,
228  const surfaceScalarField& phi,
229  const transportModel& transport,
230  const word& propertiesName = turbulenceModel::propertiesName,
231  const word& type = typeName
232  );
233 
234 
235  //- Destructor
236  virtual ~kOmegaSSTLM()
237  {}
238 
239 
240  // Member Functions
241 
242  //- Re-read model coefficients if they have changed
243  virtual bool read();
244 
245  //- Access function transition onset momentum-thickness Reynolds number
246  const volScalarField& ReThetat() const
247  {
248  return ReThetat_;
249  }
250 
251  //- Access function to intermittency
252  const volScalarField& gammaInt() const
253  {
254  return gammaInt_;
255  }
256 
257  //- Return the effective diffusivity for transition onset
258  // momentum-thickness Reynolds number
260  {
261  return tmp<volScalarField>
262  (
263  new volScalarField
264  (
265  "DReThetatEff",
266  sigmaThetat_*(this->nut_ + this->nu())
267  )
268  );
269  }
270 
271  //- Return the effective diffusivity for intermittency
273  {
274  return tmp<volScalarField>
275  (
276  new volScalarField
277  (
278  "DgammaIntEff",
279  this->nut_ + this->nu()
280  )
281  );
282  }
283 
284  //- Solve the turbulence equations and correct the turbulence viscosity
285  virtual void correct();
286 };
287 
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 } // End namespace RASModels
292 } // End namespace Foam
293 
294 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 
296 #ifdef NoRepository
297  #include "kOmegaSSTLM.C"
298 #endif
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
302 #endif
303 
304 // ************************************************************************* //
volScalarField::Internal gammaIntEff_
Effective intermittency.
Definition: kOmegaSSTLM.H:146
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
surfaceScalarField & phi
tmp< volScalarField::Internal > ReThetat0(const volScalarField::Internal &Us, const volScalarField::Internal &dUsds, const volScalarField::Internal &nu) const
Return the transition onset momentum-thickness Reynolds number.
Definition: kOmegaSSTLM.C:213
dimensionedScalar ca2_
Definition: kOmegaSSTLM.H:119
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< volScalarField > DReThetatEff() const
Return the effective diffusivity for transition onset.
Definition: kOmegaSSTLM.H:258
tmp< volScalarField > DgammaIntEff() const
Return the effective diffusivity for intermittency.
Definition: kOmegaSSTLM.H:271
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
Modified form of the k-omega SST epsilon/k.
Definition: kOmegaSSTLM.C:62
dimensionedScalar ca1_
Definition: kOmegaSSTLM.H:118
const volScalarField & gammaInt() const
Access function to intermittency.
Definition: kOmegaSSTLM.H:251
virtual ~kOmegaSSTLM()
Destructor.
Definition: kOmegaSSTLM.H:235
volScalarField ReThetat_
Transition onset momentum-thickness Reynolds number.
Definition: kOmegaSSTLM.H:140
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
tmp< volScalarField::Internal > ReThetac() const
Empirical correlation for critical Reynolds number where the.
Definition: kOmegaSSTLM.C:109
const dimensionedScalar deltaU_
Stabilization for division by the magnitude of the velocity.
Definition: kOmegaSSTLM.H:134
const volScalarField & ReThetat() const
Access function transition onset momentum-thickness Reynolds number.
Definition: kOmegaSSTLM.H:245
dimensionedScalar cThetat_
Definition: kOmegaSSTLM.H:124
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTLM.H:209
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTLM.C:607
Specialisation for RAS of the generic kOmegaSSTBase base class. For more information, see Description of kOmegaSSTBase.H.
Definition: kOmegaSST.H:57
TypeName("kOmegaSSTLM")
Runtime type information.
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTLM.C:486
dimensionedScalar sigmaThetat_
Definition: kOmegaSSTLM.H:125
tmp< volScalarField::Internal > Fonset(const volScalarField::Internal &Rev, const volScalarField::Internal &ReThetac, const volScalarField::Internal &RT) const
Transition onset location control function.
Definition: kOmegaSSTLM.C:327
tmp< volScalarField::Internal > Flength(const volScalarField::Internal &nu) const
Empirical correlation that controls the length of the.
Definition: kOmegaSSTLM.C:150
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
dimensionedScalar ce2_
Definition: kOmegaSSTLM.H:122
dimensionedScalar ce1_
Definition: kOmegaSSTLM.H:121
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTLM.H:210
U
Definition: pEqn.H:72
label maxLambdaIter_
Maximum number of iterations to converge the lambda/thetat loop.
Definition: kOmegaSSTLM.H:131
Langtry-Menter 4-equation transitional SST model based on the k-omega-SST RAS model.
Definition: kOmegaSSTLM.H:101
scalar lambdaErr_
Convergence criterion for the lambda/thetat loop.
Definition: kOmegaSSTLM.H:128
void correctReThetatGammaInt()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTLM.C:509
volScalarField gammaInt_
Intermittency.
Definition: kOmegaSSTLM.H:143
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTLM.H:211
tmp< volScalarField::Internal > Fthetat(const volScalarField::Internal &Us, const volScalarField::Internal &Omega, const volScalarField::Internal &nu) const
Freestream blending-function.
Definition: kOmegaSSTLM.C:75
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Modified form of the k-omega SST F1 function.
Definition: kOmegaSSTLM.C:39
volScalarField & nu
Namespace for OpenFOAM.
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Modified form of the k-omega SST k production rate.
Definition: kOmegaSSTLM.C:52