Go to the documentation of this file.
45 components_ =
dict.toc();
46 properties_.setSize(components_.size());
113 Tpt += X[i]*properties_[i].Tt();
131 if (
p >= pv(
p, Thi, X))
135 else if (
p < pv(
p, Tlo, X))
138 <<
"Pressure below triple point pressure: "
139 <<
"p = " <<
p <<
" < Pt = " << pv(
p, Tlo, X) <<
nl <<
endl;
144 scalar
T = (Thi + Tlo)*0.5;
146 while ((Thi - Tlo) > 1.0e-4)
148 if ((pv(
p,
T, X) -
p) <= 0.0)
170 Tpc += X[i]*properties_[i].Tc();
184 Vc += X[i]*properties_[i].Vc();
185 Zc += X[i]*properties_[i].Zc();
188 return RR*Zc*Tpc(X)/Vc;
198 omega += X[i]*properties_[i].omega();
219 scalar Ti =
min(TrMax*properties_[i].Tc(), Tl);
220 Xs[i] = properties_[i].pv(
p, Ti)*Xl[i]/
p;
233 W += X[i]*properties_[i].W();
247 Y[i] = X[i]*properties_[i].W();
264 X[i] =
Y[i]/properties_[i].W();
288 scalar Ti =
min(TrMax*properties_[i].Tc(),
T);
289 scalar
rho = properties_[i].rho(
p, Ti);
293 scalar Yi = X[i]*properties_[i].W();
318 scalar Yi = X[i]*properties_[i].W();
321 scalar Ti =
min(TrMax*properties_[i].Tc(),
T);
322 pv += Yi*properties_[i].pv(
p, Ti);
344 scalar Yi = X[i]*properties_[i].W();
347 scalar Ti =
min(TrMax*properties_[i].Tc(),
T);
348 hl += Yi*properties_[i].hl(
p, Ti);
370 scalar Yi = X[i]*properties_[i].W();
373 scalar Ti =
min(TrMax*properties_[i].Tc(),
T);
374 Cp += Yi*properties_[i].Cp(
p, Ti);
397 scalar Ti =
min(TrMax*properties_[i].Tc(),
T);
398 scalar Pvs = properties_[i].pv(
p, Ti);
410 scalar Ti =
min(TrMax*properties_[i].Tc(),
T);
411 sigma += Xs[i]*properties_[i].sigma(
p, Ti);
432 scalar Ti =
min(TrMax*properties_[i].Tc(),
T);
433 mu += X[i]*
log(properties_[i].
mu(
p, Ti));
454 scalar Ti =
min(TrMax*properties_[i].Tc(),
T);
456 scalar Vi = properties_[i].W()/properties_[i].rho(
p, Ti);
467 scalar Ti =
min(TrMax*properties_[i].Tc(),
T);
471 scalar Tj =
min(TrMax*properties_[j].Tc(),
T);
476 1.0/properties_[i].K(
p, Ti)
477 + 1.0/properties_[j].K(
p, Tj)
479 K += phii[i]*phii[j]*Kij;
501 scalar Ti =
min(TrMax*properties_[i].Tc(),
T);
502 Dinv += X[i]/properties_[i].D(
p, Ti);
const scalar RR
Universal gas constant (default in [J/(kmol K)])
scalar pvInvert(const scalar p, const scalarField &X) const
Invert the vapour pressure relationship to retrieve the boiling.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
scalarField Xs(const scalar p, const scalar Tg, const scalar Tl, const scalarField &Xg, const scalarField &Xl) const
Return the surface molar fractions.
const dimensionedScalar mu
Atomic mass unit.
#define forAll(list, i)
Loop across all elements in list.
scalar Tc(const scalarField &X) const
Calculate the critical temperature of mixture.
scalar hl(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture latent heat [J/kg].
scalarField Y(const scalarField &X) const
Returns the mass fractions corresponding to the given mole fractions.
scalar mu(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture viscosity [Pa s].
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar exp(const dimensionedScalar &ds)
scalar Tpc(const scalarField &X) const
Return pseudocritical temperature according to Kay's rule.
static autoPtr< liquidMixtureProperties > New(const dictionary &)
Select construct from dictionary.
scalar rho(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture density [kg/m^3].
liquidMixtureProperties(const dictionary &dict)
Construct from dictionary.
CGAL::Exact_predicates_exact_constructions_kernel K
scalarField X(const scalarField &Y) const
Returns the mole fractions corresponding to the given mass fractions.
static const scalar TrMax
Maximum reduced temperature.
scalar sigma(const scalar p, const scalar T, const scalarField &X) const
Estimate mixture surface tension [N/m].
scalar Tpt(const scalarField &X) const
Return pseudo triple point temperature (mole averaged formulation)
Pre-declare SubField and related Field type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
List< word > components_
The names of the liquids.
scalar Ppc(const scalarField &X) const
Return pseudocritical pressure (modified Prausnitz and Gunn)
scalar W(const scalarField &X) const
Calculate the mean molecular weight [kg/kmol].
A list of keyword definitions, which are a keyword followed by any number of values (e....
dimensionedScalar log(const dimensionedScalar &ds)
scalar D(const scalar p, const scalar T, const scalarField &X) const
Vapour diffussivity [m2/s].
scalar Cp(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture heat capacity [J/(kg K)].
scalar pv(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture vapour pressure [Pa].
PtrList< liquidProperties > properties_
The liquid properties.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
scalar omega(const scalarField &X) const
Return mixture accentric factor.
scalar K(const scalar p, const scalar T, const scalarField &X) const
Estimate thermal conductivity [W/(m K)].
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define WarningInFunction
Report a warning using Foam::Warning.
PtrList< volScalarField > & Y