30 const Foam::scalar Foam::atmBoundaryLayer::kappaDefault_ = 0.41;
32 const Foam::scalar Foam::atmBoundaryLayer::CmuDefault_ = 0.09;
37 void Foam::atmBoundaryLayer::init()
39 if (
mag(flowDir_) < small ||
mag(zDir_) < small)
42 <<
"magnitude of n or z must be greater than zero"
47 flowDir_ /=
mag(flowDir_);
50 Ustar_ = kappa_*Uref_/(
log((Zref_ + z0_)/z0_));
86 const scalar epsilonLower
101 epsilonLower_(epsilonLower)
115 kappa_(
dict.lookupOrDefault<scalar>(
"kappa",
dimless, kappaDefault_)),
116 Cmu_(
dict.lookupOrDefault<scalar>(
"Cmu",
dimless, CmuDefault_)),
127 dict.lookupOrDefault<scalar>
145 flowDir_(abl.flowDir_),
151 z0_(mapper(abl.z0_)),
152 zGround_(mapper(abl.zGround_)),
153 Ustar_(mapper(abl.Ustar_)),
154 offset_(abl.offset_),
155 Ulower_(abl.Ulower_),
156 kLower_(abl.kLower_),
157 epsilonLower_(abl.epsilonLower_)
163 flowDir_(abl.flowDir_),
170 zGround_(abl.zGround_),
172 offset_(abl.offset_),
173 Ulower_(abl.Ulower_),
174 kLower_(abl.kLower_),
175 epsilonLower_(abl.epsilonLower_)
187 mapper(z0_, blptf.z0_);
188 mapper(zGround_, blptf.zGround_);
189 mapper(Ustar_, blptf.Ustar_);
195 z0_.reset(blptf.z0_);
196 zGround_.reset(blptf.zGround_);
197 Ustar_.reset(blptf.Ustar_);
209 *
log(
max((zDir_ &
p) - zGround_ + z0_, z0_)/z0_)
214 return flowDir_*Un + flowDir_*Ulower_;
250 pow3(Ustar_)/(kappa_*((zDir_ &
p) - zGround_ + z0_))
256 tepsilon.
ref() =
pos0(z)*tepsilon() +
neg(z)*epsilonLower_;
277 writeEntry(os,
"epsilonLower", epsilonLower_);
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
This class provides functions to evaluate the velocity and turbulence distributions appropriate for a...
void map(const atmBoundaryLayer &, const fieldMapper &)
Map the given atmBoundaryLayer onto this atmBoundaryLayer.
void reset(const atmBoundaryLayer &)
Reset the atmBoundaryLayer to the given atmBoundaryLayer.
void write(Ostream &) const
Write.
tmp< vectorField > U(const vectorField &p) const
Return the velocity distribution for the ATM.
atmBoundaryLayer()
Construct null.
tmp< scalarField > k(const vectorField &p) const
Return the turbulent kinetic energy distribution for the ATM.
tmp< scalarField > epsilon(const vectorField &p) const
Return the turbulent dissipation rate distribution for the ATM.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Abstract base class for field mapping.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimEnergy
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
const dimensionSet dimless
const dimensionSet dimLength
dimensionedScalar log(const dimensionedScalar &ds)
const dimensionSet dimTime
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionSet dimMass
const dimensionSet dimVelocity