-
Notifications
You must be signed in to change notification settings - Fork 60
/
fluidDynamicLagrangian.H
195 lines (140 loc) · 5.41 KB
/
fluidDynamicLagrangian.H
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
/*---------------------------------------------------------------------------*\
Copyright (C) 2015 Cyrille Bonamy, Julien Chauchat, Tian-Jian Hsu
and contributors
License
This file is part of SedFOAM.
SedFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
SedFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with SedFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::LESModels::fluidDynamicLagrangian
Group
grpLESTurbulence
Description
Dynamic SGS model with Lagrangian averaging for the fluid phase
Reference:
\verbatim
Meneveau, C., Lund, T. S., & Cabot, W. H. (1996).
A Lagrangian dynamic subgrid-scale model of turbulence.
Journal of Fluid Mechanics, 319, 353-385.
\endverbatim
SourceFiles
fluidDynamicLagrangian.C
\*---------------------------------------------------------------------------*/
#ifndef fluidDynamicLagrangian_H
#define fluidDynamicLagrangian_H
#include "LESeddyViscosity.H"
#include "simpleFilter.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace LESModels
{
/*---------------------------------------------------------------------------*\
Class fluidDynamicLagrangian Declaration
\*---------------------------------------------------------------------------*/
template<class BasicTurbulenceModel>
class fluidDynamicLagrangian
:
public LESeddyViscosity<BasicTurbulenceModel>
{
// Private Member Functions
// Disallow default bitwise copy construct and assignment
fluidDynamicLagrangian(const fluidDynamicLagrangian&);
void operator=(const fluidDynamicLagrangian&);
protected:
// Protected data
volScalarField flmb_;
volScalarField fmmb_;
volScalarField Cyb_;
dimensionedScalar theta_;
dimensionedScalar D_;
dimensionedScalar nutMax_;
dimensionedScalar CybMax_;
simpleFilter simpleFilter_;
autoPtr<LESfilter> filterPtr_;
LESfilter& filter_;
dimensionedScalar flmb0_;
dimensionedScalar fmmb0_;
// Protected Member Functions
//- Update sub-grid eddy-viscosity
void correctNut(const tmp<volTensorField>& gradU);
volScalarField fDelta
(
const tmp<volScalarField>& delta,
const tmp<volScalarField>& alpha,
const tmp<volVectorField>& Ua,
const tmp<volVectorField>& Ub,
const tmp<volScalarField>& rho,
const tmp<volScalarField>& K
);
volScalarField hAlpha(const tmp<volScalarField>& alpha);
virtual void correctNut();
public:
typedef typename BasicTurbulenceModel::alphaField alphaField;
typedef typename BasicTurbulenceModel::rhoField rhoField;
typedef typename BasicTurbulenceModel::transportModel transportModel;
//- Runtime type information
TypeName("fluidDynamicLagrangian");
// Constructors
//- Construct from components
fluidDynamicLagrangian
(
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaRhoPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName = turbulenceModel::propertiesName,
const word& type = typeName
);
//- Destructor
virtual ~fluidDynamicLagrangian()
{}
// Member Functions
//- Read model coefficients if they have changed
virtual bool read();
//- Return SGS kinetic energy of phase a
tmp<volScalarField> k(const tmp<volTensorField>& gradU) const
{
return pow(2.0*flmb_/fmmb_, 2.0/3.0)
* pow(this->Ce_, -2.0/3.0)
* sqr(this->delta())*magSqr(dev(symm(gradU)));
}
//- Return SGS kinetic energy
virtual tmp<volScalarField> k() const
{
return k(fvc::grad(this->U_));
}
//- Return the effective diffusivity for k
tmp<volScalarField> DkEff() const
{
return tmp<volScalarField>
(
new volScalarField("DkEff", this->nut_ + this->nu())
);
}
//- Correct Eddy-Viscosity and related properties
virtual void correct();
//- Return spherical part of the SGS stress tensor
virtual tmp<volScalarField> spherSigmaSGS();
virtual tmp<volScalarField> sqrDeltaD();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace LESModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "fluidDynamicLagrangian.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //