dynamicKEqn.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "dynamicKEqn.H"
30#include "fvOptions.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace LESModels
37{
38
39// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40
41template<class BasicTurbulenceModel>
43(
44 const volSymmTensorField& D,
45 const volScalarField& KK
46) const
47{
48 const volSymmTensorField LL
49 (
50 simpleFilter_(dev(filter_(sqr(this->U_)) - (sqr(filter_(this->U_)))))
51 );
52
53 const volSymmTensorField MM
54 (
55 simpleFilter_
56 (
57 -2.0*this->delta()*sqrt
58 (
60 )*filter_(D)
61 )
62 );
63
64 const volScalarField Ck
65 (
66 simpleFilter_(0.5*(LL && MM))
67 /(
68 simpleFilter_(magSqr(MM))
69 + dimensionedScalar("small", sqr(MM.dimensions()), VSMALL)
70 )
71 );
72
73 tmp<volScalarField> tfld = 0.5*(mag(Ck) + Ck);
74 return tfld();
75}
76
77
78template<class BasicTurbulenceModel>
80(
81 const volSymmTensorField& D,
82 const volScalarField& KK
83) const
84{
85 const volScalarField Ce
86 (
87 simpleFilter_(this->nuEff()*(filter_(magSqr(D)) - magSqr(filter_(D))))
88 /simpleFilter_(pow(KK, 1.5)/(2.0*this->delta()))
89 );
90
91 tmp<volScalarField> tfld = 0.5*(mag(Ce) + Ce);
92 return tfld();
93}
94
95
96template<class BasicTurbulenceModel>
98{
99 const volSymmTensorField D(dev(symm(fvc::grad(this->U_))));
100
102 (
103 0.5*(filter_(magSqr(this->U_)) - magSqr(filter_(this->U_)))
104 );
105 KK.max(dimensionedScalar("small", KK.dimensions(), SMALL));
106
107 return Ce(D, KK);
108}
109
110
111template<class BasicTurbulenceModel>
113(
114 const volSymmTensorField& D,
115 const volScalarField& KK
116)
117{
118 this->nut_ = Ck(D, KK)*sqrt(k_)*this->delta();
119 this->nut_.correctBoundaryConditions();
120 fv::options::New(this->mesh_).correct(this->nut_);
121
122 BasicTurbulenceModel::correctNut();
123}
124
125
126template<class BasicTurbulenceModel>
128{
129 const volScalarField KK
130 (
131 0.5*(filter_(magSqr(this->U_)) - magSqr(filter_(this->U_)))
132 );
133
134 correctNut(symm(fvc::grad(this->U_)), KK);
135}
136
137
138template<class BasicTurbulenceModel>
140{
142 (
144 (
145 k_,
146 dimVolume*this->rho_.dimensions()*k_.dimensions()
147 /dimTime
148 )
149 );
150}
151
152
153// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
154
155template<class BasicTurbulenceModel>
157(
158 const alphaField& alpha,
159 const rhoField& rho,
160 const volVectorField& U,
161 const surfaceScalarField& alphaRhoPhi,
162 const surfaceScalarField& phi,
163 const transportModel& transport,
164 const word& propertiesName,
165 const word& type
166)
167:
168 LESeddyViscosity<BasicTurbulenceModel>
169 (
170 type,
171 alpha,
172 rho,
173 U,
174 alphaRhoPhi,
175 phi,
176 transport,
177 propertiesName
178 ),
179
180 k_
181 (
183 (
184 IOobject::groupName("k", this->alphaRhoPhi_.group()),
185 this->runTime_.timeName(),
186 this->mesh_,
187 IOobject::MUST_READ,
188 IOobject::AUTO_WRITE
189 ),
190 this->mesh_
191 ),
192
193 simpleFilter_(this->mesh_),
194 filterPtr_(LESfilter::New(this->mesh_, this->coeffDict())),
195 filter_(filterPtr_())
196{
197 bound(k_, this->kMin_);
198
199 if (type == typeName)
200 {
201 this->printCoeffs(type);
202 }
203}
204
205
206// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
207
208template<class BasicTurbulenceModel>
210{
212 {
213 filter_.read(this->coeffDict());
214
215 return true;
216 }
217
218 return false;
219}
220
221
222template<class BasicTurbulenceModel>
224{
225 if (!this->turbulence_)
226 {
227 return;
228 }
229
230 // Local references
231 const alphaField& alpha = this->alpha_;
232 const rhoField& rho = this->rho_;
233 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
234 const volVectorField& U = this->U_;
235 volScalarField& nut = this->nut_;
237
239
241
243 const volSymmTensorField D(dev(symm(tgradU())));
244 const volScalarField G(this->GName(), 2.0*nut*(tgradU() && D));
245 tgradU.clear();
246
247 volScalarField KK(0.5*(filter_(magSqr(U)) - magSqr(filter_(U))));
248 KK.max(dimensionedScalar("small", KK.dimensions(), SMALL));
249
251 (
252 fvm::ddt(alpha, rho, k_)
253 + fvm::div(alphaRhoPhi, k_)
254 - fvm::laplacian(alpha*rho*DkEff(), k_)
255 ==
256 alpha*rho*G
257 - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
258 - fvm::Sp(Ce(D, KK)*alpha*rho*sqrt(k_)/this->delta(), k_)
259 + kSource()
260 + fvOptions(alpha, rho, k_)
261 );
262
263 kEqn.ref().relax();
264 fvOptions.constrain(kEqn.ref());
265 solve(kEqn);
266 fvOptions.correct(k_);
267 bound(k_, this->kMin_);
268
269 correctNut(D, KK);
270}
271
272
273// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274
275} // End namespace LESModels
276} // End namespace Foam
277
278// ************************************************************************* //
scalar delta
fv::options & fvOptions
surfaceScalarField & phi
const dimensionSet & dimensions() const
Return dimensions.
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Eddy viscosity LES SGS model base class.
Dynamic one equation eddy-viscosity model.
Definition: dynamicKEqn.H:82
BasicTurbulenceModel::alphaField alphaField
Definition: dynamicKEqn.H:140
volScalarField Ce() const
Definition: dynamicKEqn.C:97
BasicTurbulenceModel::rhoField rhoField
Definition: dynamicKEqn.H:141
volScalarField Ck(const volSymmTensorField &D, const volScalarField &KK) const
Calculate Ck by filtering the velocity field U.
Definition: dynamicKEqn.C:43
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicKEqn.C:223
virtual tmp< fvScalarMatrix > kSource() const
Definition: dynamicKEqn.C:139
BasicTurbulenceModel::transportModel transportModel
Definition: dynamicKEqn.H:142
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicKEqn.C:209
One equation eddy-viscosity model.
Definition: kEqn.H:80
Abstract class for LES filters.
Definition: LESfilter.H:58
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Finite-volume options.
Definition: fvOptions.H:59
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
scalar nut
zeroField divU
Definition: alphaSuSp.H:3
word timeName
Definition: getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
volScalarField & alpha
const dimensionedScalar & D
CEqn solve()