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-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
41 template<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  (
59  max(KK, dimensionedScalar(KK.dimensions(), Zero))
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 
78 template<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 
96 template<class BasicTurbulenceModel>
98 {
99  const volSymmTensorField D(dev(symm(fvc::grad(this->U_))));
100 
101  volScalarField KK
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 
111 template<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 
126 template<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 
138 template<class BasicTurbulenceModel>
140 {
141  return tmp<fvScalarMatrix>
142  (
143  new fvScalarMatrix
144  (
145  k_,
146  dimVolume*this->rho_.dimensions()*k_.dimensions()
147  /dimTime
148  )
149  );
150 }
151 
152 
153 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
154 
155 template<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 :
169  (
170  type,
171  alpha,
172  rho,
173  U,
174  alphaRhoPhi,
175  phi,
176  transport,
177  propertiesName
178  ),
179 
180  k_
181  (
182  IOobject
183  (
184  IOobject::groupName("k", this->alphaRhoPhi_.group()),
185  this->runTime_.timeName(),
186  this->mesh_,
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 
208 template<class BasicTurbulenceModel>
210 {
212  {
213  filter_.read(this->coeffDict());
214 
215  return true;
216  }
217 
218  return false;
219 }
220 
221 
222 template<class BasicTurbulenceModel>
224 {
225  return tmp<volScalarField>
226  (
227  new volScalarField
228  (
229  IOobject
230  (
231  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
232  this->runTime_.timeName(),
233  this->mesh_,
236  ),
237  Ce()*k()*sqrt(k())/this->delta()
238  )
239  );
240 }
241 
242 
243 template<class BasicTurbulenceModel>
245 {
246  volScalarField epsilon(Ce()*k()*sqrt(k())/this->delta());
247 
249  (
250  IOobject
251  (
252  IOobject::groupName("omega", this->alphaRhoPhi_.group()),
253  this->runTime_.timeName(),
254  this->mesh_
255  ),
256  epsilon/(0.09*k())
257  );
258 }
259 
260 
261 template<class BasicTurbulenceModel>
263 {
264  if (!this->turbulence_)
265  {
266  return;
267  }
268 
269  // Local references
270  const alphaField& alpha = this->alpha_;
271  const rhoField& rho = this->rho_;
272  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
273  const volVectorField& U = this->U_;
274  volScalarField& nut = this->nut_;
275  fv::options& fvOptions(fv::options::New(this->mesh_));
276 
278 
280 
282  const volSymmTensorField D(dev(symm(tgradU())));
283  const volScalarField G(this->GName(), 2.0*nut*(tgradU() && D));
284  tgradU.clear();
285 
286  volScalarField KK(0.5*(filter_(magSqr(U)) - magSqr(filter_(U))));
287  KK.max(dimensionedScalar("small", KK.dimensions(), SMALL));
288 
290  (
291  fvm::ddt(alpha, rho, k_)
292  + fvm::div(alphaRhoPhi, k_)
293  - fvm::laplacian(alpha*rho*DkEff(), k_)
294  ==
295  alpha*rho*G
296  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
297  - fvm::Sp(Ce(D, KK)*alpha*rho*sqrt(k_)/this->delta(), k_)
298  + kSource()
299  + fvOptions(alpha, rho, k_)
300  );
301 
302  kEqn.ref().relax();
303  fvOptions.constrain(kEqn.ref());
304  solve(kEqn);
305  fvOptions.correct(k_);
306  bound(k_, this->kMin_);
307 
308  correctNut(D, KK);
309 }
310 
311 
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313 
314 } // End namespace LESModels
315 } // End namespace Foam
316 
317 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::LESModels::dynamicKEqn::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: dynamicKEqn.H:140
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:323
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:291
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:104
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::LESModels::dynamicKEqn::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: dynamicKEqn.H:142
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
dynamicKEqn.H
rho
rho
Definition: readInitialConditions.H:88
Foam::LESModels::dynamicKEqn::correctNut
virtual void correctNut()
Definition: dynamicKEqn.C:127
Foam::LESModels::kEqn
One equation eddy-viscosity model.
Definition: kEqn.H:77
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::LESModels::dynamicKEqn::read
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicKEqn.C:209
Foam::LESModels::dynamicKEqn
Dynamic one equation eddy-viscosity model.
Definition: dynamicKEqn.H:79
Foam::LESModels::dynamicKEqn::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: dynamicKEqn.H:141
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
divU
zeroField divU
Definition: alphaSuSp.H:3
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::LESModels::LESeddyViscosity
Eddy viscosity LES SGS model base class.
Definition: LESeddyViscosity.H:58
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::LESModels::dynamicKEqn::correct
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicKEqn.C:262
U
U
Definition: pEqn.H:72
Foam::LESModels::dynamicKEqn::Ce
volScalarField Ce() const
Definition: dynamicKEqn.C:97
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::LESfilter::New
static autoPtr< LESfilter > New(const fvMesh &, const dictionary &, const word &filterDictName="filter")
Return a reference to the selected LES filter.
Definition: LESfilter.C:44
Foam::LESModels::dynamicKEqn::omega
virtual tmp< volScalarField > omega() const
Return sub-grid specific dissipation rate.
Definition: dynamicKEqn.C:244
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::LESModels::dynamicKEqn::Ck
volScalarField Ck(const volSymmTensorField &D, const volScalarField &KK) const
Calculate Ck by filtering the velocity field U.
Definition: dynamicKEqn.C:43
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::LESModels::dynamicKEqn::epsilon
virtual tmp< volScalarField > epsilon() const
Return sub-grid dissipation rate.
Definition: dynamicKEqn.C:223
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::GeometricField::max
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
Definition: GeometricField.C:1143
Foam::eddyViscosity< LESModel< BasicTurbulenceModel > >::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:134
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< symmTensor, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::LESModels::dynamicKEqn::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: dynamicKEqn.C:139
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106