Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
src
turbulenceModels
compressible
LES
lowReOneEqEddy
lowReOneEqEddy.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 | Copyright (C) 1991-2010 OpenCFD Ltd.
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::compressible::LESModels::lowReOneEqEddy
26
27
Description
28
One Equation Eddy Viscosity Model for compressible flow
29
30
@verbatim
31
d/dt(rho*k) + div(rho*U*k) - div(muEff*grad(k))
32
=
33
-rho*B*L - ce*rho*k^3/2/delta
34
35
and
36
37
B = 2/3*k*I - 2*nuSgs*dev(D)
38
39
where
40
41
nuSgsHiRe = ck*sqrt(k)*delta
42
nuSgs = (nu/beta)*(1 - exp(-beta*nuSgsHiRe/nu));
43
@endverbatim
44
45
SourceFiles
46
lowReOneEqEddy.C
47
48
\*---------------------------------------------------------------------------*/
49
50
#ifndef compressibleLowReOneEqEddy_H
51
#define compressibleLowReOneEqEddy_H
52
53
#include <
compressibleLESModels/GenEddyVisc.H
>
54
55
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56
57
namespace
Foam
58
{
59
namespace
compressible
60
{
61
namespace
LESModels
62
{
63
64
/*---------------------------------------------------------------------------*\
65
Class lowReOneEqEddy Declaration
66
\*---------------------------------------------------------------------------*/
67
68
class
lowReOneEqEddy
69
:
70
public
GenEddyVisc
71
{
72
// Private data
73
74
dimensionedScalar
ck_;
75
dimensionedScalar
beta_;
76
77
// Private Member Functions
78
79
//- Update sub-grid scale fields
80
void
updateSubGridScaleFields();
81
82
// Disallow default bitwise copy construct and assignment
83
lowReOneEqEddy
(
const
lowReOneEqEddy
&);
84
lowReOneEqEddy
& operator=(
const
lowReOneEqEddy
&);
85
86
87
public
:
88
89
//- Runtime type information
90
TypeName
(
"lowReOneEqEddy"
);
91
92
93
// Constructors
94
95
//- Constructor from components
96
lowReOneEqEddy
97
(
98
const
volScalarField
&
rho
,
99
const
volVectorField
&
U
,
100
const
surfaceScalarField
&
phi
,
101
const
basicThermo
& thermoPhysicalModel
102
);
103
104
105
//- Destructor
106
virtual
~
lowReOneEqEddy
()
107
{}
108
109
110
// Member Functions
111
112
//- Return the effective diffusivity for k
113
tmp<volScalarField>
DkEff
()
const
114
{
115
return
tmp<volScalarField>
116
(
117
new
volScalarField
(
"DkEff"
,
muSgs_
+
mu
())
118
);
119
}
120
121
//- Correct Eddy-Viscosity and related properties
122
virtual
void
correct
(
const
tmp<volTensorField>
& gradU);
123
124
//- Read LESProperties dictionary
125
virtual
bool
read();
126
};
127
128
129
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130
131
}
// End namespace LESModels
132
}
// End namespace compressible
133
}
// End namespace Foam
134
135
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136
137
#endif
138
139
// ************************ vim: set sw=4 sts=4 et: ************************ //