kOmegaSSTLM.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::RASModels::kOmegaSSTLM
29
30Group
31 grpLESTurbulence
32
33Description
34 Langtry-Menter 4-equation transitional SST model
35 based on the k-omega-SST RAS model.
36
37 References:
38 \verbatim
39 Langtry, R. B., & Menter, F. R. (2009).
40 Correlation-based transition modeling for unstructured parallelized
41 computational fluid dynamics codes.
42 AIAA journal, 47(12), 2894-2906.
43
44 Menter, F. R., Langtry, R., & Volker, S. (2006).
45 Transition modelling for general purpose CFD codes.
46 Flow, turbulence and combustion, 77(1-4), 277-303.
47
48 Langtry, R. B. (2006).
49 A correlation-based transition model using local variables for
50 unstructured parallelized CFD codes.
51 Phd. Thesis, Universität Stuttgart.
52 \endverbatim
53
54 The model coefficients are
55 \verbatim
56 kOmegaSSTCoeffs
57 {
58 // Default SST coefficients
59 alphaK1 0.85;
60 alphaK2 1;
61 alphaOmega1 0.5;
62 alphaOmega2 0.856;
63 beta1 0.075;
64 beta2 0.0828;
65 betaStar 0.09;
66 gamma1 5/9;
67 gamma2 0.44;
68 a1 0.31;
69 b1 1;
70 c1 10;
71 F3 no;
72
73 // Default LM coefficients
74 ca1 2;
75 ca2 0.06;
76 ce1 1;
77 ce2 50;
78 cThetat 0.03;
79 sigmaThetat 2;
80
81 lambdaErr 1e-6;
82 maxLambdaIter 10;
83 }
84 \endverbatim
85
86SourceFiles
87 kOmegaSSTLM.C
88
89\*---------------------------------------------------------------------------*/
90
91#ifndef kOmegaSSTLM_H
92#define kOmegaSSTLM_H
93
94#include "kOmegaSST.H"
95
96// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97
98namespace Foam
99{
100namespace RASModels
101{
102
103/*---------------------------------------------------------------------------*\
104 Class kOmegaSSTLM Declaration
105\*---------------------------------------------------------------------------*/
106
107template<class BasicTurbulenceModel>
108class kOmegaSSTLM
109:
110 public kOmegaSST<BasicTurbulenceModel>
111{
112 // Private Member Functions
113
114 //- No copy construct
115 kOmegaSSTLM(const kOmegaSSTLM&) = delete;
116
117 //- No copy assignment
118 void operator=(const kOmegaSSTLM&) = delete;
119
120
121protected:
122
123 // Protected data
124
125 // Model constants
135
136 //- Convergence criterion for the lambda/thetat loop
137 scalar lambdaErr_;
138
139 //- Maximum number of iterations to converge the lambda/thetat loop
140 label maxLambdaIter_;
141
142 //- Stabilization for division by the magnitude of the velocity
144
145
146 // Fields
147
148 //- Transition onset momentum-thickness Reynolds number
150
151 //- Intermittency
153
154 //- Effective intermittency
156
157
158 // Protected Member Functions
159
160 //- Modified form of the k-omega SST F1 function
161 virtual tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
162
163 //- Modified form of the k-omega SST k production rate
165 (
167 ) const;
168
169 //- Modified form of the k-omega SST epsilon/k
171 (
172 const volScalarField& F1,
173 const volTensorField& gradU
174 ) const;
175
176 //- Freestream blending-function
178 (
180 const volScalarField::Internal& Omega,
182 ) const;
183
184 //- Empirical correlation for critical Reynolds number where the
185 // intermittency first starts to increase in the boundary layer
187
188 //- Empirical correlation that controls the length of the
189 // transition region
191 (
193 ) const;
194
195 //- Transition onset location control function
197 (
198 const volScalarField::Internal& Rev,
201 ) const;
202
203 //- Return the transition onset momentum-thickness Reynolds number
204 // (based on freestream conditions)
206 (
208 const volScalarField::Internal& dUsds,
210 ) const;
211
212 //- Solve the turbulence equations and correct the turbulence viscosity
214
215
216public:
218 typedef typename BasicTurbulenceModel::alphaField alphaField;
219 typedef typename BasicTurbulenceModel::rhoField rhoField;
220 typedef typename BasicTurbulenceModel::transportModel transportModel;
221
222
223 //- Runtime type information
224 TypeName("kOmegaSSTLM");
225
226
227 // Constructors
228
229 //- Construct from components
231 (
232 const alphaField& alpha,
233 const rhoField& rho,
234 const volVectorField& U,
235 const surfaceScalarField& alphaRhoPhi,
236 const surfaceScalarField& phi,
237 const transportModel& transport,
238 const word& propertiesName = turbulenceModel::propertiesName,
239 const word& type = typeName
240 );
241
242
243 //- Destructor
244 virtual ~kOmegaSSTLM() = default;
245
246
247 // Member Functions
248
249 //- Re-read model coefficients if they have changed
250 virtual bool read();
251
252 //- Access function transition onset momentum-thickness Reynolds number
253 const volScalarField& ReThetat() const
254 {
255 return ReThetat_;
256 }
257
258 //- Access function to intermittency
259 const volScalarField& gammaInt() const
260 {
261 return gammaInt_;
262 }
263
264 //- Return the effective diffusivity for transition onset
265 // momentum-thickness Reynolds number
267 {
269 (
271 (
272 "DReThetatEff",
273 sigmaThetat_*(this->nut_ + this->nu())
274 )
275 );
276 }
277
278 //- Return the effective diffusivity for intermittency
280 {
282 (
284 (
285 "DgammaIntEff",
286 this->nut_ + this->nu()
287 )
288 );
289 }
290
291 //- Solve the turbulence equations and correct the turbulence viscosity
292 virtual void correct();
293};
294
295
296// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
297
298} // End namespace RASModels
299} // End namespace Foam
300
301// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
302
303#ifdef NoRepository
304 #include "kOmegaSSTLM.C"
305#endif
306
307// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308
309#endif
310
311// ************************************************************************* //
#define F1(B, C, D)
Definition: SHA1.C:150
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Langtry-Menter 4-equation transitional SST model based on the k-omega-SST RAS model.
Definition: kOmegaSSTLM.H:110
label maxLambdaIter_
Maximum number of iterations to converge the lambda/thetat loop.
Definition: kOmegaSSTLM.H:139
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTLM.H:217
virtual ~kOmegaSSTLM()=default
Destructor.
tmp< volScalarField > DReThetatEff() const
Return the effective diffusivity for transition onset.
Definition: kOmegaSSTLM.H:265
volScalarField gammaInt_
Intermittency.
Definition: kOmegaSSTLM.H:151
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTLM.H:218
tmp< volScalarField::Internal > Flength(const volScalarField::Internal &nu) const
Empirical correlation that controls the length of the.
Definition: kOmegaSSTLM.C:160
volScalarField::Internal gammaIntEff_
Effective intermittency.
Definition: kOmegaSSTLM.H:154
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTLM.C:621
tmp< volScalarField::Internal > Fthetat(const volScalarField::Internal &Us, const volScalarField::Internal &Omega, const volScalarField::Internal &nu) const
Freestream blending-function.
Definition: kOmegaSSTLM.C:80
void correctReThetatGammaInt()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTLM.C:523
const dimensionedScalar deltaU_
Stabilization for division by the magnitude of the velocity.
Definition: kOmegaSSTLM.H:142
dimensionedScalar ce1_
Definition: kOmegaSSTLM.H:129
const volScalarField & gammaInt() const
Access function to intermittency.
Definition: kOmegaSSTLM.H:258
tmp< volScalarField > DgammaIntEff() const
Return the effective diffusivity for intermittency.
Definition: kOmegaSSTLM.H:278
dimensionedScalar ca2_
Definition: kOmegaSSTLM.H:127
tmp< volScalarField::Internal > ReThetac() const
Empirical correlation for critical Reynolds number where the.
Definition: kOmegaSSTLM.C:119
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTLM.H:219
volScalarField ReThetat_
Transition onset momentum-thickness Reynolds number.
Definition: kOmegaSSTLM.H:148
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Modified form of the k-omega SST epsilon/k.
Definition: kOmegaSSTLM.C:67
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Modified form of the k-omega SST k production rate.
Definition: kOmegaSSTLM.C:57
scalar lambdaErr_
Convergence criterion for the lambda/thetat loop.
Definition: kOmegaSSTLM.H:136
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:337
dimensionedScalar cThetat_
Definition: kOmegaSSTLM.H:132
dimensionedScalar sigmaThetat_
Definition: kOmegaSSTLM.H:133
dimensionedScalar ca1_
Definition: kOmegaSSTLM.H:126
dimensionedScalar ce2_
Definition: kOmegaSSTLM.H:130
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:223
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTLM.C:502
TypeName("kOmegaSSTLM")
Runtime type information.
const volScalarField & ReThetat() const
Access function transition onset momentum-thickness Reynolds number.
Definition: kOmegaSSTLM.H:252
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows.
Definition: kOmegaSST.H:132
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
volScalarField & nu
volScalarField & alpha
Us
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73