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-2023 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  (
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(devSymm(fvc::grad(this->U_)));
100 
101  volScalarField KK
102  (
103  0.5*(filter_(magSqr(this->U_)) - magSqr(filter_(this->U_)))
104  );
105  KK.clamp_min(SMALL);
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 {
142  (
143  k_,
144  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
145  );
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
150 
151 template<class BasicTurbulenceModel>
153 (
154  const alphaField& alpha,
155  const rhoField& rho,
156  const volVectorField& U,
157  const surfaceScalarField& alphaRhoPhi,
158  const surfaceScalarField& phi,
159  const transportModel& transport,
160  const word& propertiesName,
161  const word& type
162 )
163 :
164  LESeddyViscosity<BasicTurbulenceModel>
165  (
166  type,
167  alpha,
168  rho,
169  U,
170  alphaRhoPhi,
171  phi,
172  transport,
173  propertiesName
174  ),
175 
176  k_
177  (
178  IOobject
179  (
180  IOobject::groupName("k", this->alphaRhoPhi_.group()),
181  this->runTime_.timeName(),
182  this->mesh_,
183  IOobject::MUST_READ,
184  IOobject::AUTO_WRITE
185  ),
186  this->mesh_
187  ),
188 
189  simpleFilter_(this->mesh_),
190  filterPtr_(LESfilter::New(this->mesh_, this->coeffDict())),
191  filter_(filterPtr_())
192 {
193  bound(k_, this->kMin_);
194 
195  if (type == typeName)
196  {
197  this->printCoeffs(type);
198  }
199 }
200 
201 
202 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
203 
204 template<class BasicTurbulenceModel>
206 {
208  {
209  filter_.read(this->coeffDict());
210 
211  return true;
212  }
213 
214  return false;
215 }
216 
217 
218 template<class BasicTurbulenceModel>
220 {
221  if (!this->turbulence_)
222  {
223  return;
224  }
225 
226  // Local references
227  const alphaField& alpha = this->alpha_;
228  const rhoField& rho = this->rho_;
229  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
230  const volVectorField& U = this->U_;
231  volScalarField& nut = this->nut_;
232  fv::options& fvOptions(fv::options::New(this->mesh_));
233 
235 
237 
238  tmp<volTensorField> tgradU(fvc::grad(U));
239  const volSymmTensorField D(devSymm(tgradU()));
240  const volScalarField G(this->GName(), 2.0*nut*(tgradU() && D));
241  tgradU.clear();
242 
243  volScalarField KK(0.5*(filter_(magSqr(U)) - magSqr(filter_(U))));
244  KK.clamp_min(SMALL);
245 
246  tmp<fvScalarMatrix> kEqn
247  (
248  fvm::ddt(alpha, rho, k_)
249  + fvm::div(alphaRhoPhi, k_)
250  - fvm::laplacian(alpha*rho*DkEff(), k_)
251  ==
252  alpha*rho*G
253  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
254  - fvm::Sp(Ce(D, KK)*alpha*rho*sqrt(k_)/this->delta(), k_)
255  + kSource()
256  + fvOptions(alpha, rho, k_)
257  );
258 
259  kEqn.ref().relax();
260  fvOptions.constrain(kEqn.ref());
261  solve(kEqn);
262  fvOptions.correct(k_);
263  bound(k_, this->kMin_);
264 
265  correctNut(D, KK);
266 }
267 
268 
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
270 
271 } // End namespace LESModels
272 } // End namespace Foam
273 
274 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
void clamp_min(const Type &lower)
Impose lower (floor) clamp on the field values (in-place)
scalar delta
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicKEqn.C:198
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:84
Abstract class for LES filters.
Definition: LESfilter.H:53
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicKEqn.C:212
const dimensionedScalar G
Newtonian constant of gravitation.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:27
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:481
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
fv::options & fvOptions
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
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:36
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
volScalarField Ce() const
Definition: dynamicKEqn.C:90
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
virtual bool read()
Read model coefficients if they have changed.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Dynamic one equation eddy-viscosity model.
Definition: dynamicKEqn.H:75
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:183
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:29
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:41
Base-class for all transport models used by the incompressible turbulence models. ...
dimensionedScalar kMin_
Lower limit of k.
Definition: LESModel.H:92
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual tmp< fvScalarMatrix > kSource() const
Definition: dynamicKEqn.C:132
const dimensionedScalar & D
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
scalar nut
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127