kL.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) 2021-2023 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 "kL.H"
29 #include "fvOptions.H"
30 #include "bound.H"
31 #include "gravityMeshObject.H"
32 #include "wallDist.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace RASModels
39 {
40 
41 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
42 
43 template<class BasicTurbulenceModel>
44 tmp<volScalarField> kL<BasicTurbulenceModel>::Cmu() const
45 {
46  // (A:Eq. 31)
47  return (0.556 + 0.108*Rt_)/(1.0 + 0.308*Rt_ + 0.00837*sqr(Rt_));
48 }
49 
50 
51 template<class BasicTurbulenceModel>
52 tmp<volScalarField> kL<BasicTurbulenceModel>::CmuPrime() const
53 {
54  // (A:Eq. 32)
55  return 0.556/(1.0 + 0.277*Rt_);
56 }
57 
58 
59 template<class BasicTurbulenceModel>
60 tmp<volScalarField> kL<BasicTurbulenceModel>::nutPrime() const
61 {
62  // (A:Eq. 12)
63  return CmuPrime()*sqrt(k_)*L_;
64 }
65 
66 
67 template<class BasicTurbulenceModel>
68 tmp<volScalarField> kL<BasicTurbulenceModel>::epsilonCanopy() const
69 {
70  const auto* CdPtr = this->mesh_.template findObject<volScalarField>("Cd");
71  const auto* LADPtr = this->mesh_.template findObject<volScalarField>("LAD");
72  const volVectorField& U = this->U_;
73 
74  if (CdPtr && LADPtr)
75  {
76  const auto& Cd = *CdPtr;
77  const auto& LAD = *LADPtr;
78 
79  // (W:Eq. 13)
80  return Cd*LAD*mag(U)*k_;
81  }
82 
83  return volScalarField::New
84  (
85  IOobject::groupName("epsilonCanopy", this->alphaRhoPhi_.group()),
87  this->mesh_,
89  );
90 }
91 
92 
93 template<class BasicTurbulenceModel>
94 tmp<volScalarField> kL<BasicTurbulenceModel>::epsilon() const
95 {
96  // (W:Eq. 13)
97  tmp<volScalarField> tepsilonCanopy = epsilonCanopy();
98 
99  // (A:Eq. 19)
100  tmp<volScalarField> tepsilonPlain = pow3(Cmu0_)*pow(k_, 1.5)/L_;
101 
102  // (W:Eq. 13)
103  tmp<volScalarField> tepsilon = max(tepsilonPlain, tepsilonCanopy);
104  volScalarField& epsilon = tepsilon.ref();
105  bound(epsilon, this->epsilonMin_);
106 
107  return tepsilon;
108 }
109 
110 
111 template<class BasicTurbulenceModel>
112 tmp<volScalarField> kL<BasicTurbulenceModel>::canopyHeight() const
113 {
114  tmp<volScalarField> tcanopy;
115 
116  tcanopy.cref
117  (
118  this->mesh_.template cfindObject<volScalarField>("canopyHeight")
119  );
120 
121  if (tcanopy)
122  {
123  return tcanopy;
124  }
125 
126  return volScalarField::New
127  (
128  IOobject::groupName("canopyHeight", this->alphaRhoPhi_.group()),
130  this->mesh_,
132  );
133 }
134 
135 
136 template<class BasicTurbulenceModel>
137 tmp<volScalarField> kL<BasicTurbulenceModel>::L() const
138 {
139  // (A:Eq. 22)
140  const volScalarField Lplain(kappa_*y_);
141 
142  // Increase roughness for canopy (forest, vegetation etc)
143  tmp<volScalarField> tLcanopy = kappa_*canopyHeight();
144  const volScalarField& Lcanopy = tLcanopy;
145 
146  // (W:Eq. 16)
147  return max(Lcanopy, Lplain);
148 }
149 
150 
151 template<class BasicTurbulenceModel>
152 void kL<BasicTurbulenceModel>::stratification(const volScalarField& fVB)
153 {
154  tmp<volScalarField> tLg = L();
155  const volScalarField& Lg = tLg.cref();
156 
157  tmp<volScalarField> tcanopyHeight = canopyHeight();
158  const volScalarField& canopyHeight = tcanopyHeight;
159 
160  tmp<volScalarField> tLcanopy = kappa_*canopyHeight;
161  const volScalarField& Lcanopy = tLcanopy;
162 
163  const scalar Cmu0 = Cmu0_.value();
164  const scalar CbStable = CbStable_.value();
165  const scalar CbUnstable = CbUnstable_.value();
166 
167  forAll(L_, i)
168  {
169  if (y_[i] > canopyHeight[i])
170  {
171  if (fVB[i] > 0)
172  {
173  // (A:Eq. 23)
174  const scalar Lb = CbStable*sqrt(k_[i])/sqrt(fVB[i]);
175 
176  // (A:Eq. 26)
177  L_[i] = sqrt(sqr(Lg[i]*Lb)/(sqr(Lg[i]) + sqr(Lb)));
178  }
179  else
180  {
181  // For unstable/neutral boundary layer (A:p. 80)
182  // Smoothing function for turbulent Richardson
183  // number to ensure gentle transition into
184  // the regime of strong convection
185  Rt_[i] =
186  min
187  (
188  max(Rt_[i], -1.0),
189  Rt_[i] - sqr(Rt_[i] + 1.0)/(Rt_[i] - 1.0)
190  );
191 
192  // (A:Eq. 28)
193  L_[i] =
194  Lg[i]
195  *sqrt(1.0 - pow(Cmu0, 6.0)*pow(CbUnstable, -2.0)*Rt_[i]);
196  }
197  }
198  else
199  {
200  L_[i] = Lcanopy[i];
201  }
202  }
203 
204  // Limit characteristic length scale
205  L_ = min(L_, Lmax_);
206 }
207 
208 
209 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
210 
211 template<class BasicTurbulenceModel>
213 {
214  this->nut_ = Cmu()*sqrt(k_)*L_;
215  this->nut_.correctBoundaryConditions();
216  fv::options::New(this->mesh_).correct(this->nut_);
217 
218  BasicTurbulenceModel::correctNut();
219 }
220 
221 
222 template<class BasicTurbulenceModel>
224 {
226  (
227  k_,
228  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
229  );
230 }
231 
232 
233 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
234 
235 template<class BasicTurbulenceModel>
237 (
238  const alphaField& alpha,
239  const rhoField& rho,
240  const volVectorField& U,
241  const surfaceScalarField& alphaRhoPhi,
242  const surfaceScalarField& phi,
243  const transportModel& transport,
244  const word& propertiesName,
245  const word& type
246 )
247 :
248  eddyViscosity<RASModel<BasicTurbulenceModel>>
249  (
250  type,
251  alpha,
252  rho,
253  U,
254  alphaRhoPhi,
255  phi,
256  transport,
257  propertiesName
258  ),
259 
260  kappa_
261  (
262  dimensioned<scalar>::getOrAddToDict
263  (
264  "kappa",
265  this->coeffDict_,
266  0.41
267  )
268  ),
269  sigmak_
270  (
271  dimensioned<scalar>::getOrAddToDict
272  (
273  "sigmak",
274  this->coeffDict_,
275  1.0
276  )
277  ),
278  beta_
279  (
280  dimensioned<scalar>::getOrAddToDict
281  (
282  "beta",
283  this->coeffDict_,
285  3.3e-03
286  )
287  ),
288  Cmu0_
289  (
290  dimensioned<scalar>::getOrAddToDict
291  (
292  "Cmu0",
293  this->coeffDict_,
294  0.556
295  )
296  ),
297  Lmax_
298  (
299  dimensioned<scalar>::getOrAddToDict
300  (
301  "Lmax",
302  this->coeffDict_,
303  dimLength,
304  GREAT
305  )
306  ),
307  CbStable_
308  (
309  dimensioned<scalar>::getOrAddToDict
310  (
311  "CbStable",
312  this->coeffDict_,
313  0.25
314  )
315  ),
316  CbUnstable_
317  (
318  dimensioned<scalar>::getOrAddToDict
319  (
320  "CbUnstable",
321  this->coeffDict_,
322  0.35
323  )
324  ),
325 
326  k_
327  (
328  IOobject
329  (
330  IOobject::groupName("k", alphaRhoPhi.group()),
331  this->runTime_.timeName(),
332  this->mesh_,
333  IOobject::MUST_READ,
334  IOobject::AUTO_WRITE
335  ),
336  this->mesh_
337  ),
338  L_
339  (
340  IOobject
341  (
342  IOobject::groupName("L", alphaRhoPhi.group()),
343  this->runTime_.timeName(),
344  this->mesh_,
345  IOobject::READ_IF_PRESENT,
346  IOobject::AUTO_WRITE
347  ),
348  this->mesh_,
349  dimensionedScalar(dimLength, scalar(1))
350  ),
351  Rt_
352  (
353  IOobject
354  (
355  IOobject::groupName("Rt", alphaRhoPhi.group()),
356  this->runTime_.timeName(),
357  this->mesh_,
358  IOobject::READ_IF_PRESENT,
359  IOobject::AUTO_WRITE
360  ),
361  this->mesh_,
363  ),
364  g_(meshObjects::gravity::New(this->mesh_.time())),
365  y_(wallDist::New(this->mesh_).y())
366 {
367  bound(k_, this->kMin_);
368 
369  if (type == typeName)
370  {
371  this->printCoeffs(type);
372  }
373 }
374 
375 
376 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
377 
378 template<class BasicTurbulenceModel>
380 {
382  {
383  kappa_.readIfPresent(this->coeffDict());
384  sigmak_.readIfPresent(this->coeffDict());
385  beta_.readIfPresent(this->coeffDict());
386  Cmu0_.readIfPresent(this->coeffDict());
387  Lmax_.readIfPresent(this->coeffDict());
388  CbStable_.readIfPresent(this->coeffDict());
389  CbUnstable_.readIfPresent(this->coeffDict());
390 
391  return true;
392  }
393 
394  return false;
395 }
396 
397 
398 template<class BasicTurbulenceModel>
400 {
401  if (!this->turbulence_)
402  {
403  return;
404  }
405 
406  // Construct local convenience references
407  const alphaField& alpha = this->alpha_;
408  const rhoField& rho = this->rho_;
409  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
410  const volVectorField& U = this->U_;
411  const volScalarField& nut = this->nut_;
412 
413  fv::options& fvOptions(fv::options::New(this->mesh_));
414 
415  eddyViscosity<RASModel<BasicTurbulenceModel>>::correct();
416 
417  // Turbulent kinetic energy production rate
418  tmp<volTensorField> tgradU = fvc::grad(U);
420  (
421  this->GName(),
422  nut.v()*2*magSqr(devSymm(tgradU.cref().v()))
423  );
424  tgradU.clear();
425 
426  // Square of Brunt-Vaisala (buoyancy) frequency
427  const auto& T = U.mesh().lookupObject<volScalarField>("T");
428  tmp<volScalarField> tfBV = -beta_*(fvc::grad(T) & g_);
429  const volScalarField& fBV = tfBV.cref();
430 
431  // Sink or source of TKE depending on stratification type (A:Eq. 15)
432  tmp<volScalarField> tPb = -fBV*nutPrime();
433  const volScalarField& Pb = tPb.cref();
434 
435  // Turbulent kinetic energy dissipation rate due to plains and canopy
436  tmp<volScalarField> tepsilon = epsilon();
437  const volScalarField& epsilon = tepsilon.cref();
438 
439  // Divergence of velocity
440  tmp<volScalarField> tdivU = fvc::div(fvc::absolute(this->phi(), U));
441  const volScalarField::Internal& divU = tdivU.cref().v();
442 
443  // Turbulent kinetic energy equation
444  tmp<fvScalarMatrix> kEqn
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((Pb - epsilon)/k_, k_)
452  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
453  + kSource()
454  + fvOptions(alpha, rho, k_)
455  );
456 
457  tdivU.clear();
458  tPb.clear();
459 
460  kEqn.ref().relax();
461  fvOptions.constrain(kEqn.ref());
462  solve(kEqn);
463  fvOptions.correct(k_);
464  bound(k_, this->kMin_);
465 
466  // Turbulent Richardson number (A:Eq. 29)
467  Rt_ = fBV*sqr(k_/tepsilon);
468 
469  stratification(fBV);
470  tfBV.clear();
471 
472  correctNut();
473 }
474 
475 
476 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
477 
478 } // End namespace RASModels
479 } // End namespace Foam
480 
481 // ************************************************************************* //
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensionedScalar kMin_
Lower limit of k.
Definition: RASModel.H:109
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const vector L(dict.get< vector >("L"))
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)
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
Generic dimensioned Type class.
const dimensionSet dimless
Dimensionless.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:51
virtual tmp< fvScalarMatrix > kSource() const
Add explicit source for turbulent kinetic energy.
Definition: kL.C:216
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
fv::options & fvOptions
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:77
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
volScalarField k_
Turbulent kinetic energy [m2/s2].
Definition: kL.H:328
word timeName
Definition: getTimeIndex.H:3
A one-equation (turbulent kinetic energy k) turbulence closure model for incompressible and compressi...
Definition: kL.H:218
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kL.C:372
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
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
virtual void correctNut()
Correct the turbulence viscosity.
Definition: kL.C:205
A class for handling words, derived from Foam::string.
Definition: word.H:63
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
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.
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:27
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
Bound the given scalar field if it has gone unbounded.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2406/OpenFOAM-v2406/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
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
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
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. ...
scalar epsilon
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields...
Definition: wallDist.H:71
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kL.C:392
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
Do not request registration (bool: false)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
scalar nut
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127