realizableKE.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-2020 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 "realizableKE.H"
30 #include "fvOptions.H"
31 #include "bound.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
44 (
45  const volTensorField& gradU,
46  const volScalarField& S2,
47  const volScalarField& magS
48 )
49 {
50  tmp<volSymmTensorField> tS = dev(symm(gradU));
51  const volSymmTensorField& S = tS();
52 
54  (
55  (2*sqrt(2.0))*((S&S)&&S)
56  /(
57  magS*S2
58  + dimensionedScalar("small", dimensionSet(0, 0, -3, 0, 0), SMALL)
59  )
60  );
61 
62  tS.clear();
63 
65  (
66  (1.0/3.0)*acos(min(max(sqrt(6.0)*W, -scalar(1)), scalar(1)))
67  );
68  volScalarField As(sqrt(6.0)*cos(phis));
69  volScalarField Us(sqrt(S2/2.0 + magSqr(skew(gradU))));
70 
71  return 1.0/(A0_ + As*Us*k_/epsilon_);
72 }
73 
74 
75 template<class BasicTurbulenceModel>
77 (
78  const volTensorField& gradU,
79  const volScalarField& S2,
80  const volScalarField& magS
81 )
82 {
83  this->nut_ = rCmu(gradU, S2, magS)*sqr(k_)/epsilon_;
84  this->nut_.correctBoundaryConditions();
85  fv::options::New(this->mesh_).correct(this->nut_);
86 
87  BasicTurbulenceModel::correctNut();
88 }
89 
90 
91 template<class BasicTurbulenceModel>
93 {
94  tmp<volTensorField> tgradU = fvc::grad(this->U_);
95  volScalarField S2(2*magSqr(dev(symm(tgradU()))));
96  volScalarField magS(sqrt(S2));
97  correctNut(tgradU(), S2, magS);
98 }
99 
100 
101 template<class BasicTurbulenceModel>
103 {
104  return tmp<fvScalarMatrix>
105  (
106  new fvScalarMatrix
107  (
108  k_,
109  dimVolume*this->rho_.dimensions()*k_.dimensions()
111  )
112  );
113 }
114 
115 
116 template<class BasicTurbulenceModel>
118 {
119  return tmp<fvScalarMatrix>
120  (
121  new fvScalarMatrix
122  (
123  epsilon_,
124  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
125  /dimTime
126  )
127  );
128 }
129 
130 
131 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
132 
133 template<class BasicTurbulenceModel>
135 (
136  const alphaField& alpha,
137  const rhoField& rho,
138  const volVectorField& U,
139  const surfaceScalarField& alphaRhoPhi,
140  const surfaceScalarField& phi,
141  const transportModel& transport,
142  const word& propertiesName,
143  const word& type
144 )
145 :
146  eddyViscosity<RASModel<BasicTurbulenceModel>>
147  (
148  type,
149  alpha,
150  rho,
151  U,
152  alphaRhoPhi,
153  phi,
154  transport,
155  propertiesName
156  ),
157  A0_
158  (
159  dimensioned<scalar>::getOrAddToDict
160  (
161  "A0",
162  this->coeffDict_,
163  4.0
164  )
165  ),
166  C2_
167  (
168  dimensioned<scalar>::getOrAddToDict
169  (
170  "C2",
171  this->coeffDict_,
172  1.9
173  )
174  ),
175  sigmak_
176  (
177  dimensioned<scalar>::getOrAddToDict
178  (
179  "sigmak",
180  this->coeffDict_,
181  1.0
182  )
183  ),
184  sigmaEps_
185  (
186  dimensioned<scalar>::getOrAddToDict
187  (
188  "sigmaEps",
189  this->coeffDict_,
190  1.2
191  )
192  ),
193 
194  k_
195  (
196  IOobject
197  (
198  IOobject::groupName("k", U.group()),
199  this->runTime_.timeName(),
200  this->mesh_,
201  IOobject::MUST_READ,
202  IOobject::AUTO_WRITE
203  ),
204  this->mesh_
205  ),
206  epsilon_
207  (
208  IOobject
209  (
210  IOobject::groupName("epsilon", U.group()),
211  this->runTime_.timeName(),
212  this->mesh_,
213  IOobject::MUST_READ,
214  IOobject::AUTO_WRITE
215  ),
216  this->mesh_
217  )
218 {
219  bound(k_, this->kMin_);
220  bound(epsilon_, this->epsilonMin_);
221 
222  if (type == typeName)
223  {
224  this->printCoeffs(type);
225  }
226 }
227 
228 
229 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
230 
231 template<class BasicTurbulenceModel>
233 {
235  {
236  A0_.readIfPresent(this->coeffDict());
237  C2_.readIfPresent(this->coeffDict());
238  sigmak_.readIfPresent(this->coeffDict());
239  sigmaEps_.readIfPresent(this->coeffDict());
240 
241  return true;
242  }
243 
244  return false;
245 }
246 
247 
248 template<class BasicTurbulenceModel>
250 {
251  if (!this->turbulence_)
252  {
253  return;
254  }
255 
256  // Local references
257  const alphaField& alpha = this->alpha_;
258  const rhoField& rho = this->rho_;
259  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
260  const volVectorField& U = this->U_;
261  volScalarField& nut = this->nut_;
262  fv::options& fvOptions(fv::options::New(this->mesh_));
263 
264  eddyViscosity<RASModel<BasicTurbulenceModel>>::correct();
265 
267 
268  tmp<volTensorField> tgradU = fvc::grad(U);
269  volScalarField S2(2*magSqr(dev(symm(tgradU()))));
270  volScalarField magS(sqrt(S2));
271 
272  volScalarField eta(magS*k_/epsilon_);
273  volScalarField C1(max(eta/(scalar(5) + eta), scalar(0.43)));
274 
275  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
276 
277  // Update epsilon and G at the wall
278  epsilon_.boundaryFieldRef().updateCoeffs();
279 
280  // SAF: limiting thermo->nu(). If psiThermo is used rho might be < 0
281  // temporarily when p < 0 then nu < 0 which needs limiting
282  volScalarField nuLimited
283  (
284  max
285  (
286  this->nu(),
287  dimensionedScalar(this->nu()().dimensions(), Zero)
288  )
289  );
290 
291  // Dissipation equation
292  tmp<fvScalarMatrix> epsEqn
293  (
294  fvm::ddt(alpha, rho, epsilon_)
295  + fvm::div(alphaRhoPhi, epsilon_)
296  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
297  ==
298  C1*alpha*rho*magS*epsilon_
299  - fvm::Sp
300  (
301  C2_*alpha*rho*epsilon_/(k_ + sqrt(nuLimited*epsilon_)),
302  epsilon_
303  )
304  + epsilonSource()
305  + fvOptions(alpha, rho, epsilon_)
306  );
307 
308  epsEqn.ref().relax();
309  fvOptions.constrain(epsEqn.ref());
310  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
311  solve(epsEqn);
312  fvOptions.correct(epsilon_);
313  bound(epsilon_, this->epsilonMin_);
314 
315 
316  // Turbulent kinetic energy equation
317 
318  tmp<fvScalarMatrix> kEqn
319  (
320  fvm::ddt(alpha, rho, k_)
321  + fvm::div(alphaRhoPhi, k_)
322  - fvm::laplacian(alpha*rho*DkEff(), k_)
323  ==
324  alpha*rho*G
325  - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
326  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
327  + kSource()
328  + fvOptions(alpha, rho, k_)
329  );
330 
331  kEqn.ref().relax();
332  fvOptions.constrain(kEqn.ref());
333  solve(kEqn);
334  fvOptions.correct(k_);
335  bound(k_, this->kMin_);
336 
337  correctNut(tgradU(), S2, magS);
338 }
339 
340 
341 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
342 
343 } // End namespace RASModels
344 } // End namespace Foam
345 
346 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
realizableKE.H
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:77
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:190
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
Foam::RASModels::realizableKE::epsilonSource
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: realizableKE.C:110
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Definition: fvOptionListTemplates.C:348
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:280
Foam::constant::atomic::group
constexpr const char *const group
Definition: atomicConstants.H:46
Foam::skew
dimensionedTensor skew(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:131
Foam::constant::universal::G
const dimensionedScalar G
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:48
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Definition: fvOptions.C:96
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Definition: readThermalProperties.H:212
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Definition: bound.C:28
Foam::read
bool read(const char *buf, int32_t &val)
Definition: int32.H:125
Foam::RASModels::realizableKE::read
virtual bool read()
Definition: realizableKE.C:225
Foam::RASModels::realizableKE::k_
volScalarField k_
Definition: realizableKE.H:98
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:104
Foam::RASModel::epsilonMin_
dimensionedScalar epsilonMin_
Definition: RASModel.H:73
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Definition: hashSets.C:26
Foam::RASModels::realizableKE::rCmu
tmp< volScalarField > rCmu(const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS)
Definition: realizableKE.C:37
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::RASModels::realizableKE::epsilon_
volScalarField epsilon_
Definition: realizableKE.H:99
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:220
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::RASModels::realizableKE::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: realizableKE.C:95
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
rho
rho
Definition: readInitialConditions.H:88
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
divU
zeroField divU
Definition: alphaSuSp.H:3
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:36
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:51
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:94
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:36
Foam
Definition: atmBoundaryLayer.C:26
Foam::RASModels::realizableKE::correctNut
virtual void correctNut()
Definition: realizableKE.C:85
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:93
Foam::transportModel
Base-class for all transport models used by the incompressible turbulence models.
Definition: transportModel.H:49
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:95
U
U
Definition: pEqn.H:72
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:44
Foam::RASModel::kMin_
dimensionedScalar kMin_
Definition: RASModel.H:70
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:50
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:137
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Us
Us
Definition: createFaFields.H:51
Foam::RASModels::realizableKE::realizableKE
realizableKE(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)
Definition: realizableKE.C:128
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Definition: POSIX.C:717
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:64
Foam::RASModel::printCoeffs
virtual void printCoeffs(const word &type)
Definition: RASModel.C:27
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:51
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
Foam::RASModels::realizableKE::correct
virtual void correct()
Definition: realizableKE.C:242
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:88
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Definition: fvcMeshPhi.C:183
phis
edgeScalarField phis(IOobject("phis", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), linearEdgeInterpolate(Us) &aMesh.Le())
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:181
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:99
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:258