Go to the documentation of this file.
53 fixedValueFvPatchScalarField(
p, iF),
54 filmRegionName_(
"surfaceFilmProperties"),
72 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
90 fixedValueFvPatchScalarField(
p, iF,
dict),
93 dict.lookupOrDefault<
word>(
"filmRegion",
"surfaceFilmProperties")
95 B_(
dict.lookupOrDefault(
"B", 5.5)),
96 yPlusCrit_(
dict.lookupOrDefault(
"yPlusCrit", 11.05)),
97 Cmu_(
dict.lookupOrDefault(
"Cmu", 0.09)),
98 kappa_(
dict.lookupOrDefault(
"kappa", 0.41)),
99 Prt_(
dict.lookupOrDefault(
"Prt", 0.85))
109 fixedValueFvPatchScalarField(fwfpsf),
126 fixedValueFvPatchScalarField(fwfpsf, iF),
163 const modelType& filmModel =
166 const label filmPatchI = filmModel.regionPatchID(
patchi);
169 scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchI];
170 filmModel.toPrimary(filmPatchI, mDotFilmp);
191 const scalar Cmu25 =
pow(
Cmu_, 0.25);
197 label faceCellI = patch().faceCells()[faceI];
201 scalar
yPlus =
y[faceI]*
uTau/(muw[faceI]/rhow[faceI]);
203 scalar
Pr = muw[faceI]/alphaw[faceI];
206 scalar mStar = mDotFilmp[faceI]/(
y[faceI]*
uTau);
213 mStar/(expTerm*(
pow(yPlusRatio, powTerm)) - 1.0 + ROOTVSMALL);
218 factor = mStar/(expTerm - 1.0 + ROOTVSMALL);
221 scalar dx = patch().deltaCoeffs()[faceI];
231 fixedValueFvPatchScalarField::updateCoeffs();
240 writeEntryIfDifferent<word>
244 "surfaceFilmProperties",
252 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.
static word groupName(Name name, const word &group)
scalar yPlusCrit_
y+ value for laminar -> turbulent transition (default = 11.05)
A class for handling words, derived from string.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
rDeltaT dimensionedInternalField()
makePatchTypeField(fvPatchScalarField, alphatFilmWallFunctionFvPatchScalarField)
static const word propertiesName
Default name of the turbulence properties dictionary.
const char *const group
Group name for atomic constants.
dimensionedScalar exp(const dimensionedScalar &ds)
scalar kappa_
Von-Karman constant (default = 0.41)
alphatFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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.
This boundary condition provides a turbulent thermal diffusivity condition when using wall functions,...
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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.
scalar Prt_
Turbulent Prandtl number (default = 0.85)
const volScalarField & alphaEff
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
scalar B_
B Coefficient (default = 5.5)
static int & msgType()
Message tag of standard messages.
virtual void write(Ostream &) const
Write.
word filmRegionName_
Name of film region.
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.
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 tmp< volScalarField > alpha() const
Return the laminar thermal diffusivity for enthalpy [kg/m/s].
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
scalar Cmu_
Turbulent Cmu coefficient (default = 0.09)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.