-
Notifications
You must be signed in to change notification settings - Fork 1
/
kOmegaML.H
221 lines (170 loc) · 5.84 KB
/
kOmegaML.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
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 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::RASModels::kOmegaML
Group
grpRASTurbulence
Description
Standard high Reynolds-number k-omega turbulence model for
incompressible and compressible flows.
References:
\verbatim
Wilcox, D. C. (1998).
Turbulence modeling for CFD
(Vol. 2, pp. 103-217). La Canada, CA: DCW industries.
\endverbatim
The default model coefficients are
\verbatim
{
Cmu 0.09; // Equivalent to betaStar (page87)
alpha 0.52;
beta 0.072;
alphak 0.5; // Equivalent to sigmaStar (page87)
alphaOmega 0.5; // Equivalent to sigma (page87)
}
\endverbatim
SourceFiles
kOmegaML.C
\*---------------------------------------------------------------------------*/
#ifndef kOmegaML_H
#define kOmegaML_H
#include "RASModel.H"
#include "eddyViscosity.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace RASModels
{
/*---------------------------------------------------------------------------*\
Class kOmegaML Declaration
\*---------------------------------------------------------------------------*/
template<class BasicTurbulenceModel>
class kOmegaML
:
public eddyViscosity<RASModel<BasicTurbulenceModel>>
{
protected:
// Protected data
// Model coefficients
dimensionedScalar Cmu_;
dimensionedScalar beta_;
dimensionedScalar gamma_;
dimensionedScalar alphaK_;
dimensionedScalar alphaOmega_;
// Fields
volScalarField k_;
volScalarField omega_;
volScalarField MLfCorr_;
// Protected Member Functions
virtual void correctNut();
public:
typedef typename BasicTurbulenceModel::alphaField alphaField;
typedef typename BasicTurbulenceModel::rhoField rhoField;
typedef typename BasicTurbulenceModel::transportModel transportModel;
//- Runtime type information
TypeName("kOmegaML");
// Constructors
//- Construct from components
kOmegaML
(
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaRhoPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName = turbulenceModel::propertiesName,
const word& type = typeName
);
//- Destructor
virtual ~kOmegaML()
{}
// Member Functions
//- Read RASProperties dictionary
virtual bool read();
const volScalarField& MLfCorr() const
{
return MLfCorr_;
}
//- Return the effective diffusivity for k
tmp<volScalarField> DkEff() const
{
return tmp<volScalarField>
(
new volScalarField
(
"DkEff",
alphaK_*this->nut_ + this->nu()
)
);
}
//- Return the effective diffusivity for omega
tmp<volScalarField> DomegaEff() const
{
return tmp<volScalarField>
(
new volScalarField
(
"DomegaEff",
alphaOmega_*this->nut_ + this->nu()
)
);
}
//- Return the turbulence kinetic energy
virtual tmp<volScalarField> k() const
{
return k_;
}
//- Return the turbulence specific dissipation rate
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",
this->mesh_.time().timeName(),
this->mesh_
),
Cmu_*k_*omega_,
omega_.boundaryField().types()
)
);
}
//- Solve the turbulence equations and correct the turbulence viscosity
virtual void correct();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
// #include "kOmegaML.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //