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 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
36namespace Foam
37{
38namespace RASModels
39{
40
41// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
42
43template<class BasicTurbulenceModel>
44tmp<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
51template<class BasicTurbulenceModel>
52tmp<volScalarField> kL<BasicTurbulenceModel>::CmuPrime() const
53{
54 // (A:Eq. 32)
55 return 0.556/(1.0 + 0.277*Rt_);
56}
57
58
59template<class BasicTurbulenceModel>
60tmp<volScalarField> kL<BasicTurbulenceModel>::nutPrime() const
61{
62 // (A:Eq. 12)
63 return CmuPrime()*sqrt(k_)*L_;
64}
65
66
67template<class BasicTurbulenceModel>
68tmp<volScalarField> kL<BasicTurbulenceModel>::epsilonCanopy() const
69{
70 const auto* CdPtr =
71 this->mesh_.template findObject<volScalarField>("plantCd");
72 const auto* LADPtr =
73 this->mesh_.template findObject<volScalarField>("leafAreaDensity");
74 const volVectorField& U = this->U_;
75
76 if (CdPtr && LADPtr)
77 {
78 const auto& Cd = *CdPtr;
79 const auto& LAD = *LADPtr;
80
81 // (W:Eq. 13)
82 return Cd*LAD*mag(U)*k_;
83 }
84
86 (
87 IOobject
88 (
89 IOobject::groupName("epsilonCanopy", this->alphaRhoPhi_.group()),
90 this->runTime_.timeName(),
91 this->mesh_,
94 false
95 ),
96 this->mesh_,
98 );
99}
100
101
102template<class BasicTurbulenceModel>
103tmp<volScalarField> kL<BasicTurbulenceModel>::epsilon() const
104{
105 // (W:Eq. 13)
106 tmp<volScalarField> tepsilonCanopy = epsilonCanopy();
107
108 // (A:Eq. 19)
109 tmp<volScalarField> tepsilonPlain = pow3(Cmu0_)*pow(k_, 1.5)/L_;
110
111 // (W:Eq. 13)
112 tmp<volScalarField> tepsilon = max(tepsilonPlain, tepsilonCanopy);
113 volScalarField& epsilon = tepsilon.ref();
114 bound(epsilon, this->epsilonMin_);
115
116 return tepsilon;
117}
118
119
120template<class BasicTurbulenceModel>
121tmp<volScalarField> kL<BasicTurbulenceModel>::canopyHeight() const
122{
123 const auto* canopyHeightPtr =
124 this->mesh_.template findObject<volScalarField>("canopyHeight");
125
126 if (canopyHeightPtr)
127 {
128 const auto& canopyHeight = *canopyHeightPtr;
129 return canopyHeight;
130 }
131
133 (
134 IOobject
135 (
136 IOobject::groupName("canopyHeight", this->alphaRhoPhi_.group()),
137 this->runTime_.timeName(),
138 this->mesh_,
141 false
142 ),
143 this->mesh_,
145 );
146}
147
148
149template<class BasicTurbulenceModel>
150tmp<volScalarField> kL<BasicTurbulenceModel>::L() const
151{
152 // (A:Eq. 22)
153 const volScalarField Lplain(kappa_*y_);
154
155 // Increase roughness for canopy (forest, vegetation etc)
156 tmp<volScalarField> tLcanopy = kappa_*canopyHeight();
157 const volScalarField& Lcanopy = tLcanopy;
158
159 // (W:Eq. 16)
160 return max(Lcanopy, Lplain);
161}
162
163
164template<class BasicTurbulenceModel>
165void kL<BasicTurbulenceModel>::stratification(const volScalarField& fVB)
166{
167 tmp<volScalarField> tLg = L();
168 const volScalarField& Lg = tLg.cref();
169
170 tmp<volScalarField> tcanopyHeight = canopyHeight();
171 const volScalarField& canopyHeight = tcanopyHeight;
172
173 tmp<volScalarField> tLcanopy = kappa_*canopyHeight;
174 const volScalarField& Lcanopy = tLcanopy;
175
176 const scalar Cmu0 = Cmu0_.value();
177 const scalar CbStable = CbStable_.value();
178 const scalar CbUnstable = CbUnstable_.value();
179
180 forAll(L_, i)
181 {
182 if (y_[i] > canopyHeight[i])
183 {
184 if (fVB[i] > 0)
185 {
186 // (A:Eq. 23)
187 const scalar Lb = CbStable*sqrt(k_[i])/sqrt(fVB[i]);
188
189 // (A:Eq. 26)
190 L_[i] = sqrt(sqr(Lg[i]*Lb)/(sqr(Lg[i]) + sqr(Lb)));
191 }
192 else
193 {
194 // For unstable/neutral boundary layer (A:p. 80)
195 // Smoothing function for turbulent Richardson
196 // number to ensure gentle transition into
197 // the regime of strong convection
198 Rt_[i] =
199 min
200 (
201 max(Rt_[i], -1.0),
202 Rt_[i] - sqr(Rt_[i] + 1.0)/(Rt_[i] - 1.0)
203 );
204
205 // (A:Eq. 28)
206 L_[i] =
207 Lg[i]
208 *sqrt(1.0 - pow(Cmu0, 6.0)*pow(CbUnstable, -2.0)*Rt_[i]);
209 }
210 }
211 else
212 {
213 L_[i] = Lcanopy[i];
214 }
215 }
216
217 // Limit characteristic length scale
218 L_ = min(L_, Lmax_);
219}
220
221
222// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
223
224template<class BasicTurbulenceModel>
226{
227 this->nut_ = Cmu()*sqrt(k_)*L_;
228 this->nut_.correctBoundaryConditions();
229 fv::options::New(this->mesh_).correct(this->nut_);
230
231 BasicTurbulenceModel::correctNut();
232}
233
234
235template<class BasicTurbulenceModel>
237{
239 (
240 k_,
241 dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
242 );
243}
244
245
246// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
247
248template<class BasicTurbulenceModel>
250(
251 const alphaField& alpha,
252 const rhoField& rho,
253 const volVectorField& U,
254 const surfaceScalarField& alphaRhoPhi,
255 const surfaceScalarField& phi,
256 const transportModel& transport,
257 const word& propertiesName,
258 const word& type
259)
260:
261 eddyViscosity<RASModel<BasicTurbulenceModel>>
262 (
263 type,
264 alpha,
265 rho,
266 U,
267 alphaRhoPhi,
268 phi,
269 transport,
270 propertiesName
271 ),
272
273 kappa_
274 (
275 dimensioned<scalar>::getOrAddToDict
276 (
277 "kappa",
278 this->coeffDict_,
279 0.41
280 )
281 ),
282 sigmak_
283 (
284 dimensioned<scalar>::getOrAddToDict
285 (
286 "sigmak",
287 this->coeffDict_,
288 1.0
289 )
290 ),
291 beta_
292 (
293 dimensioned<scalar>::getOrAddToDict
294 (
295 "beta",
296 this->coeffDict_,
298 3.3e-03
299 )
300 ),
301 Cmu0_
302 (
303 dimensioned<scalar>::getOrAddToDict
304 (
305 "Cmu0",
306 this->coeffDict_,
307 0.556
308 )
309 ),
310 Lmax_
311 (
312 dimensioned<scalar>::getOrAddToDict
313 (
314 "Lmax",
315 this->coeffDict_,
316 dimLength,
317 GREAT
318 )
319 ),
320 CbStable_
321 (
322 dimensioned<scalar>::getOrAddToDict
323 (
324 "CbStable",
325 this->coeffDict_,
326 0.25
327 )
328 ),
329 CbUnstable_
330 (
331 dimensioned<scalar>::getOrAddToDict
332 (
333 "CbUnstable",
334 this->coeffDict_,
335 0.35
336 )
337 ),
338
339 k_
340 (
342 (
343 IOobject::groupName("k", alphaRhoPhi.group()),
344 this->runTime_.timeName(),
345 this->mesh_,
346 IOobject::MUST_READ,
347 IOobject::AUTO_WRITE
348 ),
349 this->mesh_
350 ),
351 L_
352 (
354 (
355 IOobject::groupName("L", alphaRhoPhi.group()),
356 this->runTime_.timeName(),
357 this->mesh_,
358 IOobject::READ_IF_PRESENT,
359 IOobject::AUTO_WRITE
360 ),
361 this->mesh_,
362 dimensionedScalar(dimLength, scalar(1))
363 ),
364 Rt_
365 (
367 (
368 IOobject::groupName("Rt", alphaRhoPhi.group()),
369 this->runTime_.timeName(),
370 this->mesh_,
371 IOobject::READ_IF_PRESENT,
372 IOobject::AUTO_WRITE
373 ),
374 this->mesh_,
376 ),
377 g_(meshObjects::gravity::New(this->mesh_.time())),
378 y_(wallDist::New(this->mesh_).y())
379{
380 bound(k_, this->kMin_);
381
382 if (type == typeName)
383 {
384 this->printCoeffs(type);
385 }
386}
387
388
389// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
390
391template<class BasicTurbulenceModel>
393{
395 {
396 kappa_.readIfPresent(this->coeffDict());
397 sigmak_.readIfPresent(this->coeffDict());
398 beta_.readIfPresent(this->coeffDict());
399 Cmu0_.readIfPresent(this->coeffDict());
400 Lmax_.readIfPresent(this->coeffDict());
401 CbStable_.readIfPresent(this->coeffDict());
402 CbUnstable_.readIfPresent(this->coeffDict());
403
404 return true;
405 }
406
407 return false;
408}
409
410
411template<class BasicTurbulenceModel>
413{
414 if (!this->turbulence_)
415 {
416 return;
417 }
418
419 // Construct local convenience references
420 const alphaField& alpha = this->alpha_;
421 const rhoField& rho = this->rho_;
422 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
423 const volVectorField& U = this->U_;
424 const volScalarField& nut = this->nut_;
425
427
429
430 // Turbulent kinetic energy production rate
433 (
434 this->GName(),
435 nut.v()*2*magSqr(dev(symm(tgradU.cref().v())))
436 );
437 tgradU.clear();
438
439 // Square of Brunt-Vaisala (buoyancy) frequency
440 const auto& T = U.mesh().lookupObject<volScalarField>("T");
441 tmp<volScalarField> tfBV = -beta_*(fvc::grad(T) & g_);
442 const volScalarField& fBV = tfBV.cref();
443
444 // Sink or source of TKE depending on stratification type (A:Eq. 15)
445 tmp<volScalarField> tPb = -fBV*nutPrime();
446 const volScalarField& Pb = tPb.cref();
447
448 // Turbulent kinetic energy dissipation rate due to plains and canopy
449 tmp<volScalarField> tepsilon = epsilon();
450 const volScalarField& epsilon = tepsilon.cref();
451
452 // Divergence of velocity
453 tmp<volScalarField> tdivU = fvc::div(fvc::absolute(this->phi(), U));
454 const volScalarField::Internal& divU = tdivU.cref().v();
455
456 // Turbulent kinetic energy equation
458 (
459 fvm::ddt(alpha, rho, k_)
460 + fvm::div(alphaRhoPhi, k_)
461 - fvm::laplacian(alpha*rho*DkEff(), k_)
462 ==
463 alpha()*rho()*G
464 + fvm::SuSp((Pb - epsilon)/k_, k_)
465 - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
466 + kSource()
467 + fvOptions(alpha, rho, k_)
468 );
469
470 tdivU.clear();
471 tPb.clear();
472
473 kEqn.ref().relax();
474 fvOptions.constrain(kEqn.ref());
475 solve(kEqn);
476 fvOptions.correct(k_);
477 bound(k_, this->kMin_);
478
479 // Turbulent Richardson number (A:Eq. 29)
480 Rt_ = fBV*sqr(k_/tepsilon);
481
482 stratification(fBV);
483 tfBV.clear();
484
485 correctNut();
486}
487
488
489// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
490
491} // End namespace RASModels
492} // End namespace Foam
493
494// ************************************************************************* //
scalar y
Bound the given scalar field if it has gone unbounded.
fv::options & fvOptions
surfaceScalarField & phi
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:55
A one-equation (turbulent kinetic energy k) turbulence closure model for incompressible and compressi...
Definition: kL.H:225
BasicTurbulenceModel::alphaField alphaField
Definition: kL.H:324
BasicTurbulenceModel::rhoField rhoField
Definition: kL.H:325
volScalarField k_
Turbulent kinetic energy [m2/s2].
Definition: kL.H:296
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kL.C:412
virtual void correctNut()
Correct the turbulence viscosity.
Definition: kL.C:225
virtual tmp< fvScalarMatrix > kSource() const
Add explicit source for turbulent kinetic energy.
Definition: kL.C:236
BasicTurbulenceModel::transportModel transportModel
Definition: kL.H:326
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kL.C:392
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:58
Finite-volume options.
Definition: fvOptions.H:59
A class for managing temporary objects.
Definition: tmp.H:65
const T & cref() const
Definition: tmpI.H:213
void clear() const noexcept
Definition: tmpI.H:287
T & ref() const
Definition: tmpI.H:227
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
Definition: wallDist.H:78
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
thermo correct()
const volScalarField & T
scalar epsilon
scalar nut
zeroField divU
Definition: alphaSuSp.H:3
word timeName
Definition: getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
volScalarField & alpha
CEqn solve()
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
const vector L(dict.get< vector >("L"))