Go to the documentation of this file.
67 const modelType& filmModel =
68 db().time().lookupObject<modelType>(filmRegionName_);
70 const label filmPatchI = filmModel.regionPatchID(
patchi);
73 scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchI];
74 filmModel.toPrimary(filmPatchI, mDotFilmp);
93 const scalar Cmu25 =
pow(Cmu_, 0.25);
97 label faceCellI = patch().faceCells()[faceI];
99 scalar ut = Cmu25*
sqrt(
k[faceCellI]);
101 scalar
yPlus =
y[faceI]*ut/nuw[faceI];
103 scalar mStar = mDotFilmp[faceI]/(
y[faceI]*ut);
106 if (
yPlus > yPlusCrit_)
108 scalar expTerm =
exp(
min(50.0, B_*mStar));
109 scalar powTerm =
pow(
yPlus, mStar/kappa_);
110 factor = mStar/(expTerm*powTerm - 1.0 + ROOTVSMALL);
114 scalar expTerm =
exp(
min(50.0, mStar));
115 factor = mStar/(expTerm*
yPlus - 1.0 + ROOTVSMALL);
118 uTau[faceI] =
sqrt(
max(0, magGradU[faceI]*ut*factor));
146 sqr(
calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
160 filmRegionName_(
"surfaceFilmProperties"),
191 dict.lookupOrDefault<
word>(
"filmRegion",
"surfaceFilmProperties")
193 B_(
dict.lookupOrDefault(
"B", 5.5)),
194 yPlusCrit_(
dict.lookupOrDefault(
"yPlusCrit", 11.05))
250 writeEntryIfDifferent<word>
254 "surfaceFilmProperties",
260 writeEntry(
"value", os);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
virtual void write(Ostream &os) const
Write.
A class for handling words, derived from string.
scalar yPlusCrit_
y+ value for laminar -> turbulent transition (default = 11.05)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
scalar B_
B Coefficient (default = 5.5)
rDeltaT dimensionedInternalField()
makePatchTypeField(fvPatchScalarField, alphatFilmWallFunctionFvPatchScalarField)
static const word propertiesName
Default name of the turbulence properties dictionary.
const Time & time() const
Return the reference to the time database.
const char *const group
Group name for atomic constants.
nutkFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
word filmRegionName_
Name of film region.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Macros for easy insertion into run-time selection tables.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
bool foundObject(const word &name) const
Is the named Type found?
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
label k
Boltzmann constant.
This boundary condition provides a turbulent viscosity condition when using wall functions,...
Foam::fvPatchFieldMapper.
Base class for surface film models.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...