Go to the documentation of this file.
56 fixedValueFvPatchScalarField(
p, iF),
57 filmRegionName_(
"surfaceFilmProperties"),
75 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
93 fixedValueFvPatchScalarField(
p, iF,
dict),
96 dict.getOrDefault<
word>(
"filmRegion",
"surfaceFilmProperties")
98 B_(
dict.getOrDefault(
"B", 5.5)),
99 yPlusCrit_(
dict.getOrDefault(
"yPlusCrit", 11.05)),
100 Cmu_(
dict.getOrDefault(
"Cmu", 0.09)),
101 kappa_(
dict.getOrDefault(
"kappa", 0.41)),
102 Prt_(
dict.getOrDefault(
"Prt", 0.85))
112 fixedValueFvPatchScalarField(fwfpsf),
129 fixedValueFvPatchScalarField(fwfpsf, iF),
148 const auto* filmModelPtr = db().time().findObject
158 const auto& filmModel = *filmModelPtr;
166 const label patchi =
patch().index();
169 const label filmPatchi = filmModel.regionPatchID(patchi);
172 scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchi];
173 filmModel.toPrimary(filmPatchi, mDotFilmp);
181 internalField().
group()
186 const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
194 const scalar Cmu25 =
pow(
Cmu_, 0.25);
200 label faceCelli =
patch().faceCells()[facei];
204 scalar
yPlus =
y[facei]*
uTau/(muw[facei]/rhow[facei]);
206 scalar
Pr = muw[facei]/alphaw[facei];
209 scalar mStar = mDotFilmp[facei]/(
y[facei]*
uTau);
216 mStar/(expTerm*(
pow(yPlusRatio, powTerm)) - 1.0 + ROOTVSMALL);
221 factor = mStar/(expTerm - 1.0 + ROOTVSMALL);
224 scalar dx =
patch().deltaCoeffs()[facei];
228 alphat[facei] =
max(
alphaEff - alphaw[facei], 0.0);
234 fixedValueFvPatchScalarField::updateCoeffs();
246 "surfaceFilmProperties",
254 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.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
scalar yPlusCrit_
y+ value for laminar -> turbulent transition (default = 11.05)
A class for handling words, derived from Foam::string.
constexpr const char *const group
Group name for atomic constants.
makePatchTypeField(fvPatchScalarField, alphatFilmWallFunctionFvPatchScalarField)
static const word propertiesName
Default name of the turbulence properties dictionary.
dimensionedScalar exp(const dimensionedScalar &ds)
scalar kappa_
Von-Karman constant (default = 0.41)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
#define forAll(list, i)
Loop across all elements in list.
alphatFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Base class for surface film models.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
This boundary condition provides a turbulent thermal diffusivity condition when using wall functions,...
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
scalar Prt_
Turbulent Prandtl number (default = 0.85)
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.
const std::string patch
OpenFOAM patch number as a std::string.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Foam::fvPatchFieldMapper.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
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.