55 void Foam::fv::phaseSurfaceBoiling::readCoeffs(
const dictionary&
dict)
57 saturationModelPtr_.reset
61 "saturationTemperature",
66 partitioningModel_.
reset
70 dict.subDict(
"partitioningModel")
74 nucleationSiteModel_.
reset
78 dict.subDict(
"nucleationSiteModel")
82 departureDiameterModel_.
reset
86 dict.subDict(
"departureDiameterModel")
90 departureFrequencyModel_.
reset
94 dict.subDict(
"departureFrequencyModel")
100 void Foam::fv::phaseSurfaceBoiling::correctMDot()
const
112 const rhoThermo& solidThermo = solid_.thermo();
114 const rhoThermo& liquidThermo = liquid_.thermo();
116 const rhoThermo& vapourThermo = vapour_.thermo();
121 const Pair<tmp<volScalarField>> Hs =
122 solver_.heatTransfer.Hs(liquid_, solid_, scalar(0));
127 (solidH*solidT + liquidH*liquidT + qEvaporative_ + qQuenching_)
128 /
max(solidH + liquidH, rootVSmallH)
133 saturationModelPtr_->Tsat(liquid_.fluidThermo().p()())
138 vapourThermo.ha(liquid_.fluidThermo().p()(), Tsat)
139 - liquidThermo.ha()()
144 partitioningModel_->wetFraction(liquid_()/
max(1 - solid_(), small));
148 departureDiameterModel_->dDeparture
160 departureFrequencyModel_->fDeparture
172 nucleationSiteDensity_ =
173 nucleationSiteModel_->nucleationSiteDensity
185 const tmp<volScalarField> tliquidRho(liquidThermo.rho());
187 const tmp<volScalarField> tvapourRho(vapourThermo.rho());
195 liquidRho*liquidCp*(Tsat -
min(Tsurface, Tsat))/(liquidRho*L)
199 wetFraction_*4.8*
exp(
min(-Ja/80,
log(vGreat)))
203 min(
pi*
sqr(dDeparture_)*nucleationSiteDensity_*Al/4, scalar(1))
207 min(
pi*
sqr(dDeparture_)*nucleationSiteDensity_*Al/4, scalar(5))
219 +
f*(1.0/6.0)*A2E*dDeparture_*vapourRho*fDeparture_*Av;
222 qEvaporative_ = mDot_*L;
233 /
max(fDeparture_, rootVSmallF)
234 /(
pi*(liquidkappa/liquidCp)/liquidRho)
241 +
f*A2*hQuenching*
max(Tsurface - liquidThermo.T(), rooVSmallT)*Av;
250 const word& modelType,
257 solver_(
mesh().lookupObject<solvers::multiphaseEuler>(
solver::typeName)),
258 fluid_(solver_.fluid),
259 liquid_(fluid_.phases()[phaseNames().
first()]),
260 vapour_(fluid_.phases()[phaseNames().
second()]),
261 solid_(fluid_.phases()[
dict.lookup(
"phase")]),
262 saturationModelPtr_(nullptr),
263 partitioningModel_(nullptr),
264 nucleationSiteModel_(nullptr),
265 departureDiameterModel_(nullptr),
266 departureFrequencyModel_(nullptr),
267 pressureEquationIndex_(-1),
272 name +
":wetFraction",
285 name +
":dDeparture",
298 name +
":fDeparture",
307 nucleationSiteDensity_
311 name +
":nucleationSiteDensity",
324 name +
":qQuenching",
337 name +
":qEvaporative",
370 || fieldName == solid_.thermo().he().name();
381 name() +
":Lfraction",
428 && (&
rho == &liquid_.rho() || &
rho == &vapour_.rho())
429 && &eqn.
psi() == &solver_.p_rgh
434 if (pressureEquationIndex_ % 2 == 0) correctMDot();
435 pressureEquationIndex_ ++;
453 if (&
he == &liquid_.thermo().he() || &
he == &solid_.thermo().he())
455 const scalar
sign = &
he == &solid_.thermo().he() ? -1 : +1;
457 eqn +=
sign*(qEvaporative_ + qQuenching_);
461 if (&
he == &liquid_.thermo().he())
471 pressureEquationIndex_ = 0;
482 readCoeffs(coeffs(
dict));
Macros for easy insertion into run-time selection tables.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
DimensionedField< Type, GeoMesh, PrimitiveField > Internal
Type of the internal field from which this GeometricField is derived.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Mesh data needed to do the Finite Volume discretisation.
const fvSolution & solution() const
Return the fvSolution.
Finite volume model abstract base class.
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Mix-in interface for nucleation models. Provides access to properties of the nucleation process,...
Base class for phase change models.
virtual bool read(const dictionary &dict)
Read source dictionary.
void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &heOrYi, fvMatrix< scalar > &eqn) const
Override the energy equation to add the phase change heat, or.
Model for nucleate wall boiling on the surface of a third (solid) phase.
virtual tmp< DimensionedField< scalar, volMesh > > Lfraction() const
Return the fraction of the latent heat that is transferred into.
virtual void correct()
Correct the fvModel.
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the mass transfer rate.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual tmp< DimensionedField< scalar, volMesh > > d() const
Return the diameter of nuclei.
virtual tmp< DimensionedField< scalar, volMesh > > tau() const
Return the nucleation time scale.
void addSup(const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn) const
Override the pressure equation to add the mass transfer rate.
phaseSurfaceBoiling(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
virtual tmp< DimensionedField< scalar, volMesh > > nDot() const
Return the number rate at which nuclei are generated.
static const dimensionSet dimK
Coefficient dimensions.
static autoPtr< saturationTemperatureModel > New(const dictionary &dict)
Select with dictionary.
scalar fieldRelaxationFactor(const word &name) const
Return the relaxation factor for the given field.
Abstract base class for run-time selectable region solvers.
A class for managing temporary objects.
static autoPtr< departureDiameterModel > New(const dictionary &dict)
Select null constructed.
static autoPtr< departureFrequencyModel > New(const dictionary &dict)
Select null constructed.
static autoPtr< nucleationSiteModel > New(const dictionary &dict)
Select null constructed.
static autoPtr< partitioningModel > New(const dictionary &dict)
Select null constructed.
A class for handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
const dimensionSet dimEnergy
const dimensionSet dimless
const dimensionSet dimLength
labelList second(const UList< labelPair > &p)
const dimensionSet dimTemperature
dimensionedScalar log(const dimensionedScalar &ds)
labelList first(const UList< labelPair > &p)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
const dimensionSet dimDensity
const dimensionSet dimVolume
VolField< scalar > volScalarField
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimArea
void inv(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.