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-2020 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 BasicMomentumTransportModel>
102 class kOmegaSSTLM
103 :
104  public kOmegaSST<BasicMomentumTransportModel>
105 {
106 protected:
107 
108  // Protected data
109 
110  // Model constants
120 
121  //- Convergence criterion for the lambda/thetat loop
122  scalar lambdaErr_;
123 
124  //- Maximum number of iterations to converge the lambda/thetat loop
126 
127  //- Stabilization for division by the magnitude of the velocity
129 
130 
131  // Fields
132 
133  //- Transition onset momentum-thickness Reynolds number
135 
136  //- Intermittency
138 
139  //- Effective intermittency
141 
142 
143  // Protected Member Functions
144 
145  //- Modified form of the k-omega SST F1 function
146  virtual tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
147 
148  //- Modified form of the k-omega SST k production rate
150  (
152  ) const;
153 
154  //- Modified form of the k-omega SST epsilon/k
156  (
159  ) const;
160 
161  //- Freestream blending-function
163  (
164  const volScalarField::Internal& Us,
165  const volScalarField::Internal& Omega,
166  const volScalarField::Internal& nu
167  ) const;
168 
169  //- Empirical correlation for critical Reynolds number where the
170  // intermittency first starts to increase in the boundary layer
172 
173  //- Empirical correlation that controls the length of the
174  // transition region
176  (
177  const volScalarField::Internal& nu
178  ) const;
179 
180  //- Transition onset location control function
182  (
183  const volScalarField::Internal& Rev,
185  const volScalarField::Internal& RT
186  ) const;
187 
188  //- Return the transition onset momentum-thickness Reynolds number
189  // (based on freestream conditions)
191  (
192  const volScalarField::Internal& Us,
193  const volScalarField::Internal& dUsds,
194  const volScalarField::Internal& nu
195  ) const;
196 
197  //- Solve the turbulence equations and correct the turbulence viscosity
199 
200 
201 public:
203  typedef typename BasicMomentumTransportModel::alphaField alphaField;
204  typedef typename BasicMomentumTransportModel::rhoField rhoField;
205  typedef typename BasicMomentumTransportModel::transportModel transportModel;
206 
207 
208  //- Runtime type information
209  TypeName("kOmegaSSTLM");
210 
211 
212  // Constructors
213 
214  //- Construct from components
216  (
217  const alphaField& alpha,
218  const rhoField& rho,
219  const volVectorField& U,
220  const surfaceScalarField& alphaRhoPhi,
221  const surfaceScalarField& phi,
222  const transportModel& transport,
223  const word& type = typeName
224  );
225 
226  //- Disallow default bitwise copy construction
227  kOmegaSSTLM(const kOmegaSSTLM&) = delete;
228 
229 
230  //- Destructor
231  virtual ~kOmegaSSTLM()
232  {}
233 
234 
235  // Member Functions
236 
237  //- Re-read model coefficients if they have changed
238  virtual bool read();
239 
240  //- Access function transition onset momentum-thickness Reynolds number
241  const volScalarField& ReThetat() const
242  {
243  return ReThetat_;
244  }
245 
246  //- Access function to intermittency
247  const volScalarField& gammaInt() const
248  {
249  return gammaInt_;
250  }
251 
252  //- Return the effective diffusivity for transition onset
253  // momentum-thickness Reynolds number
255  {
256  return volScalarField::New
257  (
258  "DReThetatEff",
259  sigmaThetat_*(this->nut_ + this->nu())
260  );
261  }
262 
263  //- Return the effective diffusivity for intermittency
265  {
266  return volScalarField::New
267  (
268  "DgammaIntEff",
269  this->nut_ + this->nu()
270  );
271  }
272 
273  //- Solve the turbulence equations and correct the turbulence viscosity
274  virtual void correct();
275 
276 
277  // Member Operators
278 
279  //- Disallow default bitwise assignment
280  void operator=(const kOmegaSSTLM&) = delete;
281 };
282 
283 
284 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285 
286 } // End namespace RASModels
287 } // End namespace Foam
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 #ifdef NoRepository
292  #include "kOmegaSSTLM.C"
293 #endif
294 
295 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296 
297 #endif
298 
299 // ************************************************************************* //
BasicMomentumTransportModel::alphaField alphaField
Definition: kOmegaSSTLM.H:202
volScalarField::Internal gammaIntEff_
Effective intermittency.
Definition: kOmegaSSTLM.H:139
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
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:202
dimensionedScalar ca2_
Definition: kOmegaSSTLM.H:112
tmp< volScalarField > DReThetatEff() const
Return the effective diffusivity for transition onset.
Definition: kOmegaSSTLM.H:253
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
tmp< volScalarField > DgammaIntEff() const
Return the effective diffusivity for intermittency.
Definition: kOmegaSSTLM.H:263
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:63
kOmegaSSTLM(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: kOmegaSSTLM.C:338
dimensionedScalar ca1_
Definition: kOmegaSSTLM.H:111
phi
Definition: pEqn.H:104
const volScalarField & gammaInt() const
Access function to intermittency.
Definition: kOmegaSSTLM.H:246
virtual ~kOmegaSSTLM()
Destructor.
Definition: kOmegaSSTLM.H:230
volScalarField ReThetat_
Transition onset momentum-thickness Reynolds number.
Definition: kOmegaSSTLM.H:133
tmp< volScalarField::Internal > ReThetac() const
Empirical correlation for critical Reynolds number where the.
Definition: kOmegaSSTLM.C:107
const dimensionedScalar deltaU_
Stabilization for division by the magnitude of the velocity.
Definition: kOmegaSSTLM.H:127
const volScalarField & ReThetat() const
Access function transition onset momentum-thickness Reynolds number.
Definition: kOmegaSSTLM.H:240
dimensionedScalar cThetat_
Definition: kOmegaSSTLM.H:117
A class for handling words, derived from string.
Definition: word.H:59
BasicMomentumTransportModel::rhoField rhoField
Definition: kOmegaSSTLM.H:203
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTLM.C:586
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:465
dimensionedScalar sigmaThetat_
Definition: kOmegaSSTLM.H:118
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:311
tmp< volScalarField::Internal > Flength(const volScalarField::Internal &nu) const
Empirical correlation that controls the length of the.
Definition: kOmegaSSTLM.C:143
dimensionedScalar ce2_
Definition: kOmegaSSTLM.H:115
dimensionedScalar ce1_
Definition: kOmegaSSTLM.H:114
BasicMomentumTransportModel::transportModel transportModel
Definition: kOmegaSSTLM.H:204
U
Definition: pEqn.H:72
label maxLambdaIter_
Maximum number of iterations to converge the lambda/thetat loop.
Definition: kOmegaSSTLM.H:124
Langtry-Menter 4-equation transitional SST model based on the k-omega-SST RAS model.
Definition: kOmegaSSTLM.H:101
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
scalar lambdaErr_
Convergence criterion for the lambda/thetat loop.
Definition: kOmegaSSTLM.H:121
void correctReThetatGammaInt()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTLM.C:488
volScalarField gammaInt_
Intermittency.
Definition: kOmegaSSTLM.H:136
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< volScalarField::Internal > Fthetat(const volScalarField::Internal &Us, const volScalarField::Internal &Omega, const volScalarField::Internal &nu) const
Freestream blending-function.
Definition: kOmegaSSTLM.C:76
A class for managing temporary objects.
Definition: PtrList.H:53
void operator=(const kOmegaSSTLM &)=delete
Disallow default bitwise assignment.
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Modified form of the k-omega SST F1 function.
Definition: kOmegaSSTLM.C:39
const dimensionedScalar & G
Newtonian constant of gravitation.
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