/
kOmegaSST.H
305 lines (234 loc) · 8.54 KB
/
kOmegaSST.H
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::compressible::RASModels::kOmegaSST
Description
Implementation of the k-omega-SST turbulence model for compressible flows.
Turbulence model described in:
\verbatim
Menter, F., Esch, T.
"Elements of Industrial Heat Transfer Prediction"
16th Brazilian Congress of Mechanical Engineering (COBEM),
Nov. 2001
\endverbatim
with the addition of the optional F3 term for rough walls from
\verbatim
Hellsten, A.
"Some Improvements in Menter’s k-omega-SST turbulence model"
29th AIAA Fluid Dynamics Conference,
AIAA-98-2554,
June 1998.
\endverbatim
Note that this implementation is written in terms of alpha diffusion
coefficients rather than the more traditional sigma (alpha = 1/sigma) so
that the blending can be applied to all coefficuients in a consistent
manner. The paper suggests that sigma is blended but this would not be
consistent with the blending of the k-epsilon and k-omega models.
Also note that the error in the last term of equation (2) relating to
sigma has been corrected.
Wall-functions are applied in this implementation by using equations (14)
to specify the near-wall omega as appropriate.
The blending functions (15) and (16) are not currently used because of the
uncertainty in their origin, range of applicability and that is y+ becomes
sufficiently small blending u_tau in this manner clearly becomes nonsense.
The default model coefficients correspond to the following:
\verbatim
kOmegaSSTCoeffs
{
alphaK1 0.85034;
alphaK2 1.0;
alphaOmega1 0.5;
alphaOmega2 0.85616;
Prt 1.0; // only for compressible
beta1 0.075;
beta2 0.0828;
betaStar 0.09;
gamma1 0.5532;
gamma2 0.4403;
a1 0.31;
b1 1.0;
c1 10.0;
F3 no;
}
\endverbatim
SourceFiles
kOmegaSST.C
\*---------------------------------------------------------------------------*/
#ifndef compressiblekOmegaSST_H
#define compressiblekOmegaSST_H
#include "RASModel.H"
#include "wallDist.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace compressible
{
namespace RASModels
{
/*---------------------------------------------------------------------------*\
Class kOmegaSST Declaration
\*---------------------------------------------------------------------------*/
class kOmegaSST
:
public RASModel
{
protected:
// Protected data
// Model coefficients
dimensionedScalar alphaK1_;
dimensionedScalar alphaK2_;
dimensionedScalar alphaOmega1_;
dimensionedScalar alphaOmega2_;
dimensionedScalar Prt_;
dimensionedScalar gamma1_;
dimensionedScalar gamma2_;
dimensionedScalar beta1_;
dimensionedScalar beta2_;
dimensionedScalar betaStar_;
dimensionedScalar a1_;
dimensionedScalar b1_;
dimensionedScalar c1_;
Switch F3_;
//- Wall distance
// Note: different to wall distance in parent RASModel
wallDist y_;
// Fields
volScalarField k_;
volScalarField omega_;
volScalarField mut_;
volScalarField alphat_;
// Private Member Functions
tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
tmp<volScalarField> F2() const;
tmp<volScalarField> F3() const;
tmp<volScalarField> F23() const;
tmp<volScalarField> blend
(
const volScalarField& F1,
const dimensionedScalar& psi1,
const dimensionedScalar& psi2
) const
{
return F1*(psi1 - psi2) + psi2;
}
tmp<volScalarField> alphaK(const volScalarField& F1) const
{
return blend(F1, alphaK1_, alphaK2_);
}
tmp<volScalarField> alphaOmega(const volScalarField& F1) const
{
return blend(F1, alphaOmega1_, alphaOmega2_);
}
tmp<volScalarField> beta(const volScalarField& F1) const
{
return blend(F1, beta1_, beta2_);
}
tmp<volScalarField> gamma(const volScalarField& F1) const
{
return blend(F1, gamma1_, gamma2_);
}
public:
//- Runtime type information
TypeName("kOmegaSST");
// Constructors
//- Construct from components
kOmegaSST
(
const volScalarField& rho,
const volVectorField& U,
const surfaceScalarField& phi,
const basicThermo& thermophysicalModel,
const word& turbulenceModelName = turbulenceModel::typeName,
const word& modelName = typeName
);
//- Destructor
virtual ~kOmegaSST()
{}
// Member Functions
//- Return the effective diffusivity for k
tmp<volScalarField> DkEff(const volScalarField& F1) const
{
return tmp<volScalarField>
(
new volScalarField("DkEff", alphaK(F1)*mut_ + mu())
);
}
//- Return the effective diffusivity for omega
tmp<volScalarField> DomegaEff(const volScalarField& F1) const
{
return tmp<volScalarField>
(
new volScalarField("DomegaEff", alphaOmega(F1)*mut_ + mu())
);
}
virtual tmp<volScalarField> mut() const
{
return mut_;
}
//- Return the turbulence thermal diffusivity
virtual tmp<volScalarField> alphat() const
{
return alphat_;
}
//- Return the turbulence kinetic energy
virtual tmp<volScalarField> k() const
{
return k_;
}
virtual tmp<volScalarField> omega() const
{
return omega_;
}
//- Return the turbulence kinetic energy dissipation rate
virtual tmp<volScalarField> epsilon() const
{
return tmp<volScalarField>
(
new volScalarField
(
IOobject
(
"epsilon",
mesh_.time().timeName(),
mesh_
),
betaStar_*k_*omega_,
omega_.boundaryField().types()
)
);
}
//- Return the Reynolds stress tensor
virtual tmp<volSymmTensorField> R() const;
//- Return the effective stress tensor including the laminar stress
virtual tmp<volSymmTensorField> devRhoReff() const;
//- Return the source term for the momentum equation
virtual tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U) const;
//- Solve the turbulence equations and correct the turbulence viscosity
virtual void correct();
//- Read RASProperties dictionary
virtual bool read();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace compressible
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //