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