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-2020 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 BasicMomentumTransportModel>
41 {
42  return max
43  (
44  0.5*(filter_(magSqr(this->U_)) - magSqr(filter_(this->U_))),
46  );
47 }
48 
49 
50 template<class BasicMomentumTransportModel>
52 (
53  const volSymmTensorField& D,
54  const volScalarField& KK
55 ) const
56 {
57  const volSymmTensorField LL
58  (
59  simpleFilter_(dev(filter_(sqr(this->U_)) - (sqr(filter_(this->U_)))))
60  );
61 
62  const volSymmTensorField MM
63  (
64  simpleFilter_(-2.0*this->delta()*sqrt(KK)*filter_(D))
65  );
66 
67  const volScalarField Ck
68  (
69  simpleFilter_(0.5*(LL && MM))
70  /(
71  simpleFilter_(magSqr(MM))
72  + dimensionedScalar(sqr(MM.dimensions()), vSmall)
73  )
74  );
75 
76  tmp<volScalarField> tfld = 0.5*(mag(Ck) + Ck);
77  return tfld();
78 }
79 
80 
81 template<class BasicMomentumTransportModel>
83 (
84  const volSymmTensorField& D,
85  const volScalarField& KK
86 ) const
87 {
88  const volScalarField Ce
89  (
90  simpleFilter_(this->nuEff()*(filter_(magSqr(D)) - magSqr(filter_(D))))
91  /simpleFilter_(pow(KK, 1.5)/(2.0*this->delta()))
92  );
93 
94  tmp<volScalarField> tfld = 0.5*(mag(Ce) + Ce);
95  return tfld();
96 }
97 
98 
99 template<class BasicMomentumTransportModel>
101 {
102  const volSymmTensorField D(dev(symm(fvc::grad(this->U_))));
103  return Ce(D, KK());
104 }
105 
106 
107 template<class BasicMomentumTransportModel>
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  fv::options::New(this->mesh_).correct(this->nut_);
117 }
118 
119 
120 template<class BasicMomentumTransportModel>
122 {
123  correctNut(symm(fvc::grad(this->U_)), KK());
124 }
125 
126 
127 template<class BasicMomentumTransportModel>
129 {
130  return tmp<fvScalarMatrix>
131  (
132  new fvScalarMatrix
133  (
134  k_,
135  dimVolume*this->rho_.dimensions()*k_.dimensions()
136  /dimTime
137  )
138  );
139 }
140 
141 
142 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
143 
144 template<class BasicMomentumTransportModel>
146 (
147  const alphaField& alpha,
148  const rhoField& rho,
149  const volVectorField& U,
150  const surfaceScalarField& alphaRhoPhi,
151  const surfaceScalarField& phi,
152  const transportModel& transport,
153  const word& type
154 )
155 :
157  (
158  type,
159  alpha,
160  rho,
161  U,
162  alphaRhoPhi,
163  phi,
164  transport
165  ),
166 
167  k_
168  (
169  IOobject
170  (
171  IOobject::groupName("k", this->alphaRhoPhi_.group()),
172  this->runTime_.timeName(),
173  this->mesh_,
176  ),
177  this->mesh_
178  ),
179 
180  simpleFilter_(this->mesh_),
181  filterPtr_(LESfilter::New(this->mesh_, this->coeffDict())),
182  filter_(filterPtr_())
183 {
184  bound(k_, this->kMin_);
185 
186  if (type == typeName)
187  {
188  this->printCoeffs(type);
189  }
190 }
191 
192 
193 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
194 
195 template<class BasicMomentumTransportModel>
197 {
199  {
200  filter_.read(this->coeffDict());
201 
202  return true;
203  }
204  else
205  {
206  return false;
207  }
208 }
209 
210 
211 template<class BasicMomentumTransportModel>
213 {
214  return volScalarField::New
215  (
216  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
217  Ce()*k()*sqrt(k())/this->delta()
218  );
219 }
220 
221 
222 template<class BasicMomentumTransportModel>
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_;
236  fv::options& fvOptions(fv::options::New(this->mesh_));
237 
239 
241 
242  tmp<volTensorField> tgradU(fvc::grad(U));
243  const volSymmTensorField D(dev(symm(tgradU())));
244  const volScalarField G(this->GName(), 2.0*nut*(tgradU() && D));
245  tgradU.clear();
246 
247  const volScalarField KK(this->KK());
248 
250  (
251  fvm::ddt(alpha, rho, k_)
252  + fvm::div(alphaRhoPhi, k_)
253  - fvm::laplacian(alpha*rho*DkEff(), k_)
254  ==
255  alpha*rho*G
256  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
257  - fvm::Sp(Ce(D, KK)*alpha*rho*sqrt(k_)/this->delta(), k_)
258  + kSource()
259  + fvOptions(alpha, rho, k_)
260  );
261 
262  kEqn.ref().relax();
263  fvOptions.constrain(kEqn.ref());
264  solve(kEqn);
265  fvOptions.correct(k_);
266  bound(k_, this->kMin_);
267 
268  correctNut(D, KK);
269 }
270 
271 
272 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
273 
274 } // End namespace LESModels
275 } // End namespace Foam
276 
277 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
scalar delta
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicKEqn.C:196
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
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicKEqn.C:223
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
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
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Finite-volume options.
Definition: fvOptions.H:52
One equation eddy-viscosity model.
Definition: kEqn.H:71
BasicMomentumTransportModel::transportModel transportModel
Definition: dynamicKEqn.H:130
phi
Definition: pEqn.H:104
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:52
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
volScalarField Ce() const
Definition: dynamicKEqn.C:100
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:212
virtual tmp< fvScalarMatrix > kSource() const
Definition: dynamicKEqn.C:128
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:40
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
BasicMomentumTransportModel::alphaField alphaField
Definition: dynamicKEqn.H:128
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
dynamicKEqn(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: dynamicKEqn.C:146
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionedScalar & G
Newtonian constant of gravitation.
scalar nut
Namespace for OpenFOAM.
BasicMomentumTransportModel::rhoField rhoField
Definition: dynamicKEqn.H:129
const dimensionSet dimVelocity