Go to the documentation of this file.
40 fixedValueFvPatchScalarField(
p, iF),
41 filmRegionName_(
"surfaceFilmProperties"),
57 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
60 a_(ptf.
a_().clone().ptr()),
61 omega_(ptf.
omega_().clone().ptr())
73 fixedValueFvPatchScalarField(
p, iF),
76 dict.lookupOrDefault<
word>(
"filmRegion",
"surfaceFilmProperties")
92 fixedValueFvPatchScalarField(wmfrhpsf),
94 GammaMean_(wmfrhpsf.
GammaMean_().clone().ptr()),
95 a_(wmfrhpsf.
a_().clone().ptr()),
96 omega_(wmfrhpsf.
omega_().clone().ptr())
107 fixedValueFvPatchScalarField(wmfrhpsf, iF),
109 GammaMean_(wmfrhpsf.
GammaMean_().clone().ptr()),
110 a_(wmfrhpsf.
a_().clone().ptr()),
111 omega_(wmfrhpsf.
omega_().clone().ptr())
124 const label patchI = patch().index();
144 if (patch().size() && (
max(
mag(gTan)) < SMALL))
147 <<
"is designed to operate on patches inclined with respect to "
157 nTan /=
mag(nTan) + ROOTVSMALL;
166 const scalar t = db().time().timeOutputValue();
169 const scalar a =
a_->value(t);
170 const scalar omega =
omega_->value(t);
175 const scalarField mup(
mu.boundaryField()[patchI].patchInternalField());
178 const scalarField rhop(
rho.boundaryField()[patchI].patchInternalField());
184 pow(3.0*
sqr(mup/rhop)/(gTan + ROOTVSMALL), 0.333)*
pow(
Re, 0.333)
187 fixedValueFvPatchScalarField::updateCoeffs();
197 writeEntryIfDifferent<word>
201 "surfaceFilmProperties",
204 GammaMean_->writeData(os);
206 omega_->writeData(os);
207 writeEntry(
"value", os);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
inclinedFilmNusseltHeightFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for handling words, derived from string.
const dimensionedScalar mu
Atomic mass unit.
const dimensionedScalar G
Newtonian constant of gravitation.
dimensionedScalar sin(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cf)
autoPtr< DataEntry< scalar > > GammaMean_
Mean mass flow rate per unit length [kg/s/m].
const Time & time() const
Return the reference to the time database.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const scalar twoPi(2 *pi)
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
autoPtr< DataEntry< scalar > > a_
Perturbation amplitude [m].
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.
virtual void write(Ostream &) const
Write.
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.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField scalarField(fieldObject, mesh)
autoPtr< DataEntry< scalar > > omega_
Perturbation frequency [rad/s/m].
Foam::fvPatchFieldMapper.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Film height boundary condition for inclined films that imposes a sinusoidal perturbation on top of a ...
Generic GeometricField class.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
#define WarningInFunction
Report a warning using Foam::Warning.
word filmRegionName_
Name of film region.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...