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