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
SpalartAllmaras
SpalartAllmaras.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::SpalartAllmaras
26
27
Description
28
SpalartAllmaras for compressible flows
29
30
SourceFiles
31
SpalartAllmaras.C
32
33
\*---------------------------------------------------------------------------*/
34
35
#ifndef compressibleSpalartAllmaras_H
36
#define compressibleSpalartAllmaras_H
37
38
#include <
compressibleLESModels/LESModel.H
>
39
#include <
finiteVolume/volFields.H
>
40
41
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
43
namespace
Foam
44
{
45
namespace
compressible
46
{
47
namespace
LESModels
48
{
49
50
/*---------------------------------------------------------------------------*\
51
Class SpalartAllmaras Declaration
52
\*---------------------------------------------------------------------------*/
53
54
class
SpalartAllmaras
55
:
56
public
LESModel
57
{
58
// Private data
59
60
// Model coefficients
61
62
dimensionedScalar
sigmaNut_;
63
dimensionedScalar
Prt_;
64
65
dimensionedScalar
Cb1_;
66
dimensionedScalar
Cb2_;
67
dimensionedScalar
Cv1_;
68
dimensionedScalar
Cv2_;
69
dimensionedScalar
CDES_;
70
dimensionedScalar
ck_;
71
dimensionedScalar
kappa_;
72
dimensionedScalar
Cw1_;
73
dimensionedScalar
Cw2_;
74
dimensionedScalar
Cw3_;
75
76
77
// Fields
78
79
volScalarField
nuTilda_;
80
volScalarField
dTilda_;
81
volScalarField
muSgs_;
82
volScalarField
alphaSgs_;
83
84
85
// Private member functions
86
87
//- Update sub-grid scale fields
88
void
updateSubGridScaleFields();
89
90
tmp<volScalarField>
fv1()
const
;
91
tmp<volScalarField>
fv2()
const
;
92
tmp<volScalarField>
fv3()
const
;
93
tmp<volScalarField>
fw(
const
volScalarField
& Stilda)
const
;
94
95
// Disallow default bitwise copy construct and assignment
96
SpalartAllmaras
(
const
SpalartAllmaras
&);
97
SpalartAllmaras
& operator=(
const
SpalartAllmaras
&);
98
99
100
public
:
101
102
//- Runtime type information
103
TypeName
(
"SpalartAllmaras"
);
104
105
106
// Constructors
107
108
//- Constructor from components
109
SpalartAllmaras
110
(
111
const
volScalarField
&
rho
,
112
const
volVectorField
&
U
,
113
const
surfaceScalarField
&
phi
,
114
const
basicThermo
& thermoPhysicalModel
115
);
116
117
118
//- Destructor
119
virtual
~
SpalartAllmaras
()
120
{}
121
122
123
// Member Functions
124
125
tmp<volScalarField>
nuTilda
()
const
126
{
127
return
nuTilda_;
128
}
129
130
//- Return SGS kinetic energy
131
virtual
tmp<volScalarField>
k
()
const
132
{
133
return
sqr
(
muSgs
()/
rho
()/ck_/dTilda_);
134
}
135
136
//- Return sub-grid disipation rate
137
virtual
tmp<volScalarField>
epsilon()
const
;
138
139
//- Return SGS viscosity
140
virtual
tmp<volScalarField>
muSgs
()
const
141
{
142
return
muSgs_;
143
}
144
145
//- Return SGS thermal diffusivity
146
virtual
tmp<volScalarField>
alphaSgs
()
const
147
{
148
return
alphaSgs_;
149
}
150
151
//- Return thermal conductivity
152
virtual
tmp<volScalarField>
alphaEff
()
const
153
{
154
return
tmp<volScalarField>
155
(
156
new
volScalarField
(
"alphaEff"
, alphaSgs_ +
alpha
())
157
);
158
}
159
160
//- Return the sub-grid stress tensor.
161
virtual
tmp<volSymmTensorField>
B()
const
;
162
163
//- Return the deviatoric part of the effective sub-grid
164
// turbulence stress tensor including the laminar stress
165
virtual
tmp<volSymmTensorField>
devRhoBeff()
const
;
166
167
//- Returns div(rho*dev(B)).
168
// This is the additional term due to the filtering of the NSE.
169
virtual
tmp<fvVectorMatrix>
divDevRhoBeff(
volVectorField
& U)
const
;
170
171
//- Correct nuTilda and related properties
172
virtual
void
correct
(
const
tmp<volTensorField>
& gradU);
173
174
//- Read LESProperties dictionary
175
virtual
bool
read();
176
};
177
178
179
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180
181
}
// End namespace LESModels
182
}
// End namespace compressible
183
}
// End namespace Foam
184
185
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186
187
#endif
188
189
// ************************ vim: set sw=4 sts=4 et: ************************ //