dynamicKEqn.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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  return max
42  (
43  0.5*(filter_(magSqr(this->U_)) - magSqr(filter_(this->U_))),
45  );
46 }
47 
48 
49 template<class BasicTurbulenceModel>
51 (
52  const volSymmTensorField& D,
53  const volScalarField& KK
54 ) const
55 {
56  const volSymmTensorField LL
57  (
58  simpleFilter_(dev(filter_(sqr(this->U_)) - (sqr(filter_(this->U_)))))
59  );
60 
61  const volSymmTensorField MM
62  (
63  simpleFilter_(-2.0*this->delta()*sqrt(KK)*filter_(D))
64  );
65 
66  const volScalarField Ck
67  (
68  simpleFilter_(0.5*(LL && MM))
69  /(
70  simpleFilter_(magSqr(MM))
71  + dimensionedScalar(sqr(MM.dimensions()), vSmall)
72  )
73  );
74 
75  tmp<volScalarField> tfld = 0.5*(mag(Ck) + Ck);
76  return tfld();
77 }
78 
79 
80 template<class BasicTurbulenceModel>
82 (
83  const volSymmTensorField& D,
84  const volScalarField& KK
85 ) const
86 {
87  const volScalarField Ce
88  (
89  simpleFilter_(this->nuEff()*(filter_(magSqr(D)) - magSqr(filter_(D))))
90  /simpleFilter_(pow(KK, 1.5)/(2.0*this->delta()))
91  );
92 
93  tmp<volScalarField> tfld = 0.5*(mag(Ce) + Ce);
94  return tfld();
95 }
96 
97 
98 template<class BasicTurbulenceModel>
100 {
101  const volSymmTensorField D(dev(symm(fvc::grad(this->U_))));
102  return Ce(D, KK());
103 }
104 
105 
106 template<class BasicTurbulenceModel>
108 (
109  const volSymmTensorField& D,
110  const volScalarField& KK
111 )
112 {
113  this->nut_ = Ck(D, KK)*sqrt(k_)*this->delta();
114  this->nut_.correctBoundaryConditions();
115  fv::options::New(this->mesh_).correct(this->nut_);
116 
117  BasicTurbulenceModel::correctNut();
118 }
119 
120 
121 template<class BasicTurbulenceModel>
123 {
124  correctNut(symm(fvc::grad(this->U_)), KK());
125 }
126 
127 
128 template<class BasicTurbulenceModel>
130 {
131  return tmp<fvScalarMatrix>
132  (
133  new fvScalarMatrix
134  (
135  k_,
136  dimVolume*this->rho_.dimensions()*k_.dimensions()
137  /dimTime
138  )
139  );
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
144 
145 template<class BasicTurbulenceModel>
147 (
148  const alphaField& alpha,
149  const rhoField& rho,
150  const volVectorField& U,
151  const surfaceScalarField& alphaRhoPhi,
152  const surfaceScalarField& phi,
153  const transportModel& transport,
154  const word& propertiesName,
155  const word& type
156 )
157 :
159  (
160  type,
161  alpha,
162  rho,
163  U,
164  alphaRhoPhi,
165  phi,
166  transport,
167  propertiesName
168  ),
169 
170  k_
171  (
172  IOobject
173  (
174  IOobject::groupName("k", this->alphaRhoPhi_.group()),
175  this->runTime_.timeName(),
176  this->mesh_,
179  ),
180  this->mesh_
181  ),
182 
183  simpleFilter_(this->mesh_),
184  filterPtr_(LESfilter::New(this->mesh_, this->coeffDict())),
185  filter_(filterPtr_())
186 {
187  bound(k_, this->kMin_);
188 
189  if (type == typeName)
190  {
191  this->printCoeffs(type);
192  }
193 }
194 
195 
196 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
197 
198 template<class BasicTurbulenceModel>
200 {
202  {
203  filter_.read(this->coeffDict());
204 
205  return true;
206  }
207  else
208  {
209  return false;
210  }
211 }
212 
213 
214 template<class BasicTurbulenceModel>
216 {
217  return volScalarField::New
218  (
219  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
220  Ce()*k()*sqrt(k())/this->delta()
221  );
222 }
223 
224 
225 template<class BasicTurbulenceModel>
227 {
228  if (!this->turbulence_)
229  {
230  return;
231  }
232 
233  // Local references
234  const alphaField& alpha = this->alpha_;
235  const rhoField& rho = this->rho_;
236  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
237  const volVectorField& U = this->U_;
238  volScalarField& nut = this->nut_;
239  fv::options& fvOptions(fv::options::New(this->mesh_));
240 
242 
244 
245  tmp<volTensorField> tgradU(fvc::grad(U));
246  const volSymmTensorField D(dev(symm(tgradU())));
247  const volScalarField G(this->GName(), 2.0*nut*(tgradU() && D));
248  tgradU.clear();
249 
250  const volScalarField KK(this->KK());
251 
253  (
254  fvm::ddt(alpha, rho, k_)
255  + fvm::div(alphaRhoPhi, k_)
256  - fvm::laplacian(alpha*rho*DkEff(), k_)
257  ==
258  alpha*rho*G
259  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
260  - fvm::Sp(Ce(D, KK)*alpha*rho*sqrt(k_)/this->delta(), k_)
261  + kSource()
262  + fvOptions(alpha, rho, k_)
263  );
264 
265  kEqn.ref().relax();
266  fvOptions.constrain(kEqn.ref());
267  solve(kEqn);
268  fvOptions.correct(k_);
269  bound(k_, this->kMin_);
270 
271  correctNut(D, KK);
272 }
273 
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 } // End namespace LESModels
278 } // End namespace Foam
279 
280 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
scalar delta
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicKEqn.C:199
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
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
fv::options & fvOptions
surfaceScalarField & phi
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicKEqn.C:226
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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:129
Finite-volume options.
Definition: fvOptions.H:52
One equation eddy-viscosity model.
Definition: kEqn.H:71
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
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:51
A class for handling words, derived from string.
Definition: word.H:59
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:130
BasicTurbulenceModel::alphaField alphaField
Definition: dynamicKEqn.H:128
volScalarField Ce() const
Definition: dynamicKEqn.C:99
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
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
U
Definition: pEqn.H:72
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
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:215
dynamicKEqn(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: dynamicKEqn.C:147
virtual tmp< fvScalarMatrix > kSource() const
Definition: dynamicKEqn.C:129
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< volScalarField > KK() const
Return the KK field obtained by filtering the velocity field U.
Definition: dynamicKEqn.C:39
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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:101
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.
const dimensionSet dimVelocity