Go to the documentation of this file.
56 fixedValueFvPatchScalarField(
p, iF),
57 filmRegionName_(
"surfaceFilmProperties"),
75 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
76 filmRegionName_(ptf.filmRegionName_),
78 yPlusCrit_(ptf.yPlusCrit_),
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),
113 filmRegionName_(fwfpsf.filmRegionName_),
115 yPlusCrit_(fwfpsf.yPlusCrit_),
117 kappa_(fwfpsf.kappa_),
129 fixedValueFvPatchScalarField(fwfpsf, iF),
130 filmRegionName_(fwfpsf.filmRegionName_),
132 yPlusCrit_(fwfpsf.yPlusCrit_),
134 kappa_(fwfpsf.kappa_),
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();
243 os.writeEntryIfDifferent<
word>
246 "surfaceFilmProperties",
249 os.writeEntry(
"B",
B_);
251 os.writeEntry(
"Cmu",
Cmu_);
253 os.writeEntry(
"Prt",
Prt_);
254 writeEntry(
"value",
os);
263 alphatFilmWallFunctionFvPatchScalarField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
A class for handling words, derived from Foam::string.
constexpr const char *const group
A class for managing temporary objects.
makePatchTypeField(fvPatchScalarField, alphatFilmWallFunctionFvPatchScalarField)
static const word propertiesName
virtual tmp< volScalarField > mu() const =0
dimensionedScalar exp(const dimensionedScalar &ds)
const nearWallDist & y() const
label min(const labelHashSet &set, label minValue=labelMax)
alphatFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Generic templated field type.
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)
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > k() const =0
virtual void write(Ostream &) const
static int & msgType() noexcept
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
dimensionedScalar sqrt(const dimensionedScalar &ds)
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
Foam::fvPatchFieldMapper.
static word groupName(StringType base, const word &group)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
virtual tmp< volScalarField > alpha() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void updateCoeffs()