Index
Clapeyron.AlyLeeIdeal
Clapeyron.BasicIdeal
Clapeyron.JobackIdeal
Clapeyron.LJRefIdeal
Clapeyron.MonomerIdeal
Clapeyron.PPDSIdeal
Clapeyron.ReidIdeal
Clapeyron.WalkerIdeal
Correlation Models
Correlation models are, as their name says, fitted equations that express one property of a compound. They meant to be used in conjunction with other models (like Activity models that require a saturated liquid volume), or via a CompositeModel
. Because they only overload one property, the way to define a correlation is different than normal EoSModel
s.
Saturation Correlations
Saturation Correlations are any EoSModel
that are subtypes of SaturationModel
return psat(T)
and the upper limit (Tc,Pc)
pair. To define saturation correlations, you need to overload:
function crit_pure(model::MySaturationModel <: SaturationModel)
...
return (Tc,Pc,NaN)
end
function Clapeyron.saturation_pressure_impl(model::MySaturationModel <: SaturationModel,T,::SaturationCorrelation)
...
return (psat,NaN,NaN)
end
Saturation Models and Types
Clapeyron.SaturationModel
— TypeSaturationModel <: EoSModel
Abstract type for Saturation correlation models.
Clapeyron.SaturationCorrelation
— TypeSaturationCorrelation <: SaturationMethod
Saturation method used for dispatch on saturation correlations.
Clapeyron.LeeKeslerSat
— TypeLeeKeslerSat <: SaturationModel
LeeKeslerSat(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
acentricfactor
: Single Parameter (Float64
) - acentric factor
Description
Lee-Kesler correlation for saturation pressure:
psat(T) = f₀ + ω•f₁
Tr = T/Tc
f₀ = 5.92714 - 6.09648/Tr - 1.28862•log(Tr) + 0.169347•Tr⁶
f₁ = 15.2518 - 15.6875/Tr - 13.4721•log(Tr) + 0.43577•Tr⁶
Model Construction Examples
# Using the default database
sat = LeeKeslerSat("water") #single input
sat = LeeKeslerSat(["water","ethanol"]) #multiple components
# User-provided parameters, passing files or folders
sat = LeeKeslerSat(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical.csv"])
# User-provided parameters, passing parameters directly
sat = LeeKeslerSat(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
acentricfactor = [-0.03,-0.21])
)
References
- Lee, B. I., & Kesler, M. G. (1975). A generalized thermodynamic correlation based on three-parameter corresponding states. AIChE journal. American Institute of Chemical Engineers, 21(3), 510–527. doi:10.1002/aic.690210313
Clapeyron.DIPPR101Sat
— TypeDIPPR101Sat <: SaturationModel
DIPPR101Sat(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
A
: Single Parameter (Float64
)B
: Single Parameter (Float64
)C
: Single Parameter (Float64
)D
: Single Parameter (Float64
)E
: Single Parameter (Float64
)Tmin
: Single Parameter (Float64
) - mininum Temperature range[K]
Tmax
: Single Parameter (Float64
) - maximum Temperature range[K]
Model Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
A
: Single Parameter (Float64
)B
: Single Parameter (Float64
)C
: Single Parameter (Float64
)D
: Single Parameter (Float64
)E
: Single Parameter (Float64
)Tmin
: Single Parameter (Float64
) - mininum Temperature range[K]
Tmax
: Single Parameter (Float64
) - maximum Temperature range[K]
Description
DIPPR 101 Equation for saturation pressure:
psat(T) = exp(A + B/T + C•log(T) + D•T^E)
References
- Design Institute for Physical Properties, 1996. DIPPR Project 801 DIPPR/AIChE
Clapeyron.AntoineEqSat
— TypeAntoineEqSat <: SaturationModel
AntoineEqSat(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
A
: Single Parameter (Float64
) - First coefficient[dimensionless]
B
: Single Parameter (Float64
) - Second coefficent[°C]
C
: Single Parameter (Float64
) - Third coefficent[°C]
Tmin
: Single Parameter (Float64
) - Mininum Temperature range[K]
Tmax
: Single Parameter (Float64
) - Maximum Temperature range[K]
Model Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
A
: Single Parameter (Float64
) - First coefficient[dimensionless]
B
: Single Parameter (Float64
) - Second coefficent[°C]
C
: Single Parameter (Float64
) - Third coefficent[°C]
Tmin
: Single Parameter (Float64
) - Mininum Temperature range[K]
Tmax
: Single Parameter (Float64
) - Maximum Temperature range[K]
Description
Antoine Equation for saturation pressure:
psat(T) = 10^(A + B/(T + C))
Returns saturation pressure of a pure substance, given temperature T
.
- Temperature
T
[K]
- Saturation Pressure
psat
[Pa]
References
- Antoine, C. (1888), «Tensions des vapeurs; nouvelle relation entre les tensions et les températures», Comptes Rendus des Séances de l'Académie des Sciences 107: 681-684, 778-780, 836-837.
Liquid Volume Correlations
Liquid Volume Correlations are any EoSModel
that are subtypes of LiquidVolumeModel
.
Most Liquid volume correlations are only defined by a volume equation, but some models are defined in terms of the gibbs free energy.
Clapeyron.RackettLiquid
— TypeRackettLiquid(components;
userlocations::Vector{String}=String[],
verbose::Bool=false)
Input parameters
Tc
: Single Parameter (Float64) - Critical Temperature[K]
Pc
: Single Parameter (Float64) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Critical Volume[m³·mol⁻¹]
Model Parameters
Tc
: Single Parameter (Float64) - Critical Temperature[K]
Pc
: Single Parameter (Float64) - Critical Pressure[Pa]
Zc
: Single Parameter (Float64) - Critical Compressibility Factor
Description
Rackett Equation of State for saturated liquids. It is independent of the pressure.
Tr = T/Tc
V = (R̄Tc/Pc)Zc^(1+(1-Tr)^(2/7))
Model Construction Examples
# Using the default database
model = RackettLiquid("water") #single input
model = RackettLiquid(["water","ethanol"]) #multiple components
# User-provided parameters, passing files or folders
model = RackettLiquid(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical.csv"])
# User-provided parameters, passing parameters directly
model = RackettLiquid(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Vc = [4.25e-5, 6.43e-5],
Pc = [2679000, 1296400])
)
References
- Rackett, H. G. (1970). Equation of state for saturated liquids. Journal of Chemical and Engineering Data, 15(4), 514–517. doi:10.1021/je60047a012
Clapeyron.YamadaGunnLiquid
— FunctionYamadaGunnLiquid(components;
userlocations::Vector{String}=String[],
verbose::Bool=false)::RackettLiquid
Input parameters
Tc
: Single Parameter (Float64) - Critical Temperature[K]
Pc
: Single Parameter (Float64) - Critical Pressure[Pa]
acentricfactor
: Single Parameter (Float64
) - Acentric Factor
Model Parameters
Tc
: Single Parameter (Float64) - Critical Temperature[K]
Pc
: Single Parameter (Float64) - Critical Pressure[Pa]
Zc
: Single Parameter (Float64) - Critical Compressibility Factor
Description
The Yamada-Gunn equation of state is a modification of the Rackett equation of state that uses a different approach to calculate the compressibility factor Zc
:
Tr = T/Tc
Zc = 0.29056 - 0.08775ω
V = (R̄Tc/Pc)Zc^(1+(1-Tr)^(2/7))
It can be used as a substitute of RackettLiquid
when Vc
is not known.
Model Construction Examples
# Using the default database
model = YamadaGunnLiquid("water") #single input
model = YamadaGunnLiquid(["water","ethanol"]) #multiple components
# User-provided parameters, passing files or folders
model = YamadaGunnLiquid(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical.csv"])
# User-provided parameters, passing parameters directly
model = YamadaGunnLiquid(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
acentricfactor = [-0.03,-0.21])
)
References
- Rackett, H. G. (1970). Equation of state for saturated liquids. Journal of Chemical and Engineering Data, 15(4), 514–517. doi:10.1021/je60047a012
- Gunn, R. D., & Yamada, T. (1971). A corresponding states correlation of saturated liquid volumes. AIChE Journal. American Institute of Chemical Engineers, 17(6), 1341–1345. doi:10.1002/aic.690170613
Clapeyron.COSTALD
— TypeCOSTALD(components;
userlocations::Vector{String}=String[],
verbose::Bool=false)
Input parameters
Tc
: Single Parameter (Float64) - Critical Temperature[K]
Vc
: Single Parameter (Float64
) - Critical Volume[m³·mol⁻¹]
acentricfactor
: Single Parameter (Float64
) - Acentric Factor
Description
COSTALD Equation of State for saturated liquids. It is independent of the pressure.
Model Construction Examples
# Using the default database
model = COSTALD("water") #single input
model = COSTALD(["water","ethanol"]) #multiple components
# User-provided parameters, passing files or folders
model = COSTALD(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical.csv"])
# User-provided parameters, passing parameters directly
model = COSTALD(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Vc = [4.25e-5, 6.43e-5],
acentricfactor = [-0.03,-0.21])
)
References
Hankinson, R. W., & Thomson, G. H. (1979). A new correlation for saturated densities of liquids and their mixtures. AIChE Journal. American Institute of Chemical Engineers, 25(4), 653–663. doi:10.1002/aic.690250412
Missing docstring for GrenkeElliotWater
. Check Documenter's build log for details.
Clapeyron.HoltenWater
— TypeGrenkeElliottWater <: GibbsBasedModel
GrenkeElliottWater()
Input parameters
None
Description
Holten's model for liquid water at low temperatures (homogeneous ice nucleation temperature up to 300 K) and high pressures (up to 400 MPa, but can be extrapolated up to 1000 MPa).
g(P,T) = R*Tc*ĝ
ĝ = ĝᴬ + T̂*(x*L + x*log(x) + (1 - x)*log(1 - x) + ω*x*(1-x))
x: solution of L + log(x) - log(1 - x) + ω*(1 - 2*x) = 0
L = L₀*K₂*(1 + k₀*k₂ + k₁*(p + k₂*t) - K₁)/(2*k₁*k₂)
K₁ = sqrt((1 + k₀*k₂ + k₁*(p - k₂*t))^2 - 4*k₀*k₁*k₂*(p - k₂*t))
K₂ = sqrt(1 + k₂*k₂)
t = (T - Tc)/Tc
p = P/(R*Tc*ρ0)
ω = 2 + ω₀*p
ĝᴬ = ∑(cᵢ * τ^aᵢ * π^bᵢ * exp(-dᵢ*π))
Where x is the fraction of water molecules in with low-density structure.
References
- Holten, V., Sengers, J. V., & Anisimov, M. A. (2014). Equation of state for supercooled water at pressures up to 400 MPa. Journal of Physical and Chemical Reference Data, 43(4), 043101. doi:10.1063/1.4895593
Virial Models
Virial models are defined in terms of the second virial coefficient, B(T,z)
. The reduced residual Helmholtz energy is defined as:
$\frac{A_\mathrm{res}}{Nk_\mathrm{B}T} = \frac{B}{V}$,
To implement a virial model, it is necessary to overload Clapeyron.second_virial_coefficient_impl(model::<:SecondVirialModel,T,z)
.
Clapeyron.AbbottVirial
— TypeAbbottVirial <: SecondVirialModel
AbbottVirial(components;
idealmodel = BasicIdeal,
userlocations = String[],
ideal_userlocations = String[],
verbose = false)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
acentricfactor
: Single Parameter (Float64
) - Acentric FactorMw
: Single Parameter (Float64
) - Molecular Weight[g·mol⁻¹]
Input models
idealmodel
: Ideal Model
Description
Virial model using Corresponding State Principles:
B = ∑xᵢxⱼBᵢⱼ
Bᵢⱼ = BrᵢⱼRTcᵢⱼ/Pcᵢⱼ
Brᵢⱼ = B₀ + ωᵢⱼB₁
B₀ = 0.083 + 0.422/Trᵢⱼ^1.6
B₁ = 0.139 - 0.172/Trᵢⱼ^4.2
Trᵢⱼ = T/Tcᵢⱼ
Tcᵢⱼ = √TcᵢTcⱼ
Pcᵢⱼ = (Pcᵢ + Pcⱼ)/2
ωᵢⱼ = (ωᵢ + ωⱼ)/2
Model Construction Examples
# Using the default database
model = AbbottVirial("water") #single input
model = AbbottVirial(["water","ethanol"]) #multiple components
model = AbbottVirial(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
# Passing a prebuilt model
my_idealmodel = MonomerIdeal(["neon","hydrogen"];userlocations = (;Mw = [20.17, 2.]))
model = AbbottVirial(["neon","hydrogen"],idealmodel = my_idealmodel)
# User-provided parameters, passing files or folders
model = AbbottVirial(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical.csv"])
# User-provided parameters, passing parameters directly
model = AbbottVirial(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
Mw = [20.17, 2.],
acentricfactor = [-0.03,-0.21])
)
References
- Smith, H. C. Van Ness Joseph M. Introduction to Chemical Engineering Thermodynamics 4E 1987.
Clapeyron.TsonopoulosVirial
— TypeTsonopoulosVirial <: SecondVirialModel
TsonopoulosVirial(components;
idealmodel = BasicIdeal,
userlocations = String[],
ideal_userlocations = String[],
verbose = false)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
acentricfactor
: Single Parameter (Float64
) - Acentric FactorMw
: Single Parameter (Float64
) - Molecular Weight[g·mol⁻¹]
Input models
idealmodel
: Ideal Model
Description
Virial model using Corresponding State Principles:
B = ∑xᵢxⱼBᵢⱼ
Bᵢⱼ = BrᵢⱼRTcᵢⱼ/Pcᵢⱼ
Brᵢⱼ = B₀ + ωᵢⱼB₁
B₀ = 0.1445 - 0.330/Trᵢⱼ - 0.1385/Trᵢⱼ^2 - 0.0121/Trᵢⱼ^3 - 0.000607/Trᵢⱼ^8
B₁ = 0.0637 + 0.331/Trᵢⱼ - 0.423/Trᵢⱼ^2 - 0.423/Trᵢⱼ^3 - 0.008/Trᵢⱼ^8
Trᵢⱼ = T/Tcᵢⱼ
Tcᵢⱼ = √TcᵢTcⱼ
Pcᵢⱼ = (Pcᵢ + Pcⱼ)/2
ωᵢⱼ = (ωᵢ + ωⱼ)/2
Model Construction Examples
# Using the default database
model = TsonopoulosVirial("water") #single input
model = TsonopoulosVirial(["water","ethanol"]) #multiple components
model = TsonopoulosVirial(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
# Passing a prebuilt model
my_idealmodel = MonomerIdeal(["neon","hydrogen"];userlocations = (;Mw = [20.17, 2.]))
model = TsonopoulosVirial(["neon","hydrogen"],idealmodel = my_idealmodel)
# User-provided parameters, passing files or folders
model = TsonopoulosVirial(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical.csv"])
# User-provided parameters, passing parameters directly
model = TsonopoulosVirial(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
Mw = [20.17, 2.],
acentricfactor = [-0.03,-0.21])
)
References
- Tsonopoulos, C. (1974). An empirical correlation of second virial coefficients. AIChE Journal. American Institute of Chemical Engineers, 20(2), 263–272. doi:10.1002/aic.690200209
Clapeyron.EoSVirial2
— TypeEoSVirial2 <: SecondVirialModel
EoSVirial2(model;idealmodel = idealmodel(model))
Input models
model
: Model providing second virial coefficient is obtainedidealmodel
: Ideal Model
Description
Virial model, that just calls the second virial coefficient of the underlying model.
B(T,z) = B(model,T,z)
Solid Models
Solid models provide simple approximations to the excess chemical potential in the solid phase. Intended to be used in conjunction with a liquid model within a CompositeModel.
Clapeyron.SolidHfus
— TypeSolidHfusModel <: EoSModel
SolidHfus(components;
userlocations = String[],
verbose::Bool=false)
Parameters
Hfus
: Single Parameter (Float64
) - Enthalpy of Fusion at 1 bar[J·mol⁻¹]
Tm
: Single Parameter (Float64
) - Melting Temperature[K]
CpSL
: Single Parameter (Float64
) (optional) - Heat Capacity of the Solid-Liquid Phase Transition[J·mol⁻¹·K⁻¹]
Description
Approximation of the excess chemical potential in the solid phase (CpSL
is not necessary by default):
ln(xᵢγᵢ) = Hfusᵢ*T*(1/Tmᵢ-1/T)-CpSLᵢ/R̄*(Tmᵢ/T-1-log(Tmᵢ/T))
Clapeyron.SolidKs
— TypeSolidKsModel <: EoSModel
SolidKs(components;
userlocations = String[],
verbose::Bool=false)
Parameters
Gform
: Single Parameter (Float64
) - Formation Gibbs energy in water at infinite dilution 1 bar and the reference temperature[J·mol⁻¹]
Hform
: Single Parameter (Float64
) -Formation enthalpy in water at infinite dilution 1 bar and the reference temperature[J·mol⁻¹]
Tref
: Single Parameter (Float64
) - Reference temperature[K]
Description
Approximation of the excess chemical potential in the solid phase, using enthalpies and gibbs energies of formation where the liquid reference is at infinite dilution in water:
ln(xᵢγᵢ) = -Gformᵢ*T/Trefᵢ - Hformᵢ*(1 - T/Trefᵢ)
Clapeyron.IAPWS06
— TypeIAPWS06 <: EmpiricHelmholtzModel
IAPWS06()
Input parameters
None
Description
IAPWS (International Association for the Properties of Water and Steam) Pure water Model for ice Ih, 2009 update.
References
- Feistel, R., & Wagner, W. (2006). A new equation of state for H2O ice Ih. Journal of Physical and Chemical Reference Data, 35(2), 1021–1047. doi:10.1063/1.2183324
- IAPWS R10-06 (2009). Revised Release on the Equation of State 2006 for H2O Ice Ih