-
Notifications
You must be signed in to change notification settings - Fork 2
/
dynamicSmagorinsky.C
211 lines (162 loc) · 5.21 KB
/
dynamicSmagorinsky.C
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
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM 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.
OpenFOAM 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 OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "dynamicSmagorinsky.H"
#include "fvOptions.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace LESModels
{
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
template<class BasicTurbulenceModel>
void dynamicSmagorinsky<BasicTurbulenceModel>::correctNut
(
const tmp<volTensorField>& gradU
)
{
this->nut_ = max(Cs_*sqr(this->delta())*mag(dev(symm(gradU))),-1.0*this->nu());
this->nut_.correctBoundaryConditions();
fv::options::New(this->mesh_).correct(this->nut_);
BasicTurbulenceModel::correctNut();
}
template<class BasicTurbulenceModel>
void dynamicSmagorinsky<BasicTurbulenceModel>::correctNut()
{
correctNut(fvc::grad(this->U_));
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
dynamicSmagorinsky<BasicTurbulenceModel>::dynamicSmagorinsky
(
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaRhoPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName,
const word& type
)
:
LESeddyViscosity<BasicTurbulenceModel>
(
type,
alpha,
rho,
U,
alphaRhoPhi,
phi,
transport,
propertiesName
),
k_
(
IOobject
(
IOobject::groupName("k", this->U_.group()),
this->runTime_.timeName(),
this->mesh_,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
this->mesh_
),
Cs_
(
IOobject
(
IOobject::groupName("Cs", this->U_.group()),
this->runTime_.timeName(),
this->mesh_,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
this->mesh_,
dimensionedScalar ("Cs", dimless,SMALL)
),
simpleFilter_(U.mesh()),
filterPtr_(LESfilter::New(U.mesh(), this->coeffDict())),
filter_(filterPtr_())
{
// bound(k_, this->kMin_);
if (type == typeName)
{
this->printCoeffs(type);
}
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
bool dynamicSmagorinsky<BasicTurbulenceModel>::read()
{
if (LESeddyViscosity<BasicTurbulenceModel>::read())
{
filter_.read(this->coeffDict());
return true;
}
else
{
return false;
}
}
template<class BasicTurbulenceModel>
void dynamicSmagorinsky<BasicTurbulenceModel>::correct()
{
if (!this->turbulence_)
{
return;
}
// Local references
const surfaceScalarField& phi = this->phi_;
const volVectorField& U = this->U_;
fv::options& fvOptions(fv::options::New(this->mesh_));
LESeddyViscosity<BasicTurbulenceModel>::correct();
tmp<volTensorField> tgradU(fvc::grad(U));
const volTensorField& gradU = tgradU();
volSymmTensorField S(dev(symm(gradU)));
volScalarField magS(mag(S));
volVectorField Uf(filter_(U));
volSymmTensorField Sf(filter_(S));
// volSymmTensorField Sf(dev(symm(fvc::grad(Uf))));
volScalarField magSf(mag(Sf));
volSymmTensorField LL =
(dev(filter_(sqr(U)) - (sqr(filter_(U)))));
volSymmTensorField MM
(
sqr(this->delta())*(filter_(magS*S) - 4.0*magSf*Sf)
);
volScalarField MMMM = fvc::average(magSqr(MM));
MMMM.max(VSMALL);
Cs_= 0.5* fvc::average(LL && MM)/MMMM;
volScalarField KK =
0.5*(filter_(magSqr(U)) - magSqr(filter_(U)));
volScalarField mm
(
sqr(this->delta())*(4.0*sqr(mag(Sf)) - filter_(sqr(magS)))
);
volScalarField mmmm = fvc::average(magSqr(mm));
mmmm.max(VSMALL);
k_ = fvc::average(KK*mm)/mmmm * sqr(this->delta())*magSqr(S);
correctNut(gradU);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace LESModels
} // End namespace Foam
// ************************************************************************* //