kEpsilonPhitF.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) 2019-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "kEpsilonPhitF.H"
29 #include "fvOptions.H"
30 #include "bound.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace RASModels
37 {
38 
39 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
40 
41 template<class BasicTurbulenceModel>
43 {
44  // (LUU:p. 173)
45  this->nut_ = Cmu_*phit_*k_*T_;
46  this->nut_.correctBoundaryConditions();
47  fv::options::New(this->mesh_).correct(this->nut_);
48 
49  BasicTurbulenceModel::correctNut();
50 }
51 
52 
53 template<class BasicTurbulenceModel>
55 {
56  // (LUU:Eq. 7)
57  return
58  max
59  (
60  k_/epsilon_,
61  CT_*sqrt
62  (
63  max
64  (
65  this->nu(),
66  dimensionedScalar(this->nu()().dimensions(), Zero)
67  )/epsilon_
68  )
69  );
70 }
71 
72 
73 template<class BasicTurbulenceModel>
75 {
76  // (LUU:Eq. 7)
77  return
78  CL_*max
79  (
80  pow(k_, 1.5)/epsilon_,
81  Ceta_*pow025
82  (
83  pow3
84  (
85  max
86  (
87  this->nu(),
88  dimensionedScalar(this->nu()().dimensions(), Zero)
89  )
90  )/epsilon_
91  )
92  );
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97 
98 template<class BasicTurbulenceModel>
100 (
101  const alphaField& alpha,
102  const rhoField& rho,
103  const volVectorField& U,
104  const surfaceScalarField& alphaRhoPhi,
105  const surfaceScalarField& phi,
106  const transportModel& transport,
107  const word& propertiesName,
108  const word& type
109 )
110 :
112  (
113  type,
114  alpha,
115  rho,
116  U,
117  alphaRhoPhi,
118  phi,
119  transport,
120  propertiesName
121  ),
122 
123  includeNu_
124  (
126  (
127  "includeNu",
128  this->coeffDict_,
129  true
130  )
131  ),
132  Cmu_
133  (
135  (
136  "Cmu",
137  this->coeffDict_,
138  0.22
139  )
140  ),
141  Ceps1a_
142  (
144  (
145  "Ceps1a",
146  this->coeffDict_,
147  1.4
148  )
149  ),
150  Ceps1b_
151  (
153  (
154  "Ceps1b",
155  this->coeffDict_,
156  1.0
157  )
158  ),
159  Ceps1c_
160  (
162  (
163  "Ceps1c",
164  this->coeffDict_,
165  0.05
166  )
167  ),
168  Ceps2_
169  (
171  (
172  "Ceps2",
173  this->coeffDict_,
174  1.9
175  )
176  ),
177  Cf1_
178  (
180  (
181  "Cf1",
182  this->coeffDict_,
183  1.4
184  )
185  ),
186  Cf2_
187  (
189  (
190  "Cf2",
191  this->coeffDict_,
192  0.3
193  )
194  ),
195  CL_
196  (
198  (
199  "CL",
200  this->coeffDict_,
201  0.25
202  )
203  ),
204  Ceta_
205  (
207  (
208  "Ceta",
209  this->coeffDict_,
210  110.0
211  )
212  ),
213  CT_
214  (
216  (
217  "CT",
218  this->coeffDict_,
219  6.0
220  )
221  ),
222  sigmaK_
223  (
225  (
226  "sigmaK",
227  this->coeffDict_,
228  1.0
229  )
230  ),
231  sigmaEps_
232  (
234  (
235  "sigmaEps",
236  this->coeffDict_,
237  1.3
238  )
239  ),
240  sigmaPhit_
241  (
243  (
244  "sigmaPhit",
245  this->coeffDict_,
246  1.0
247  )
248  ),
249 
250  k_
251  (
252  IOobject
253  (
254  IOobject::groupName("k", alphaRhoPhi.group()),
255  this->runTime_.timeName(),
256  this->mesh_,
259  ),
260  this->mesh_
261  ),
262  epsilon_
263  (
264  IOobject
265  (
266  IOobject::groupName("epsilon", alphaRhoPhi.group()),
267  this->runTime_.timeName(),
268  this->mesh_,
271  ),
272  this->mesh_
273  ),
274  phit_
275  (
276  IOobject
277  (
278  IOobject::groupName("phit", alphaRhoPhi.group()),
279  this->runTime_.timeName(),
280  this->mesh_,
283  ),
284  this->mesh_
285  ),
286  f_
287  (
288  IOobject
289  (
290  IOobject::groupName("f", alphaRhoPhi.group()),
291  this->runTime_.timeName(),
292  this->mesh_,
295  ),
296  this->mesh_
297  ),
298  T_
299  (
300  IOobject
301  (
302  "T",
303  this->runTime_.timeName(),
304  this->mesh_,
307  false
308  ),
309  this->mesh_,
311  ),
312 
313  phitMin_(dimensionedScalar("phitMin", phit_.dimensions(), SMALL)),
314  fMin_(dimensionedScalar("fMin", f_.dimensions(), SMALL)),
315  TMin_(dimensionedScalar("TMin", dimTime, SMALL)),
316  L2Min_(dimensionedScalar("L2Min", sqr(dimLength), SMALL))
317 {
318  bound(k_, this->kMin_);
319  bound(epsilon_, this->epsilonMin_);
320  bound(phit_, phitMin_);
321  bound(f_, fMin_);
322 
323  if (type == typeName)
324  {
325  this->printCoeffs(type);
326  }
327 
328  if
329  (
330  mag(sigmaK_.value()) < VSMALL
331  || mag(sigmaEps_.value()) < VSMALL
332  || mag(sigmaPhit_.value()) < VSMALL
333  )
334  {
336  << "Non-zero values are required for the model constants:" << nl
337  << "sigmaK = " << sigmaK_ << nl
338  << "sigmaEps = " << sigmaEps_ << nl
339  << "sigmaPhit = " << sigmaPhit_ << nl
340  << exit(FatalError);
341  }
342 }
343 
344 
345 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
346 
347 template<class BasicTurbulenceModel>
349 {
351  {
352  includeNu_.readIfPresent("includeNu", this->coeffDict());
353  Cmu_.readIfPresent(this->coeffDict());
354  Ceps1a_.readIfPresent(this->coeffDict());
355  Ceps1b_.readIfPresent(this->coeffDict());
356  Ceps1c_.readIfPresent(this->coeffDict());
357  Ceps2_.readIfPresent(this->coeffDict());
358  Cf1_.readIfPresent(this->coeffDict());
359  Cf2_.readIfPresent(this->coeffDict());
360  CL_.readIfPresent(this->coeffDict());
361  Ceta_.readIfPresent(this->coeffDict());
362  CT_.readIfPresent(this->coeffDict());
363  sigmaK_.readIfPresent(this->coeffDict());
364  sigmaEps_.readIfPresent(this->coeffDict());
365  sigmaPhit_.readIfPresent(this->coeffDict());
366 
367  return true;
368  }
369 
370  return false;
371 }
372 
373 
374 template<class BasicTurbulenceModel>
376 {
377  if (!this->turbulence_)
378  {
379  return;
380  }
381 
382  // Construct local convenience references
383  const alphaField& alpha = this->alpha_;
384  const rhoField& rho = this->rho_;
385  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
386  const volVectorField& U = this->U_;
387  volScalarField& nut = this->nut_;
388  fv::options& fvOptions(fv::options::New(this->mesh_));
389 
391 
393  (
394  fvc::div(fvc::absolute(this->phi(), U))
395  );
396 
398  volScalarField G(this->GName(), nut*(2.0*(dev(tS()) && tS())));
399  tS.clear();
400 
401  T_ = Ts();
402  bound(T_, TMin_);
403 
404  const volScalarField L2(type() + "L2", sqr(Ls()) + L2Min_);
405 
406  const volScalarField::Internal Ceps1
407  (
408  "Ceps1",
409  Ceps1a_*(Ceps1b_ + Ceps1c_*sqrt(1.0/phit_()))
410  );
411 
412 
413  // Update epsilon and G at the wall
414  epsilon_.boundaryFieldRef().updateCoeffs();
415 
416  // Turbulent kinetic energy dissipation rate equation (LUU:Eq. 4)
417  // k/T ~ epsilon
418  tmp<fvScalarMatrix> epsEqn
419  (
420  fvm::ddt(alpha, rho, epsilon_)
421  + fvm::div(alphaRhoPhi, epsilon_)
422  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
423  ==
424  alpha()*rho()*Ceps1*G()/T_()
425  - fvm::SuSp
426  (
427  (2.0/3.0*Ceps1)*(alpha()*rho()*divU),
428  epsilon_
429  )
430  - fvm::Sp(alpha()*rho()*Ceps2_/T_(), epsilon_)
431  + fvOptions(alpha, rho, epsilon_)
432  );
433 
434  epsEqn.ref().relax();
435  fvOptions.constrain(epsEqn.ref());
436  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
437  solve(epsEqn);
438  fvOptions.correct(epsilon_);
439  bound(epsilon_, this->epsilonMin_);
440 
441 
442  // Turbulent kinetic energy equation (LUU:Eq. 3)
443  // epsilon/k ~ 1/Ts
445  (
446  fvm::ddt(alpha, rho, k_)
447  + fvm::div(alphaRhoPhi, k_)
448  - fvm::laplacian(alpha*rho*DkEff(), k_)
449  ==
450  alpha()*rho()*G()
451  - fvm::SuSp(2.0/3.0*alpha()*rho()*divU, k_)
452  - fvm::Sp(alpha()*rho()*(1.0/T_()), k_)
453  + fvOptions(alpha, rho, k_)
454  );
455 
456  kEqn.ref().relax();
457  fvOptions.constrain(kEqn.ref());
458  solve(kEqn);
459  fvOptions.correct(k_);
460  bound(k_, this->kMin_);
461 
462 
463  // Elliptic relaxation function equation (LUU:Eq. 18)
464  // All source terms are non-negative functions (LUU:p. 176)
466  (
467  - fvm::laplacian(f_)
468  ==
469  - fvm::Sp(1.0/L2(), f_)
470  - (
471  (Cf1_ - 1.0)*(phit_() - 2.0/3.0)/T_()
472  -(Cf2_*G())/k_()
473  +(Cf2_*(2.0/3.0)*divU)
474  -(2.0*this->nu()*(fvc::grad(phit_) & fvc::grad(k_)))()/k_()
475  -(this->nu()*fvc::laplacian(phit_))()
476  )/L2()
477  );
478 
479  fEqn.ref().relax();
480  solve(fEqn);
481  bound(f_, fMin_);
482 
483 
484  // Normalised wall-normal fluctuating velocity scale equation (LUU:Eq. 17)
485  // All source terms are non-negative functions (LUU:p. 176)
486  tmp<fvScalarMatrix> phitEqn
487  (
488  fvm::ddt(alpha, rho, phit_)
489  + fvm::div(alphaRhoPhi, phit_)
490  - fvm::laplacian(alpha*rho*DphitEff(), phit_)
491  ==
492  alpha()*rho()*f_()
493  - fvm::SuSp
494  (
495  alpha()*rho()*
496  (
497  G()/k_()
498  - (2.0/3.0)*divU
499  - (2.0*nut*(fvc::grad(phit_) & fvc::grad(k_)))()
500  /(k_()*sigmaPhit_*phit_())
501  )
502  , phit_
503  )
504  + fvOptions(alpha, rho, phit_)
505  );
506 
507  phitEqn.ref().relax();
508  fvOptions.constrain(phitEqn.ref());
509  solve(phitEqn);
510  fvOptions.correct(phit_);
511  bound(phit_, phitMin_);
512 
513  correctNut();
514 }
515 
516 
517 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
518 
519 } // End namespace RASModels
520 } // End namespace Foam
521 
522 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::RASModels::kEpsilonPhitF::correct
virtual void correct()
Solve the transport equations and correct the turbulent viscosity.
Definition: kEpsilonPhitF.C:375
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:287
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:52
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::RASModels::kEpsilonPhitF::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kEpsilonPhitF.H:207
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::RASModels::kEpsilonPhitF::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kEpsilonPhitF.H:208
rho
rho
Definition: readInitialConditions.H:88
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::RASModels::kEpsilonPhitF::Ls
tmp< volScalarField > Ls() const
Return the turbulent length scale, L.
Definition: kEpsilonPhitF.C:74
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:47
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::RASModels::kEpsilonPhitF::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kEpsilonPhitF.H:209
Foam::RASModels::kEpsilonPhitF::correctNut
virtual void correctNut()
Update nut with the latest available k, phit, and T.
Definition: kEpsilonPhitF.C:42
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
divU
zeroField divU
Definition: alphaSuSp.H:3
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::FatalError
error FatalError
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimensioned::getOrAddToDict
static dimensioned< Type > getOrAddToDict(const word &name, dictionary &dict, const dimensionSet &dims=dimless, const Type &deflt=Type(Zero))
Construct dimensioned from dictionary, with default value.
Definition: dimensionedType.C:374
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
kEpsilonPhitF.H
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:55
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::RASModels::kEpsilonPhitF
The k-epsilon-phit-f turbulence closure model for incompressible and compressible flows.
Definition: kEpsilonPhitF.H:131
Foam::RASModels::kEpsilonPhitF::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kEpsilonPhitF.C:348
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::RASModels::kEpsilonPhitF::Ts
tmp< volScalarField > Ts() const
Return the turbulent time scale, T.
Definition: kEpsilonPhitF.C:54
Foam::Switch::getOrAddToDict
static Switch getOrAddToDict(const word &key, dictionary &dict, const Switch deflt=switchType::FALSE)
Definition: Switch.C:164
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106