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
70 const auto* CdPtr = this->mesh_.template findObject<volScalarField>(
"Cd");
71 const auto* LADPtr = this->mesh_.template findObject<volScalarField>(
"LAD");
76 const auto& Cd = *CdPtr;
77 const auto& LAD = *LADPtr;
80 return Cd*LAD*
mag(
U)*k_;
93template<
class BasicTurbulenceModel>
94tmp<volScalarField> kL<BasicTurbulenceModel>::epsilon()
const
97 tmp<volScalarField> tepsilonCanopy = epsilonCanopy();
100 tmp<volScalarField> tepsilonPlain =
pow3(Cmu0_)*
pow(k_, 1.5)/L_;
103 tmp<volScalarField> tepsilon =
max(tepsilonPlain, tepsilonCanopy);
111template<
class BasicTurbulenceModel>
112tmp<volScalarField> kL<BasicTurbulenceModel>::canopyHeight()
const
114 tmp<volScalarField> tcanopy;
118 this->mesh_.template cfindObject<volScalarField>(
"canopyHeight")
136template<
class BasicTurbulenceModel>
137tmp<volScalarField> kL<BasicTurbulenceModel>::L()
const
143 tmp<volScalarField> tLcanopy = kappa_*canopyHeight();
147 return max(Lcanopy, Lplain);
151template<
class BasicTurbulenceModel>
152void kL<BasicTurbulenceModel>::stratification(
const volScalarField& fVB)
154 tmp<volScalarField> tLg =
L();
157 tmp<volScalarField> tcanopyHeight = canopyHeight();
160 tmp<volScalarField> tLcanopy = kappa_*canopyHeight;
163 const scalar Cmu0 = Cmu0_.value();
164 const scalar CbStable = CbStable_.value();
165 const scalar CbUnstable = CbUnstable_.value();
169 if (y_[i] > canopyHeight[i])
174 const scalar Lb = CbStable*
sqrt(k_[i])/
sqrt(fVB[i]);
189 Rt_[i] -
sqr(Rt_[i] + 1.0)/(Rt_[i] - 1.0)
195 *
sqrt(1.0 -
pow(Cmu0, 6.0)*
pow(CbUnstable, -2.0)*Rt_[i]);
211template<
class BasicTurbulenceModel>
214 this->nut_ = Cmu()*
sqrt(k_)*L_;
215 this->nut_.correctBoundaryConditions();
218 BasicTurbulenceModel::correctNut();
222template<
class BasicTurbulenceModel>
235template<
class BasicTurbulenceModel>
236kL<BasicTurbulenceModel>::kL
238 const alphaField&
alpha,
244 const word& propertiesName,
330 IOobject::groupName(
"k", alphaRhoPhi.group()),
342 IOobject::groupName(
"L", alphaRhoPhi.group()),
355 IOobject::groupName(
"Rt", alphaRhoPhi.group()),
371 this->printCoeffs(
type);
378template<
class BasicTurbulenceModel>
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());
398template<
class BasicTurbulenceModel>
401 if (!this->turbulence_)
407 const alphaField&
alpha = this->alpha_;
408 const rhoField&
rho = this->rho_;
415 eddyViscosity<RASModel<BasicTurbulenceModel>>
::correct();
428 tmp<volScalarField> tfBV = -beta_*(
fvc::grad(
T) & g_);
432 tmp<volScalarField> tPb = -fBV*nutPrime();
436 tmp<volScalarField> tepsilon =
epsilon();
444 tmp<fvScalarMatrix> kEqn
464 bound(k_, this->kMin_);
467 Rt_ = fBV*
sqr(k_/tepsilon);
Bound the given scalar field if it has gone unbounded.
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())
DimensionedField< scalar, volMesh > Internal
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
@ NO_REGISTER
Do not request registration (bool: false).
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.
BasicTurbulenceModel::alphaField alphaField
dimensionedScalar sigmak_
Empirical model coefficient.
BasicTurbulenceModel::rhoField rhoField
volScalarField k_
Turbulent kinetic energy [m2/s2].
volScalarField L_
Characteristic length scale [m].
const uniformDimensionedVectorField & g_
Gravitational acceleration [m2/s2].
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedScalar kappa_
von Karman constant
dimensionedScalar CbUnstable_
Unstable stratification constant.
dimensionedScalar Cmu0_
Empirical model coefficient.
volScalarField Rt_
Turbulent Richardson number [-].
const volScalarField & y_
Wall distance.
dimensionedScalar CbStable_
Stable stratification constant.
virtual void correctNut()
Correct the turbulence viscosity.
virtual tmp< fvScalarMatrix > kSource() const
Add explicit source for turbulent kinetic energy.
BasicTurbulenceModel::transportModel transportModel
dimensionedScalar Lmax_
Maximum mixing-length scalar [m].
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
dimensionedScalar beta_
Thermal expansion coefficient [1/K].
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
eddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
virtual tmp< volScalarField > nut() const
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
A class for managing temporary objects.
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Base-class for all transport models used by the incompressible turbulence models.
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.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
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.
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
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
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).
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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"))