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
DeardorffDiffStress
DeardorffDiffStress.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::DeardorffDiffStress
26
27
Description
28
Differential SGS Stress Equation Model for compressible flows
29
30
The DSEM uses a model version of the full balance equation for the SGS
31
stress tensor to simulate the behaviour of B.
32
Thus,
33
@verbatim
34
d/dt(rho*B) + div(rho*U*B) - div(muSgs*grad(B))
35
=
36
P - c1*rho*epsilon/k*B - 0.667*(1 - c1)*rho*epsilon*I - c2*(P - 0.333*trP*I)
37
38
where
39
40
k = 0.5*trB,
41
epsilon = ce*k^3/2/delta,
42
epsilon/k = ce*k^1/2/delta
43
P = -rho*(B'L + L'B)
44
muSgs = ck*rho*sqrt(k)*delta
45
muEff = muSgs + mu
46
@endverbatim
47
48
SourceFiles
49
DeardorffDiffStress.C
50
51
\*---------------------------------------------------------------------------*/
52
53
#ifndef compressibleDeardorffDiffStress_H
54
#define compressibleDeardorffDiffStress_H
55
56
#include <
compressibleLESModels/GenSGSStress.H
>
57
58
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59
60
namespace
Foam
61
{
62
namespace
compressible
63
{
64
namespace
LESModels
65
{
66
67
/*---------------------------------------------------------------------------*\
68
Class DeardorffDiffStress Declaration
69
\*---------------------------------------------------------------------------*/
70
71
class
DeardorffDiffStress
72
:
73
public
GenSGSStress
74
{
75
// Private data
76
77
dimensionedScalar
ck_;
78
dimensionedScalar
cm_;
79
80
81
// Private Member Functions
82
83
//- Update sub-grid scale fields
84
void
updateSubGridScaleFields(
const
volScalarField
&
K
);
85
86
// Disallow default bitwise copy construct and assignment
87
DeardorffDiffStress
(
const
DeardorffDiffStress
&);
88
DeardorffDiffStress
& operator=(
const
DeardorffDiffStress
&);
89
90
91
public
:
92
93
//- Runtime type information
94
TypeName
(
"DeardorffDiffStress"
);
95
96
// Constructors
97
98
//- Constructor from components
99
DeardorffDiffStress
100
(
101
const
volScalarField
&
rho
,
102
const
volVectorField
&
U
,
103
const
surfaceScalarField
&
phi
,
104
const
basicThermo
& thermoPhysicalModel
105
);
106
107
108
//- Destructor
109
virtual
~
DeardorffDiffStress
()
110
{}
111
112
113
// Member Functions
114
115
//- Return the effective diffusivity for B
116
tmp<volScalarField>
DBEff
()
const
117
{
118
return
tmp<volScalarField>
119
(
120
new
volScalarField
(
"DBEff"
,
muSgs_
+
mu
())
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 compressible
136
}
// End namespace Foam
137
138
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139
140
#endif
141
142
// ************************ vim: set sw=4 sts=4 et: ************************ //