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 | Copyright (C) 2011-2015 OpenFOAM Foundation
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 \*---------------------------------------------------------------------------*/
25 
26 #include "dynamicKEqn.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 template<class BasicTurbulenceModel>
39 (
40  const volSymmTensorField& D,
41  const volScalarField& KK
42 ) const
43 {
44  const volSymmTensorField LL
45  (
46  simpleFilter_(dev(filter_(sqr(this->U_)) - (sqr(filter_(this->U_)))))
47  );
48 
49  const volSymmTensorField MM
50  (
51  simpleFilter_
52  (
53  -2.0*this->delta()*sqrt
54  (
55  max(KK, dimensionedScalar("zero", KK.dimensions(), 0.0))
56  )*filter_(D)
57  )
58  );
59 
60  const volScalarField Ck
61  (
62  simpleFilter_(0.5*(LL && MM))
63  /(
64  simpleFilter_(magSqr(MM))
65  + dimensionedScalar("small", sqr(MM.dimensions()), VSMALL)
66  )
67  );
68 
69  tmp<volScalarField> tfld = 0.5*(mag(Ck) + Ck);
70  return tfld();
71 }
72 
73 
74 template<class BasicTurbulenceModel>
76 (
77  const volSymmTensorField& D,
78  const volScalarField& KK
79 ) const
80 {
81  const volScalarField Ce
82  (
83  simpleFilter_(this->nuEff()*(filter_(magSqr(D)) - magSqr(filter_(D))))
84  /simpleFilter_(pow(KK, 1.5)/(2.0*this->delta()))
85  );
86 
87  tmp<volScalarField> tfld = 0.5*(mag(Ce) + Ce);
88  return tfld();
89 }
90 
91 
92 template<class BasicTurbulenceModel>
94 {
95  const volSymmTensorField D(dev(symm(fvc::grad(this->U_))));
96 
98  (
99  0.5*(filter_(magSqr(this->U_)) - magSqr(filter_(this->U_)))
100  );
101  KK.max(dimensionedScalar("small", KK.dimensions(), SMALL));
102 
103  return Ce(D, KK);
104 }
105 
106 
107 template<class BasicTurbulenceModel>
109 (
110  const volSymmTensorField& D,
111  const volScalarField& KK
112 )
113 {
114  this->nut_ = Ck(D, KK)*sqrt(k_)*this->delta();
115  this->nut_.correctBoundaryConditions();
116 }
117 
118 
119 template<class BasicTurbulenceModel>
121 {
122  const volScalarField KK
123  (
124  0.5*(filter_(magSqr(this->U_)) - magSqr(filter_(this->U_)))
125  );
126 
127  correctNut(symm(fvc::grad(this->U_)), KK);
128 }
129 
130 
131 template<class BasicTurbulenceModel>
133 {
134  return tmp<fvScalarMatrix>
135  (
136  new fvScalarMatrix
137  (
138  k_,
139  dimVolume*this->rho_.dimensions()*k_.dimensions()
140  /dimTime
141  )
142  );
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
147 
148 template<class BasicTurbulenceModel>
150 (
151  const alphaField& alpha,
152  const rhoField& rho,
153  const volVectorField& U,
154  const surfaceScalarField& alphaRhoPhi,
155  const surfaceScalarField& phi,
156  const transportModel& transport,
157  const word& propertiesName,
158  const word& type
159 )
160 :
162  (
163  type,
164  alpha,
165  rho,
166  U,
167  alphaRhoPhi,
168  phi,
169  transport,
170  propertiesName
171  ),
172 
173  k_
174  (
175  IOobject
176  (
177  IOobject::groupName("k", this->U_.group()),
178  this->runTime_.timeName(),
179  this->mesh_,
182  ),
183  this->mesh_
184  ),
185 
186  simpleFilter_(this->mesh_),
187  filterPtr_(LESfilter::New(this->mesh_, this->coeffDict())),
188  filter_(filterPtr_())
189 {
190  bound(k_, this->kMin_);
191 
192  if (type == typeName)
193  {
194  this->printCoeffs(type);
195  }
196 }
197 
198 
199 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 
201 template<class BasicTurbulenceModel>
203 {
205  {
206  filter_.read(this->coeffDict());
207 
208  return true;
209  }
210  else
211  {
212  return false;
213  }
214 }
215 
216 
217 template<class BasicTurbulenceModel>
219 {
220  return tmp<volScalarField>
221  (
222  new volScalarField
223  (
224  IOobject
225  (
226  IOobject::groupName("epsilon", this->U_.group()),
227  this->runTime_.timeName(),
228  this->mesh_,
231  ),
232  Ce()*k()*sqrt(k())/this->delta()
233  )
234  );
235 }
236 
237 
238 template<class BasicTurbulenceModel>
240 {
241  if (!this->turbulence_)
242  {
243  return;
244  }
245 
246  // Local references
247  const alphaField& alpha = this->alpha_;
248  const rhoField& rho = this->rho_;
249  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
250  const volVectorField& U = this->U_;
251  volScalarField& nut = this->nut_;
252 
254 
256 
258  const volSymmTensorField D(dev(symm(tgradU())));
259  const volScalarField G(this->GName(), 2.0*nut*(tgradU() && D));
260  tgradU.clear();
261 
262  volScalarField KK(0.5*(filter_(magSqr(U)) - magSqr(filter_(U))));
263  KK.max(dimensionedScalar("small", KK.dimensions(), SMALL));
264 
266  (
267  fvm::ddt(alpha, rho, k_)
268  + fvm::div(alphaRhoPhi, k_)
269  - fvm::laplacian(alpha*rho*DkEff(), k_)
270  ==
271  alpha*rho*G
272  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
273  - fvm::Sp(Ce(D, KK)*alpha*rho*sqrt(k_)/this->delta(), k_)
274  + kSource()
275  );
276 
277  kEqn().relax();
278  kEqn().solve();
279  bound(k_, this->kMin_);
280 
281  correctNut(D, KK);
282 }
283 
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 } // End namespace LESModels
288 } // End namespace Foam
289 
290 // ************************************************************************* //
Foam::LESModels::dynamicKEqn::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: dynamicKEqn.H:135
Foam::LESModels::dynamicKEqn::dynamicKEqn
dynamicKEqn(const dynamicKEqn &)
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:82
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Foam::LESModels::dynamicKEqn::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: dynamicKEqn.H:137
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
dynamicKEqn.H
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
U
U
Definition: pEqn.H:46
Foam::LESModels::dynamicKEqn::correctNut
virtual void correctNut()
Definition: dynamicKEqn.C:120
Foam::LESModels::kEqn
One equation eddy-viscosity model.
Definition: kEqn.H:74
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::LESModels::dynamicKEqn::read
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicKEqn.C:202
Foam::LESModels::dynamicKEqn::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: dynamicKEqn.H:136
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
nut
nut
Definition: createTDFields.H:71
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::LESModels::LESeddyViscosity
Eddy viscosity LES SGS model base class.
Definition: LESeddyViscosity.H:55
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::LESModels::dynamicKEqn::correct
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicKEqn.C:239
rho
rho
Definition: pEqn.H:3
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::GeometricField::max
void max(const dimensioned< Type > &)
Foam::LESModels::dynamicKEqn::Ce
volScalarField Ce() const
Definition: dynamicKEqn.C:93
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
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:41
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::LESModels::dynamicKEqn::Ck
volScalarField Ck(const volSymmTensorField &D, const volScalarField &KK) const
Calculate Ck by filtering the velocity field U.
Definition: dynamicKEqn.C:39
Foam::LESModels::dynamicKEqn::epsilon
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: dynamicKEqn.C:218
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::eddyViscosity< LESModel< BasicTurbulenceModel > >::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:132
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::tmp::clear
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:172
Foam::LESModels::dynamicKEqn::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: dynamicKEqn.C:132
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104