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 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  ),
123 
124  Cmu_
125  (
127  (
128  "Cmu",
129  this->coeffDict_,
130  0.22
131  )
132  ),
133  Ceps1a_
134  (
136  (
137  "Ceps1a",
138  this->coeffDict_,
139  1.4
140  )
141  ),
142  Ceps1b_
143  (
145  (
146  "Ceps1b",
147  this->coeffDict_,
148  1.0
149  )
150  ),
151  Ceps1c_
152  (
154  (
155  "Ceps1c",
156  this->coeffDict_,
157  0.05
158  )
159  ),
160  Ceps2_
161  (
163  (
164  "Ceps2",
165  this->coeffDict_,
166  1.9
167  )
168  ),
169  Cf1_
170  (
172  (
173  "Cf1",
174  this->coeffDict_,
175  1.4
176  )
177  ),
178  Cf2_
179  (
181  (
182  "Cf2",
183  this->coeffDict_,
184  0.3
185  )
186  ),
187  CL_
188  (
190  (
191  "CL",
192  this->coeffDict_,
193  0.25
194  )
195  ),
196  Ceta_
197  (
199  (
200  "Ceta",
201  this->coeffDict_,
202  110.0
203  )
204  ),
205  CT_
206  (
208  (
209  "CT",
210  this->coeffDict_,
211  6.0
212  )
213  ),
214  sigmaK_
215  (
217  (
218  "sigmaK",
219  this->coeffDict_,
220  1.0
221  )
222  ),
223  sigmaEps_
224  (
226  (
227  "sigmaEps",
228  this->coeffDict_,
229  1.3
230  )
231  ),
232  sigmaPhit_
233  (
235  (
236  "sigmaPhit",
237  this->coeffDict_,
238  1.0
239  )
240  ),
241 
242  k_
243  (
244  IOobject
245  (
246  IOobject::groupName("k", alphaRhoPhi.group()),
247  this->runTime_.timeName(),
248  this->mesh_,
251  ),
252  this->mesh_
253  ),
254  epsilon_
255  (
256  IOobject
257  (
258  IOobject::groupName("epsilon", alphaRhoPhi.group()),
259  this->runTime_.timeName(),
260  this->mesh_,
263  ),
264  this->mesh_
265  ),
266  phit_
267  (
268  IOobject
269  (
270  IOobject::groupName("phit", alphaRhoPhi.group()),
271  this->runTime_.timeName(),
272  this->mesh_,
275  ),
276  this->mesh_
277  ),
278  f_
279  (
280  IOobject
281  (
282  IOobject::groupName("f", alphaRhoPhi.group()),
283  this->runTime_.timeName(),
284  this->mesh_,
287  ),
288  this->mesh_
289  ),
290  T_
291  (
292  IOobject
293  (
294  "T",
295  this->runTime_.timeName(),
296  this->mesh_,
299  false
300  ),
301  this->mesh_,
303  ),
304 
305  phitMin_(dimensionedScalar("phitMin", phit_.dimensions(), SMALL)),
306  fMin_(dimensionedScalar("fMin", f_.dimensions(), SMALL)),
307  TMin_(dimensionedScalar("TMin", dimTime, SMALL)),
308  L2Min_(dimensionedScalar("L2Min", sqr(dimLength), SMALL))
309 {
310  bound(k_, this->kMin_);
311  bound(epsilon_, this->epsilonMin_);
312  bound(phit_, phitMin_);
313  bound(f_, fMin_);
314 
315  if (type == typeName)
316  {
317  this->printCoeffs(type);
318  }
319 
320  if
321  (
322  mag(sigmaK_.value()) < VSMALL
323  || mag(sigmaEps_.value()) < VSMALL
324  || mag(sigmaPhit_.value()) < VSMALL
325  )
326  {
328  << "Non-zero values are required for the model constants:" << nl
329  << "sigmaK = " << sigmaK_ << nl
330  << "sigmaEps = " << sigmaEps_ << nl
331  << "sigmaPhit = " << sigmaPhit_ << nl
332  << exit(FatalError);
333  }
334 }
335 
336 
337 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
338 
339 template<class BasicTurbulenceModel>
341 {
343  {
344  Cmu_.readIfPresent(this->coeffDict());
345  Ceps1a_.readIfPresent(this->coeffDict());
346  Ceps1b_.readIfPresent(this->coeffDict());
347  Ceps1c_.readIfPresent(this->coeffDict());
348  Ceps2_.readIfPresent(this->coeffDict());
349  Cf1_.readIfPresent(this->coeffDict());
350  Cf2_.readIfPresent(this->coeffDict());
351  CL_.readIfPresent(this->coeffDict());
352  Ceta_.readIfPresent(this->coeffDict());
353  CT_.readIfPresent(this->coeffDict());
354  sigmaK_.readIfPresent(this->coeffDict());
355  sigmaEps_.readIfPresent(this->coeffDict());
356  sigmaPhit_.readIfPresent(this->coeffDict());
357 
358  return true;
359  }
360 
361  return false;
362 }
363 
364 
365 template<class BasicTurbulenceModel>
367 {
368  if (!this->turbulence_)
369  {
370  return;
371  }
372 
373  // Construct local convenience references
374  const alphaField& alpha = this->alpha_;
375  const rhoField& rho = this->rho_;
376  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
377  const volVectorField& U = this->U_;
378  volScalarField& nut = this->nut_;
379  fv::options& fvOptions(fv::options::New(this->mesh_));
380 
382 
384  (
385  fvc::div(fvc::absolute(this->phi(), U))
386  );
387 
389  volScalarField G(this->GName(), nut*(2.0*(dev(tS()) && tS())));
390  tS.clear();
391 
392  T_ = Ts();
393  bound(T_, TMin_);
394 
395  const volScalarField L2(type() + "L2", sqr(Ls()) + L2Min_);
396 
397  const volScalarField::Internal Ceps1
398  (
399  "Ceps1",
400  Ceps1a_*(Ceps1b_ + Ceps1c_*sqrt(1.0/phit_()))
401  );
402 
403 
404  // Update epsilon and G at the wall
405  epsilon_.boundaryFieldRef().updateCoeffs();
406 
407  // Turbulent kinetic energy dissipation rate equation (LUU:Eq. 4)
408  // k/T ~ epsilon
409  tmp<fvScalarMatrix> epsEqn
410  (
411  fvm::ddt(alpha, rho, epsilon_)
412  + fvm::div(alphaRhoPhi, epsilon_)
413  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
414  ==
415  alpha()*rho()*Ceps1*G()/T_()
416  - fvm::SuSp
417  (
418  (2.0/3.0*Ceps1)*(alpha()*rho()*divU),
419  epsilon_
420  )
421  - fvm::Sp(alpha()*rho()*Ceps2_/T_(), epsilon_)
422  + fvOptions(alpha, rho, epsilon_)
423  );
424 
425  epsEqn.ref().relax();
426  fvOptions.constrain(epsEqn.ref());
427  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
428  solve(epsEqn);
429  fvOptions.correct(epsilon_);
430  bound(epsilon_, this->epsilonMin_);
431 
432 
433  // Turbulent kinetic energy equation (LUU:Eq. 3)
434  // epsilon/k ~ 1/Ts
436  (
437  fvm::ddt(alpha, rho, k_)
438  + fvm::div(alphaRhoPhi, k_)
439  - fvm::laplacian(alpha*rho*DkEff(), k_)
440  ==
441  alpha()*rho()*G()
442  - fvm::SuSp(2.0/3.0*alpha()*rho()*divU, k_)
443  - fvm::Sp(alpha()*rho()*(1.0/T_()), k_)
444  + fvOptions(alpha, rho, k_)
445  );
446 
447  kEqn.ref().relax();
448  fvOptions.constrain(kEqn.ref());
449  solve(kEqn);
450  fvOptions.correct(k_);
451  bound(k_, this->kMin_);
452 
453 
454  // Elliptic relaxation function equation (LUU:Eq. 18)
455  // All source terms are non-negative functions (LUU:p. 176)
457  (
458  - fvm::laplacian(f_)
459  ==
460  - fvm::Sp(1.0/L2(), f_)
461  - (
462  (Cf1_ - 1.0)*(phit_() - 2.0/3.0)/T_()
463  -(Cf2_*G())/k_()
464  +(Cf2_*(2.0/3.0)*divU)
465  -(2.0*this->nu()*(fvc::grad(phit_) & fvc::grad(k_)))()/k_()
466  -(this->nu()*fvc::laplacian(phit_))()
467  )/L2()
468  );
469 
470  fEqn.ref().relax();
471  solve(fEqn);
472  bound(f_, fMin_);
473 
474 
475  // Normalised wall-normal fluctuating velocity scale equation (LUU:Eq. 17)
476  // All source terms are non-negative functions (LUU:p. 176)
477  tmp<fvScalarMatrix> phitEqn
478  (
479  fvm::ddt(alpha, rho, phit_)
480  + fvm::div(alphaRhoPhi, phit_)
481  - fvm::laplacian(alpha*rho*DphitEff(), phit_)
482  ==
483  alpha()*rho()*f_()
484  - fvm::SuSp
485  (
486  alpha()*rho()*
487  (
488  G()/k_()
489  - (2.0/3.0)*divU
490  - (2.0*nut*(fvc::grad(phit_) & fvc::grad(k_)))()
491  /(k_()*sigmaPhit_*phit_())
492  )
493  , phit_
494  )
495  + fvOptions(alpha, rho, phit_)
496  );
497 
498  phitEqn.ref().relax();
499  fvOptions.constrain(phitEqn.ref());
500  solve(phitEqn);
501  fvOptions.correct(phit_);
502  bound(phit_, phitMin_);
503 
504  correctNut();
505 }
506 
507 
508 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
509 
510 } // End namespace RASModels
511 } // End namespace Foam
512 
513 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
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:366
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:321
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:325
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp< volScalarField >
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:51
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:104
Foam::RASModels::kEpsilonPhitF::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kEpsilonPhitF.H:178
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:179
rho
rho
Definition: readInitialConditions.H:96
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:54
Foam::RASModels::kEpsilonPhitF::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kEpsilonPhitF.H:180
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:258
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:43
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::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModels::kEpsilonPhitFBase
Abstract base-class for the k-epsilon-phit-f model to provide boundary condition access to the phit a...
Definition: kEpsilonPhitFBase.H:59
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:355
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:372
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::IOobject::groupName
static word groupName(StringType name, const word &group)
Create dot-delimited name.group.
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
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
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:103
Foam::RASModels::kEpsilonPhitF::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kEpsilonPhitF.C:340
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::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106