Contents
Index
Clapeyron.T_scale
Clapeyron.activity
Clapeyron.activity_coefficient
Clapeyron.adiabatic_index
Clapeyron.antoine_coef
Clapeyron.aqueous_activity
Clapeyron.chemical_potential
Clapeyron.chemical_potential_res
Clapeyron.compressibility_factor
Clapeyron.enthalpy
Clapeyron.entropy
Clapeyron.entropy_res
Clapeyron.fugacity_coefficient
Clapeyron.gibbs_free_energy
Clapeyron.gibbs_free_energy_res
Clapeyron.helmholtz_free_energy
Clapeyron.helmholtz_free_energy_res
Clapeyron.internal_energy
Clapeyron.internal_energy_res
Clapeyron.isentropic_compressibility
Clapeyron.isobaric_expansivity
Clapeyron.isobaric_heat_capacity
Clapeyron.isochoric_heat_capacity
Clapeyron.isothermal_compressibility
Clapeyron.joule_thomson_coefficient
Clapeyron.lb_volume
Clapeyron.mass_density
Clapeyron.mixing
Clapeyron.molar_density
Clapeyron.p_scale
Clapeyron.pip
Clapeyron.pressure
Clapeyron.reference_chemical_potential
Clapeyron.reference_chemical_potential_type
Clapeyron.second_virial_coefficient
Clapeyron.speed_of_sound
Clapeyron.volume
Clapeyron.volume_virial
Clapeyron.x0_crit_pure
Clapeyron.x0_psat
Clapeyron.x0_sat_pure
Clapeyron.x0_saturation_temperature
Clapeyron.x0_volume
Clapeyron.x0_volume_gas
Clapeyron.x0_volume_liquid
Clapeyron.x0_volume_solid
Volume–Temperature Based Properties
Clapeyron.pressure
— Functionpressure(model::EoSModel, V, T, z=SA[1.])
default units: [Pa]
Returns the pressure of the model at a given volume, temperature and composition, defined as:
p = -∂A/∂V
Clapeyron.second_virial_coefficient
— Functionsecond_virial_coefficient(model::EoSModel, T, z=SA[1.])
Default units: [m^3]
Calculates the second virial coefficient B
, defined as:
B = lim(ρ->0)[∂Aᵣ/∂ρ]
where Aᵣ
is the residual helmholtz energy.
Clapeyron.pip
— Functionpip(model::EoSModel,V,T,z=[1.0])
Phase identification parameter Π
. as described in 1. If Π > 1
, then the phase is clasified as a liquid or a liquid-like vapor, being a vapor or vapor-like liquid otherwise.
This identification parameter fails at temperatures and pressures well aboVe the critical point.
Calculated as:
Π = V*((∂²p/∂V∂T)/(∂p/∂T) - (∂²p/∂V²)/(∂p/∂V))
- G. Venkatarathnama, L.R. Oellrich, Identification of the phase of a fluid using partial derivatives of pressure, volume,and temperature without reference to saturation properties: Applications in phase equilibria calculations, Fluid Phase Equilibria 301 (2011) 225–233
Pressure–Temperature Based Bulk Properties
In general almost all bulk properties follow the pattern:
function property(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
V = volume(model, p, T, z; phase, threaded, vol0)
return VT_property(model,V,T,z)
end
So, you can calculate the property with Volume–Temperature variables by calling VT_property(model,V,T,z)
. Another way to do this is by using units, provided by Unitful.jl
:
using Unitful
r = 18u"kg/m^3"
T = 373.15"K"
prop = helmholtz_free_energy(model,r,T,z,output = u"kJ")
Where r
could be any molar or mass density, molar or mass volume, total volume or pressure. It also supports mass and mol amounts defined as units for the composition (z
). If no units are provided for the composition, they will be considered moles.
Methods that require first order VT derivatives
Clapeyron.volume
— Functionvolume(model::EoSModel, p, T, z=SA[1.0]; phase=:unknown, threaded=true, vol0=nothing)
Calculates the volume (m³) of the compound modelled by model
at a certain pressure, temperature and moles. phase
is a Symbol that determines the initial volume root to look for:
- If
phase =:unknown
(Default), it will return the physically correct volume root with the least gibbs energy. - If
phase =:liquid
, it will return the volume of the phase using a liquid initial point. - If
phase =:vapor
, it will return the volume of the phase using a gas initial point. - If
phase =:solid
, it will return the volume of the phase using a solid initial point (only supported for EoS that support a solid phase) - If
phase =:stable
, it will return the physically correct volume root with the least gibbs energy, and perform a stability test on the result.
All volume calculations are checked for mechanical stability, that is: dP/dV <= 0
.
The calculation of both volume roots can be calculated in serial (threaded=false
) or in parallel (threaded=true
).
An initial estimate of the volume vol0
can be optionally be provided.
The volume computation may fail and return NaN
because the default initial point is too far from the actual volume. Providing a value for vol0
may help in these situations. Such a starting point can be found from physical knowledge, or by computing the volume using a different model for example.
The stability check is disabled by default. that means that the volume obtained just follows the the relation P = pressure(model,V,T,z)
. For single component models, this is alright, but phase splits (with different compositions that the input) can and will occur, meaning that the volume solution does not correspond to an existing phase. For unknown multicomponent mixtures, it is recommended to use a phase equilibrium procedure (like tp_flash
) to obtain a list of valid compositions, and then perform a volume calculation over those compositions. You can also pass phase=:stable
to perform the stability test inside the volume solver. Finally, you can perform the stability test after the volume solver:
v = volume(model,p,T,z)
isstable(model,v,T,z)
Clapeyron.helmholtz_free_energy
— Functionhelmholtz_free_energy(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
helmholtz_energy(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J]
Calculates the helmholtz free energy, defined as:
A = eos(model,V(p),T,z)
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_helmholtz_free_energy(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.helmholtz_free_energy_res
— Functionhelmholtz_free_energy_res(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
helmholtz_energy_res(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J]
Calculates the residual helmholtz free energy, defined as:
A = eos_res(model,V(p),T,z)
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_helmholtz_free_energy_res(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.molar_density
— Functionmolar_density(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
default units: [mol/m^3]
Calculates the molar density, defined as:
ρₙ = ∑nᵢ/V
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_molar_density(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.mass_density
— Functionmass_density(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true)
default units: [kg/m^3]
Calculates the mass density, defined as:
ρₙ = Mr/V
Where Mr
is the molecular weight of the model at the input composition.
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_mass_density(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.compressibility_factor
— Functioncompressibility_factor(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Calculates the compressibility factor Z
, defined as:
Z = p*V(p)/R*T
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.gibbs_free_energy
— Functiongibbs_free_energy(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
gibbs_energy(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J]
Calculates the gibbs free energy, defined as:
G = A + p*V
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_gibbs_free_energy(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.gibbs_free_energy_res
— Functiongibbs_free_energy_res(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
gibbs_energy_res(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J]
Calculates the residual gibbs free energy, defined as:
G = Ar - V*∂Ar/∂V
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_gibbs_free_energy_res(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.entropy
— Functionentropy(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J/K]
Calculates entropy, defined as:
S = -∂A/∂T
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_entropy(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.entropy_res
— Functionentropy_res(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J/K]
Calculates residual entropy, defined as:
S = -∂Ares/∂T
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_entropy_res(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.enthalpy
— Functionenthalpy(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J]
Calculates the enthalpy, defined as:
H = A - T * ∂A/∂T - V * ∂A/∂V
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_enthalpy(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.internal_energy
— Functioninternal_energy(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J]
Calculates the internal energy, defined as:
U = A - T * ∂A/∂T
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_internal_energy(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.internal_energy_res
— Functioninternal_energy_res(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J]
Calculates the residual internal energy, defined as:
U = Ar - T * ∂Ar/∂T
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_internal_energy_res(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Methods that require second order VT derivatives
Clapeyron.isochoric_heat_capacity
— Functionisochoric_heat_capacity(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J/K]
Calculates the isochoric heat capacity, defined as:
Cv = -T * ∂²A/∂T²
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_isochoric_heat_capacity(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
This property requires at least second order ideal model temperature derivatives. If you are computing these properties, consider using a different ideal model than the BasicIdeal
default (e.g. EoS(["species"];idealmodel = ReidIdeal)
).
Clapeyron.isobaric_heat_capacity
— Functionisobaric_heat_capacity(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J/K]
Calculates the isobaric heat capacity, defined as:
Cp = -T*(∂²A/∂T² - (∂²A/∂V∂T)^2 / ∂²A/∂V²)
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_isobaric_heat_capacity(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
This property requires at least second order ideal model temperature derivatives. If you are computing these properties, consider using a different ideal model than the BasicIdeal
default (e.g. EoS(["species"];idealmodel = ReidIdeal)
).
Clapeyron.adiabatic_index
— Functionadiabatic_index(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J/K]
Calculates the isobaric heat capacity, defined as:
γ = Cp/Cv
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_adiabatic_index(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
This property requires at least second order ideal model temperature derivatives. If you are computing these properties, consider using a different ideal model than the BasicIdeal
default (e.g. EoS(["species"];idealmodel = ReidIdeal)
).
Clapeyron.isothermal_compressibility
— Functionisothermal_compressibility(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
default units: [Pa^-1]
Calculates the isothermal compressibility, defined as:
κT = (V*∂p/∂V)^-1
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_isothermal_compressibility(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.isentropic_compressibility
— Functionisentropic_compressibility(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
default units: [Pa^-1]
Calculates the isentropic compressibility, defined as:
κS = (V*( ∂²A/∂V² - ∂²A/∂V∂T^2 / ∂²A/∂T² ))^-1
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_isentropic_compressibility(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
This property requires at least second order ideal model temperature derivatives. If you are computing these properties, consider using a different ideal model than the BasicIdeal
default (e.g. EoS(["species"];idealmodel = ReidIdeal)
).
Clapeyron.speed_of_sound
— Functionspeed_of_sound(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
default units: [m/s]
Calculates the speed of sound, defined as:
c = V * √(∂²A/∂V² - ∂²A/∂V∂T^2 / ∂²A/∂T²)/Mr)
Where Mr
is the molecular weight of the model at the input composition.
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_speed_of_sound(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
This property requires at least second order ideal model temperature derivatives. If you are computing these properties, consider using a different ideal model than the BasicIdeal
default (e.g. EoS(["species"];idealmodel = ReidIdeal)
).
Clapeyron.isobaric_expansivity
— Functionisobaric_expansivity(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
default units: [K^-1]
Calculates the isobaric expansivity, defined as:
α = -∂²A/∂V∂T / (V*∂²A/∂V²)
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_isobaric_expansivity(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.joule_thomson_coefficient
— Functionjoule_thomson_coefficient(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
default units: [K/Pa]
Calculates the joule thomson coefficient, defined as:
μⱼₜ = -(∂²A/∂V∂T - ∂²A/∂V² * ((T*∂²A/∂T² + V*∂²A/∂V∂T) / (T*∂²A/∂V∂T + V*∂²A/∂V²)))^-1
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_joule_thomson_coefficient(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
This property requires at least second order ideal model temperature derivatives. If you are computing these properties, consider using a different ideal model than the BasicIdeal
default (e.g. EoS(["species"];idealmodel = ReidIdeal)
).
Methods that require first order composition derivatives
Clapeyron.chemical_potential
— Functionchemical_potential(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J/mol]
Calculates the chemical potential, defined as:
μᵢ = ∂A/∂nᵢ
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_chemical_potential(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.chemical_potential_res
— Functionchemical_potential_res(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Default units: [J/mol]
Calculates the residual chemical potential, defined as:
μresᵢ = ∂Ares/∂nᵢ
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_chemical_potential_res(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.fugacity_coefficient
— Functionfugacity_coefficient(model::EoSModel, p, T, z=SA[1.]; phase=:unknown, threaded=true, vol0=nothing)
Calculates the fugacity coefficient φᵢ, defined as:
log(φᵢ) = μresᵢ/RT - log(Z)
Where μresᵢ
is the vector of residual chemical potentials and Z
is the compressibility factor.
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_fugacity_coefficient(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Activity Coefficient
Clapeyron.reference_chemical_potential
— Functionreference_chemical_potential(model::EoSModel,p,T,reference; phase=:unknown, threaded=true, vol0=nothing)
Returns a reference chemical potential. used in calculation of activity
and actitivy_coefficient. there are two available references:
:pure
: the reference potential is a pure component at specifiedT
,p
andphase
:aqueous
: the chemical potential of the pure components at specifiedT
,p
andphase
:sat_pure_T
: the reference potential is the pure saturated liquid phase at specifiedT
.:zero
: the reference potential is equal to zero for all components (used forActivityModel
)
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Clapeyron.reference_chemical_potential_type
— Functionreference_chemical_potential_type(model)::Symbol
Returns a symbol with the type of reference chemical potential used by the input model
.
Clapeyron.activity_coefficient
— Functionactivity_coefficient(model::EoSModel,p,T,z=SA[1.0];reference = :pure, phase=:unknown, threaded=true, vol0=nothing)
Calculates the activity, defined as:
log(γ*z) = (μ_mixt - μ_ref) / R̄ / T
where μ_mixt
is the chemical potential of the mixture and μ_ref
is the reference chemical potential for the model at p
,T
conditions, calculated via Clapeyron.reference_chemical_potential
. Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_fugacity_coefficient(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
If the μ_ref
keyword argument is not provided, the reference
keyword is used to specify the reference chemical potential..
Clapeyron.activity
— Functionactivity(model::EoSModel,p,T,z=SA[1.0];reference = :pure, phase=:unknown, threaded=true, vol0=nothing)
Calculates the activity, defined as:
log(a) = (μ_mixt - μ_ref) / R̄ / T
where μ_mixt
is the chemical potential of the mixture and μ_ref
is the reference chemical potential for the model at p
,T
conditions, calculated via Clapeyron.reference_chemical_potential
. Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_fugacity_coefficient(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
If the μ_ref
keyword argument is not provided, the reference
keyword is used to specify the reference chemical potential..
Clapeyron.aqueous_activity
— Functionaqueous_activity(model::EoSModel,p,T,z=SA[1.0]; phase=:unknown, threaded=true, vol0=nothing)
Calculates the activity with the reference being infinite dilution in water, defined as:
log(a) = (μ_mixt - μ_inf) / R̄ / T
where μ_mixt
is the chemical potential of the mixture and μ_inf
is the chemical potential of the components at infinite dilution in water.
Internally, it calls Clapeyron.volume
to obtain V
and calculates the property via VT_fugacity_coefficient(model,V,T,z)
.
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Mixing
Clapeyron.mixing
— Functionmixing(model::EoSModel, p, T, z=SA[1.], property; phase=:unknown, threaded=true, vol0=nothing)
Calculates the mixing function for a specified property as:
f_mix = f(p,T,z) - ∑zᵢ*f_pureᵢ(p,T)
The keywords phase
, threaded
and vol0
are passed to the Clapeyron.volume
solver.
Initial guess functions
These methods are considered internal, they don't support Symbolics.jl
or Unitful.jl
overloads.
Clapeyron.lb_volume
— Functionlb_volume(model::EoSModel)
lb_volume(model::EoSModel,z)
lb_volume(model::EoSModel,T,z)
Returns the lower bound volume. It has different meanings depending on the Equation of State, but symbolizes the minimum allowable volume at a certain composition:
- SAFT EoS: the packing volume
- Cubic EoS, covolume (b) parameter
On empiric equations of state, the value is chosen to match the volume of the conditions at maximum pressure and minimum temperature , but the equation itself normally can be evaluated at lower volumes. On SAFT and Cubic EoS, volumes lower than lb_volume
will likely error. The lower bound volume is used for guesses of liquid volumes at a certain pressure, saturated liquid volumes and critical volumes.
In most cases, the lower bound volume is independent of temperature. Some notable exceptions are the Quantum-Corrected Peng-Robinson cubic (QCPR
) and Cubic-plus-Chain (CPC) models. For those, it is better to define the three-argument variant lb_volume(model,T,z)
Clapeyron.T_scale
— FunctionT_scale(model::EoSModel,z)
Represents a temperature scaling factor.
On any EoS based on Critical parameters (Cubic or Empiric EoS), the temperature scaling factor is chosen to be the critical temperature. On SAFT or other molecular EoS, the temperature scaling factor is chosen to be a function of the potential depth ϵ. Used as scaling factors in saturation_pressure
and as input for solving crit_pure
Clapeyron.p_scale
— Functionp_scale(model::EoSModel,z)
Represents a pressure scaling factor.
On any EoS based on Critical parameters (Cubic or Empiric EoS), the pressure scaling factor is chosen to be a function of the critical pressure. On SAFT or other molecular EoS, the pressure scaling factor is chosen to a function of ∑(zᵢϵᵢ(σᵢᵢ)³) Used as scaling factors in saturation_pressure
and as input for solving crit_pure
.
By default, it can be defined as a function of Clapeyron.lb_volume
and Clapeyron.T_scale
Clapeyron.x0_volume
— Functionx0_volume(model,p,T,z; phase = :unknown)
Returns an initial guess of the volume at a pressure, temperature, composition and suggested phase. If the suggested phase is :unknown
or :liquid
, calls x0_volume_liquid
. If the suggested phase is :gas
, calls x0_volume_gas
. If the suggested phase is solid
, calls x0_volume_solid
. Returns NaN
otherwise
Clapeyron.x0_volume_solid
— Functionx0_volume_solid(model,T,z)
x0_volume_solid(model,p,T,z)
Returns an initial guess to the solid volume, dependent on temperature and composition. needs to be defined for EoS that support solid phase. by default returns NaN. can be overrided if the EoS defines is_solid(::EoSModel) = true
Clapeyron.x0_volume_liquid
— Functionx0_volume_liquid(model,T,z)
x0_volume_liquid(model,p,T,z)
Returns an initial guess to the liquid volume, dependent on temperature and composition. by default is 1.25 times lb_volume
.
Clapeyron.x0_volume_gas
— Functionx0_volume_gas(model,p,T,z)
Returns an initial guess to the gas volume, depending of pressure, temperature and composition. by default uses volume_virial
Clapeyron.volume_virial
— Functionvolume_virial(model::EoSModel,p,T,z=SA[1.0])
volume_virial(B::Real,p,T,z=SA[1.0])
Calculates an aproximation to the gas volume at specified pressure, volume and composition, by aproximating:
Z(v) ≈ 1 + B(T)/v
where Z
is the compressibility factor and B
is the second virial coefficient. If B>0
, (over the inversion temperature) returns NaN
. If the solution to the problem is complex (Z = 1 + B/v
implies solving a quadratic polynomial), returns -2*B
. If you pass an EoSModel
as the first argument, B
will be calculated from the EoS at the input T
. You can provide your own second virial coefficient instead of a model.
Clapeyron.x0_sat_pure
— Functionx0_sat_pure(model::EoSModel,T)
x0_sat_pure(model,T,crit)
Returns a 2-tuple corresponding to (Vₗ,Vᵥ)
, where Vₗ
and Vᵥ
are the liquid and vapor initial guesses. Used in saturation_pressure
methods that require initial volume guesses. It can be overloaded to provide more accurate estimates if necessary. If an EoS model provides a fast method for crit_pure
, overloading has_fast_crit_pure
will provide x0_sat_pure
with additional information to improve its accuracy.
Clapeyron.x0_psat
— Functionx0_psat(model::EoSModel, T,crit = nothing)
Initial point for saturation pressure, given the temperature and V,T critical coordinates. On moderate pressures it will use a Zero Pressure initialization. On pressures near the critical point it will switch to spinodal finding. Used in saturation_pressure
methods that require initial pressure guesses. if the initial temperature is over the critical point, it returns NaN
. It can be overloaded to provide more accurate estimates if necessary.
Clapeyron.x0_saturation_temperature
— Functionx0_saturation_temperature(model::EoSModel,p)
Returns a 3-tuple corresponding to (T,Vₗ,Vᵥ)
, T
is the initial guess for temperature and Vₗ
and Vᵥ
are the liquid and vapor initial guesses. Used in saturation_temperature
with AntoineSaturation
.
Clapeyron.antoine_coef
— Functionantoine_coef(model)
should return a 3-Tuple containing reduced Antoine Coefficients. The Coefficients follow the correlation:
lnp̄ = log(p / p_scale(model))
T̃ = T/T_scale(model)
lnp̄ = A - B/(T̄ + C))
By default returns nothing
. This is to use alternative methods in case Antoine coefficients aren't available. Used mainly in single and multicomponent temperature calculations.
Clapeyron.x0_crit_pure
— Functionx0_crit_pure(model::EoSModel)
Returns a 2-tuple corresponding to (k,log10(Vc0))
, where k
is Tc0/T_scale(model,z)