Go to the documentation of this file.
36 namespace regionModels
38 namespace surfaceFilmModels
57 if (!isA<thermoSingleLayer>(
owner_))
60 <<
"Thermo model requires a " << thermoSingleLayer::typeName
61 <<
" film to supply the pressure and temperature, but "
62 <<
owner_.type() <<
" film model selected. "
63 <<
"Use the 'useReferenceValues' flag to employ reference "
67 return refCast<const thermoSingleLayer>(
owner_);
109 name_(
"unknown_liquid"),
112 useReferenceValues_(
readBool(coeffDict_.lookup(
"useReferenceValues"))),
116 initLiquid(coeffDict_);
118 if (useReferenceValues_)
120 coeffDict_.lookup(
"pRef") >> pRef_;
121 coeffDict_.lookup(
"TRef") >> TRef_;
151 return liquidPtr_->rho(
p,
T);
161 return liquidPtr_->mu(
p,
T);
171 return liquidPtr_->sigma(
p,
T);
181 return liquidPtr_->Cp(
p,
T);
191 return liquidPtr_->K(
p,
T);
201 return liquidPtr_->D(
p,
T);
211 return liquidPtr_->hl(
p,
T);
221 return liquidPtr_->pv(
p,
T);
247 owner().regionMesh(),
251 owner().regionMesh(),
253 zeroGradientFvPatchScalarField::typeName
272 rho[cellI] = this->
rho(p[cellI],
T[cellI]);
276 trho().correctBoundaryConditions();
292 owner().regionMesh(),
296 owner().regionMesh(),
298 zeroGradientFvPatchScalarField::typeName
317 mu[cellI] = this->
mu(p[cellI],
T[cellI]);
321 tmu().correctBoundaryConditions();
337 owner().regionMesh(),
341 owner().regionMesh(),
343 zeroGradientFvPatchScalarField::typeName
366 tsigma().correctBoundaryConditions();
382 owner().regionMesh(),
386 owner().regionMesh(),
388 zeroGradientFvPatchScalarField::typeName
407 Cp[cellI] = this->
Cp(p[cellI],
T[cellI]);
411 tCp().correctBoundaryConditions();
427 owner().regionMesh(),
431 owner().regionMesh(),
433 zeroGradientFvPatchScalarField::typeName
456 tkappa().correctBoundaryConditions();
const dimensionSet dimPressure
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const thermoSingleLayer & thermoFilm() const
Return a reference to a thermo film.
A class for handling words, derived from string.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.
scalar pRef_
Reference pressure [pa].
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
A class for managing temporary objects.
const dimensionSet dimEnergy
const dimensionSet dimDensity
virtual const word & name() const
Return the specie name.
virtual tmp< volScalarField > kappa() const
Return thermal conductivity [W/m/K].
Template functions to aid in the implementation of demand driven data.
virtual scalar D(const scalar p, const scalar T) const
Return diffusivity [m2/s].
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
virtual tmp< volScalarField > sigma() const
Return surface tension [kg/s2].
const surfaceFilmModel & owner() const
Return const access to the owner surface film model.
virtual tmp< volScalarField > rho() const
Return density [kg/m3].
scalar W() const
Molecular weight [kg/kmol].
liquidFilmThermo(const liquidFilmThermo &)
Disallow default bitwise copy construct.
bool ownLiquid_
Flag to indicate that model owns the liquid object.
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
void deleteDemandDrivenData(DataPtr &dataPtr)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
The thermophysical properties of a liquidProperties.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
void initLiquid(const dictionary &dict)
Initialise the liquid pointer.
virtual tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
Pre-declare SubField and related Field type.
const dictionary & dict() const
Return const access to the cloud dictionary.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual tmp< volScalarField > mu() const
Return dynamic viscosity [Pa.s].
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const liquidProperties * liquidPtr_
Pointer to the liquid properties.
const volScalarField & pPrimary() const
Pressure / [Pa].
virtual scalar pvInvert(scalar p) const
Invert the vapour pressure relationship to retrieve the.
virtual const volScalarField & T() const
Return the film mean temperature [K].
const dimensionSet dimPower
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.
virtual scalar hl(const scalar p, const scalar T) const
Return latent heat [J/kg].
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
virtual scalar Tb(const scalar p) const
Return boiling temperature [K].
defineTypeNameAndDebug(kinematicSingleLayer, 0)
bool foundObject(const word &name) const
Is the named Type found?
bool useReferenceValues_
Flag to indicate that reference values of p and T should be used.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
tmp< volScalarField > trho
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Base class for film thermo models.
scalar TRef_
Reference temperature [K].
surfaceFilmModel & owner_
Reference to the owner surface film model.
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Base class for surface film models.
virtual scalar W() const
Return molecular weight [kg/kmol].
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Generic GeometricField class.
virtual ~liquidFilmThermo()
Destructor.
virtual scalar pv(const scalar p, const scalar T) const
Return vapour pressure [Pa].