43template<
class BasicTurbulenceModel>
44tmp<volScalarField> kL<BasicTurbulenceModel>::Cmu()
const
47 return (0.556 + 0.108*Rt_)/(1.0 + 0.308*Rt_ + 0.00837*
sqr(Rt_));
51template<
class BasicTurbulenceModel>
52tmp<volScalarField> kL<BasicTurbulenceModel>::CmuPrime()
const
55 return 0.556/(1.0 + 0.277*Rt_);
59template<
class BasicTurbulenceModel>
60tmp<volScalarField> kL<BasicTurbulenceModel>::nutPrime()
const
63 return CmuPrime()*
sqrt(k_)*L_;
67template<
class BasicTurbulenceModel>
68tmp<volScalarField> kL<BasicTurbulenceModel>::epsilonCanopy()
const
71 this->mesh_.template findObject<volScalarField>(
"plantCd");
73 this->mesh_.template findObject<volScalarField>(
"leafAreaDensity");
78 const auto& Cd = *CdPtr;
79 const auto& LAD = *LADPtr;
82 return Cd*LAD*
mag(
U)*k_;
90 this->runTime_.timeName(),
102template<
class BasicTurbulenceModel>
103tmp<volScalarField> kL<BasicTurbulenceModel>::epsilon()
const
106 tmp<volScalarField> tepsilonCanopy = epsilonCanopy();
109 tmp<volScalarField> tepsilonPlain =
pow3(Cmu0_)*
pow(k_, 1.5)/L_;
112 tmp<volScalarField> tepsilon =
max(tepsilonPlain, tepsilonCanopy);
120template<
class BasicTurbulenceModel>
121tmp<volScalarField> kL<BasicTurbulenceModel>::canopyHeight()
const
123 const auto* canopyHeightPtr =
124 this->mesh_.template findObject<volScalarField>(
"canopyHeight");
128 const auto& canopyHeight = *canopyHeightPtr;
137 this->runTime_.timeName(),
149template<
class BasicTurbulenceModel>
150tmp<volScalarField> kL<BasicTurbulenceModel>::L()
const
156 tmp<volScalarField> tLcanopy = kappa_*canopyHeight();
160 return max(Lcanopy, Lplain);
164template<
class BasicTurbulenceModel>
165void kL<BasicTurbulenceModel>::stratification(
const volScalarField& fVB)
167 tmp<volScalarField> tLg =
L();
170 tmp<volScalarField> tcanopyHeight = canopyHeight();
173 tmp<volScalarField> tLcanopy = kappa_*canopyHeight;
176 const scalar Cmu0 = Cmu0_.value();
177 const scalar CbStable = CbStable_.value();
178 const scalar CbUnstable = CbUnstable_.value();
182 if (y_[i] > canopyHeight[i])
187 const scalar Lb = CbStable*
sqrt(k_[i])/
sqrt(fVB[i]);
202 Rt_[i] -
sqr(Rt_[i] + 1.0)/(Rt_[i] - 1.0)
208 *
sqrt(1.0 -
pow(Cmu0, 6.0)*
pow(CbUnstable, -2.0)*Rt_[i]);
224template<
class BasicTurbulenceModel>
227 this->nut_ = Cmu()*
sqrt(k_)*L_;
228 this->nut_.correctBoundaryConditions();
231 BasicTurbulenceModel::correctNut();
235template<
class BasicTurbulenceModel>
248template<
class BasicTurbulenceModel>
257 const word& propertiesName,
343 IOobject::groupName(
"k", alphaRhoPhi.group()),
355 IOobject::groupName(
"L", alphaRhoPhi.group()),
368 IOobject::groupName(
"Rt", alphaRhoPhi.group()),
377 g_(meshObjects::gravity::
New(this->mesh_.time())),
382 if (
type == typeName)
384 this->printCoeffs(
type);
391template<
class BasicTurbulenceModel>
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());
411template<
class BasicTurbulenceModel>
414 if (!this->turbulence_)
477 bound(k_, this->kMin_);
480 Rt_ = fBV*
sqr(k_/tepsilon);
Bound the given scalar field if it has gone unbounded.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Templated abstract base class for RAS turbulence models.
A one-equation (turbulent kinetic energy k) turbulence closure model for incompressible and compressi...
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
volScalarField k_
Turbulent kinetic energy [m2/s2].
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual void correctNut()
Correct the turbulence viscosity.
virtual tmp< fvScalarMatrix > kSource() const
Add explicit source for turbulent kinetic energy.
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
A class for managing temporary objects.
void clear() const noexcept
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
A class for handling words, derived from Foam::string.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
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.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
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.
static constexpr const zero Zero
Global zero (0)
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))
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
#define forAll(list, i)
Loop across all elements in list.
const vector L(dict.get< vector >("L"))