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
GenEddyVisc
GenEddyVisc.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::GenEddyVisc
26
27
Description
28
General base class for all compressible models that can be implemented as
29
an eddy viscosity, i.e. algebraic and one-equation models.
30
31
Contains fields for k (SGS turbulent kinetic energy), gamma
32
(modelled viscosity) and epsilon (SGS dissipation).
33
34
SourceFiles
35
GenEddyVisc.C
36
37
\*---------------------------------------------------------------------------*/
38
39
#ifndef compressibleGenEddyVisc_H
40
#define compressibleGenEddyVisc_H
41
42
#include <
compressibleLESModels/LESModel.H
>
43
44
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46
namespace
Foam
47
{
48
namespace
compressible
49
{
50
namespace
LESModels
51
{
52
53
/*---------------------------------------------------------------------------*\
54
Class GenEddyVisc Declaration
55
\*---------------------------------------------------------------------------*/
56
57
class
GenEddyVisc
58
:
59
virtual
public
LESModel
60
{
61
// Private Member Functions
62
63
// Disallow default bitwise copy construct and assignment
64
GenEddyVisc
(
const
GenEddyVisc
&);
65
GenEddyVisc
& operator=(
const
GenEddyVisc
&);
66
67
68
protected
:
69
70
// Model coefficients
71
72
dimensionedScalar
ce_
;
73
dimensionedScalar
Prt_
;
74
75
76
// Fields
77
78
volScalarField
k_
;
79
volScalarField
muSgs_
;
80
volScalarField
alphaSgs_
;
81
82
83
public
:
84
85
// Constructors
86
87
//- Construct from components
88
GenEddyVisc
89
(
90
const
volScalarField
&
rho
,
91
const
volVectorField
&
U
,
92
const
surfaceScalarField
&
phi
,
93
const
basicThermo
& thermoPhysicalModel
94
);
95
96
97
//- Destructor
98
virtual
~
GenEddyVisc
()
99
{}
100
101
102
// Member Functions
103
104
//- Return SGS kinetic energy
105
virtual
tmp<volScalarField>
k
()
const
106
{
107
return
k_
;
108
}
109
110
//- Return sub-grid disipation rate
111
virtual
tmp<volScalarField>
epsilon
()
const
112
{
113
return
ce_
*
k_
*
sqrt
(
k_
)/
delta
();
114
}
115
116
//- Return viscosity
117
virtual
tmp<volScalarField>
muSgs
()
const
118
{
119
return
muSgs_
;
120
}
121
122
//- Return thermal diffusivity
123
virtual
tmp<volScalarField>
alphaSgs
()
const
124
{
125
return
alphaSgs_
;
126
}
127
128
//- Return thermal diffusivity
129
virtual
tmp<volScalarField>
alphaEff
()
const
130
{
131
return
tmp<volScalarField>
132
(
133
new
volScalarField
(
"alphaEff"
,
alphaSgs_
+
alpha
())
134
);
135
}
136
137
//- Return the sub-grid stress tensor.
138
virtual
tmp<volSymmTensorField>
B()
const
;
139
140
//- Return the deviatoric part of the effective sub-grid
141
// turbulence stress tensor including the laminar stress
142
virtual
tmp<volSymmTensorField>
devRhoBeff()
const
;
143
144
//- Returns div(rho*dev(B)).
145
// This is the additional term due to the filtering of the NSE.
146
virtual
tmp<fvVectorMatrix>
divDevRhoBeff(
volVectorField
&
U
)
const
;
147
148
//- Correct Eddy-Viscosity and related properties
149
virtual
void
correct
(
const
tmp<volTensorField>
& gradU);
150
151
//- Read LESProperties dictionary
152
virtual
bool
read();
153
};
154
155
156
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157
158
}
// End namespace LESModels
159
}
// End namespace compressible
160
}
// End namespace Foam
161
162
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163
164
#endif
165
166
// ************************ vim: set sw=4 sts=4 et: ************************ //