Go to the documentation of this file.
64 initABL_(
dict.getOrDefault<
bool>(
"initABL", true)),
70 C1_(
dict.getOrDefault(
"C1", 0.0)),
71 C2_(
dict.getOrDefault(
"C2", 1.0)),
91 initABL_(abl.initABL_),
99 flowDir_(abl.flowDir_.clone()),
100 zDir_(abl.zDir_.clone()),
101 Uref_(abl.Uref_.clone()),
102 Zref_(abl.Zref_.clone()),
103 z0_(abl.z0_.clone(patch_)),
104 d_(abl.d_.clone(patch_))
110 initABL_(abl.initABL_),
118 flowDir_(abl.flowDir_),
122 z0_(abl.z0_.clone(patch_)),
123 d_(abl.d_.clone(patch_))
132 const vector dir(flowDir_->value(t));
133 const scalar magDir =
mag(dir);
138 <<
"magnitude of " << flowDir_->name() <<
" = " << magDir
139 <<
" vector must be greater than zero"
150 const vector dir(zDir_->value(t));
151 const scalar magDir =
mag(dir);
156 <<
"magnitude of " << zDir_->name() <<
" = " << magDir
157 <<
" vector must be greater than zero"
168 const scalar Uref = Uref_->value(t);
169 const scalar Zref = Zref_->value(t);
174 <<
"Negative entry in " << Zref_->name() <<
" = " << Zref
179 return kappa_*Uref/
log((Zref + z0)/z0);
187 z0_->autoMap(mapper);
198 const atmBoundaryLayer& abl,
204 z0_->rmap(abl.z0_(), addr);
208 d_->rmap(abl.d_(), addr);
218 const scalar groundMin =
zDir() & ppMin_;
235 const scalar groundMin =
zDir() & ppMin_;
240 *
sqrt(C1_*
log(((
zDir() & pCf) - groundMin - d + z0)/z0) + C2_);
249 const scalar groundMin =
zDir() & ppMin_;
253 pow3(
Ustar(z0))/(kappa_*((
zDir() & pCf) - groundMin - d + z0))
254 *
sqrt(C1_*
log(((
zDir() & pCf) - groundMin - d + z0)/z0) + C2_);
263 const scalar groundMin =
zDir() & ppMin_;
266 return Ustar(z0)/(kappa_*
sqrt(Cmu_)*((
zDir() & pCf) - groundMin - d + z0));
273 os.writeEntry(
"kappa", kappa_);
274 os.writeEntry(
"Cmu", Cmu_);
275 os.writeEntry(
"C1", C1_);
276 os.writeEntry(
"C2", C2_);
277 flowDir_->writeData(
os);
278 zDir_->writeData(
os);
279 Uref_->writeData(
os);
280 Zref_->writeData(
os);
const Field< point_type > & points() const noexcept
List< label > labelList
A List of labels.
scalar timeOutputValue() const
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< scalarField > Ustar(const scalarField &z0) const
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for managing temporary objects.
Base class to set log-law type ground-normal inlet boundary conditions for wind velocity and turbulen...
label min(const labelHashSet &set, label minValue=labelMax)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Field< vector > vectorField
Specialisation of Field<T> for vector.
tmp< scalarField > k(const vectorField &pCf) const
tmp< vectorField > U(const vectorField &pCf) const
void write(Ostream &) const
Generic templated field type.
dimensionedScalar pow3(const dimensionedScalar &ds)
A patch is a list of labels that address the faces in the global face list.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
label max(const labelHashSet &set, label maxValue=labelMin)
tmp< scalarField > omega(const vectorField &pCf) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
dimensionedScalar log(const dimensionedScalar &ds)
OBJstream os(runTime.globalPath()/outputName)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
errorManip< error > abort(error &err)
Vector< scalar > vector
A scalar version of the templated Vector.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
#define FatalErrorInFunction
atmBoundaryLayer(const Time &time, const polyPatch &pp)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< scalarField > epsilon(const vectorField &pCf) const
static MinMax< T > ge(const T &minVal)
void autoMap(const fvPatchFieldMapper &)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void rmap(const atmBoundaryLayer &, const labelList &)
A bounding box defined in terms of min/max extrema points.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
A min/max value pair with additional methods. In addition to conveniently storing values,...