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
incompressible
LES
Smagorinsky
Smagorinsky.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::incompressible::LESModels::Smagorinsky
26
27
Description
28
The Isochoric Smagorinsky Model for incompressible flows.
29
30
Algebraic eddy viscosity SGS model founded on the assumption that
31
local equilibrium prevails.
32
Thus,
33
@verbatim
34
B = 2/3*k*I - 2*nuSgs*dev(D)
35
Beff = 2/3*k*I - 2*nuEff*dev(D)
36
37
where
38
39
D = symm(grad(U));
40
k = (2*ck/ce)*delta^2*||D||^2
41
nuSgs = ck*sqrt(k)*delta
42
nuEff = nuSgs + nu
43
@endverbatim
44
45
SourceFiles
46
Smagorinsky.C
47
48
\*---------------------------------------------------------------------------*/
49
50
#ifndef Smagorinsky_H
51
#define Smagorinsky_H
52
53
#include <
incompressibleLESModels/GenEddyVisc.H
>
54
55
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56
57
namespace
Foam
58
{
59
namespace
incompressible
60
{
61
namespace
LESModels
62
{
63
64
/*---------------------------------------------------------------------------*\
65
Class Smagorinsky Declaration
66
\*---------------------------------------------------------------------------*/
67
68
class
Smagorinsky
69
:
70
public
GenEddyVisc
71
{
72
// Private data
73
74
dimensionedScalar
ck_;
75
76
77
// Private Member Functions
78
79
//- Update sub-grid scale fields
80
void
updateSubGridScaleFields(
const
volTensorField
& gradU);
81
82
// Disallow default bitwise copy construct and assignment
83
Smagorinsky
(
const
Smagorinsky
&);
84
Smagorinsky
& operator=(
const
Smagorinsky
&);
85
86
87
public
:
88
89
//- Runtime type information
90
TypeName
(
"Smagorinsky"
);
91
92
// Constructors
93
94
//- Construct from components
95
Smagorinsky
96
(
97
const
volVectorField
&
U
,
98
const
surfaceScalarField
&
phi
,
99
transportModel
&
transport
100
);
101
102
103
//- Destructor
104
virtual
~Smagorinsky
()
105
{}
106
107
108
// Member Functions
109
110
//- Return SGS kinetic energy
111
// calculated from the given velocity gradient
112
tmp<volScalarField>
k
(
const
tmp<volTensorField>
& gradU)
const
113
{
114
return
(2.0*ck_/
ce_
)*
sqr
(
delta
())*
magSqr
(
dev
(
symm
(gradU)));
115
}
116
117
//- Return SGS kinetic energy
118
virtual
tmp<volScalarField>
k
()
const
119
{
120
return
k
(
fvc::grad
(
U
()));
121
}
122
123
124
//- Correct Eddy-Viscosity and related properties
125
virtual
void
correct
(
const
tmp<volTensorField>
& gradU);
126
127
//- Read LESProperties dictionary
128
virtual
bool
read
();
129
};
130
131
132
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133
134
}
// End namespace LESModels
135
}
// End namespace incompressible
136
}
// End namespace Foam
137
138
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139
140
#endif
141
142
// ************************ vim: set sw=4 sts=4 et: ************************ //