Index
Clapeyron.BMAlpha
Clapeyron.Berthelot
Clapeyron.CPAAlpha
Clapeyron.Clausius
Clapeyron.ClausiusAlpha
Clapeyron.ConstantTranslation
Clapeyron.HVRule
Clapeyron.KU
Clapeyron.KUAlpha
Clapeyron.KayRule
Clapeyron.LCVMRule
Clapeyron.MHV1Rule
Clapeyron.MHV2Rule
Clapeyron.MTAlpha
Clapeyron.MTTranslation
Clapeyron.NoAlpha
Clapeyron.NoTranslation
Clapeyron.PPR78Rule
Clapeyron.PR
Clapeyron.PR78Alpha
Clapeyron.PRAlpha
Clapeyron.PSRKRule
Clapeyron.PTV
Clapeyron.PatelTeja
Clapeyron.PatelTejaAlpha
Clapeyron.PenelouxTranslation
Clapeyron.QCPRRule
Clapeyron.RK
Clapeyron.RKAlpha
Clapeyron.RKPR
Clapeyron.RKPRAlpha
Clapeyron.RackettTranslation
Clapeyron.Soave2019Alpha
Clapeyron.SoaveAlpha
Clapeyron.TwuAlpha
Clapeyron.UMRRule
Clapeyron.VTPRRule
Clapeyron.WSRule
Clapeyron.gErRule
Clapeyron.modWSRule
Clapeyron.sCPAAlpha
Clapeyron.vdW
Clapeyron.vdW1fRule
Clapeyron.EPPR78
Clapeyron.PR78
Clapeyron.PSRK
Clapeyron.QCPR
Clapeyron.SRK
Clapeyron.Twu88Alpha
Clapeyron.UMRPR
Clapeyron.VTPR
Clapeyron.ab_premixing
Clapeyron.cPR
Clapeyron.mixing_rule
Clapeyron.tcPR
Clapeyron.tcPRW
Clapeyron.tcRK
Clapeyron.translation
Clapeyron.α_function
Cubic Models
All cubic models in Clapeyron.jl
follow a common evaluation order:
function CubicModel(args...)
# get params for database, initialize other models, etc.
recombine!(model) # we calculate the mixing rules, caches for the translation models if necessary, etc.
end
function cubic_ab(model::CubicModel,V,T,z=SA[1.0])
invn2 = (one(n)/n)^2
a = model.params.a.values
b = model.params.b.values
α = α_function(model,V,T,z,model.alpha)
c = translation(model,V,T,z,model.translation)
ā,b̄,c̄ = mixing_rule(model,V,T,z,model.mixing,α,a,b,c)
return ā, b̄, c̄
end
function a_res(model::CubicModel,V,T,z,data = (sum(z),cubic_ab(model,V,T,z)))
n, ā, b̄, c̄ = data
# depends on the specific EoS
return result
end
- A Mixing Rule Model creates
aᵢⱼ
andbᵢⱼ
from the critical temperature, critical pressure and a matrix of pair coefficients. - An Alpha Model creates a vector of
αᵢ(T)
values. - A Translation Model creates a vector of
cᵢ
values. - The same Mixing rule, given
aᵢⱼ
,bᵢⱼ
,αᵢ(T)
andcᵢ
returns the the mixture values ofā
,b̄
andc̄
that are then used by the corresponding cubic model. A Mixing Rule can contain activity models to participate in the mixing (for example, Huron–Vidal rules).
Common Definitions
Clapeyron.ab_premixing
— Functionab_premixing(model,mixing,kij = nothing,lij = nothing)
given a model::CubicModel, that has a::PairParam
, b::PairParam
, a mixing::MixingRule and kij
,lij
matrices, ab_premixing
will perform an implace calculation to obtain the values of a
and b
, containing values aᵢⱼ and bᵢⱼ. by default, it performs the van der Wals One-Fluid mixing rule. that is:
aᵢⱼ = sqrt(aᵢ*aⱼ)*(1-kᵢⱼ)
bᵢⱼ = (bᵢ + bⱼ)/2
Clapeyron.mixing_rule
— Functionmixing_rule(model::CubicModel,V,T,z,mixing_model::MixingRule,α,a,b,c)
Interface function used by cubic models. with matrices a
and b
, vectors α
and c
, a model::CubicModel
and mixing_model::MixingRule
, returns the scalars ā
,b̄
and c̄
, corresponding to the values mixed by the amount of components and the specifics of the mixing rule.
Example
function mixing_rule(model::CubicModel,V,T,z,mixing_model::vdW1fRule,α,a,b,c)
∑z = sum(z)
ā = dot(z .* sqrt(α),a,z .* sqrt(α))/(∑z*∑z) #∑∑aᵢⱼxᵢxⱼ√(αᵢαⱼ)
b̄ = dot(z,b,z)/(∑z*∑z) #∑∑bᵢⱼxᵢxⱼ
c̄ = dot(z,c)/∑z ∑cᵢxᵢ
return ā,b̄,c̄
end
Main Models
Clapeyron.vdW
— TypevdW(components;
idealmodel = BasicIdeal,
alpha = NoAlpha,
mixing = vdW1fRule,
activity = nothing,
translation = NoTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
k
: Pair Parameter (Float64
) (optional)l
: Pair Parameter (Float64
) (optional)
Model Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
a
: Pair Parameter (Float64
)b
: Pair Parameter (Float64
)
Input models
idealmodel
: Ideal Modelalpha
: Alpha modelmixing
: Mixing modelactivity
: Activity Model, used in the creation of the mixing model.translation
: Translation Model
Description
van der Wals Equation of state.
P = RT/(V-Nb) + a•α(T)/V²
Model Construction Examples
# Using the default database
model = vdW("water") #single input
model = vdW(["water","ethanol"]) #multiple components
model = vdW(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
model = vdW(["water","ethanol"],alpha = SoaveAlpha) #modifying alpha function
model = vdW(["water","ethanol"],translation = RackettTranslation) #modifying translation
model = vdW(["water","ethanol"],mixing = KayRule) #using another mixing rule
model = vdW(["water","ethanol"],mixing = WSRule, activity = NRTL) #using advanced EoS+gᴱ mixing rule
# Passing a prebuilt model
my_alpha = SoaveAlpha(["ethane","butane"],userlocations = Dict(:acentricfactor => [0.1,0.2]))
model = vdW(["ethane","butane"],alpha = my_alpha)
# User-provided parameters, passing files or folders
model = vdW(["neon","hydrogen"]; userlocations = ["path/to/my/db","cubic/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = vdW(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
Mw = [20.17, 2.],
acentricfactor = [-0.03,-0.21]
k = [0. 0.18; 0.18 0.], #k,l can be ommited in single-component models.
l = [0. 0.01; 0.01 0.])
)
References
- van der Waals JD. Over de Continuiteit van den Gasen Vloeistoftoestand. PhD thesis, University of Leiden; 1873
Clapeyron.Clausius
— TypeClausius(components;
idealmodel = BasicIdeal,
alpha = NoAlpha,
mixing = vdW1fRule,
activity = nothing,
translation = ClausiusTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Molar Volume[m^3/mol]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
k
: Pair Parameter (Float64
) (optional)l
: Pair Parameter (Float64
) (optional)
Model Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Molar Volume[m^3/mol]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
a
: Pair Parameter (Float64
)b
: Pair Parameter (Float64
)
Description
Clausius Equation of state.
P = RT/(v-b) + a•α(T)/((v - Δ₀b)^2)
aᵢᵢ =27/64 * (RTcᵢ)²/Pcᵢ
bᵢᵢ = Vcᵢ - 1/4 * RTcᵢ/Pcᵢ
cᵢ = 3/8 * RTcᵢ/Pcᵢ - Vcᵢ
Δ₀ = ∑cᵢxᵢ/∑bᵢxᵢ
References
- Clausius, R. (1880). Ueber das Verhalten der Kohlensäure in Bezug auf Druck, Volumen und Temperatur. Annalen der Physik, 245(3), 337–357. doi:10.1002/andp.18802450302
Clapeyron.RK
— TypeRK(components;
idealmodel = BasicIdeal,
alpha = PRAlpha,
mixing = vdW1fRule,
activity = nothing,
translation = NoTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
verbose = false,
reference_state = nothing)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
k
: Pair Parameter (Float64
) (optional)l
: Pair Parameter (Float64
) (optional)
Model Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
a
: Pair Parameter (Float64
)b
: Pair Parameter (Float64
)
Input models
idealmodel
: Ideal Modelalpha
: Alpha modelmixing
: Mixing modelactivity
: Activity Model, used in the creation of the mixing model.translation
: Translation Model
Description
Redlich-Kwong Equation of state.
P = RT/(V-Nb) + a•α(T)/(V(V+Nb))
Model Construction Examples
# Using the default database
model = RK("water") #single input
model = RK(["water","ethanol"]) #multiple components
model = RK(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
model = RK(["water","ethanol"],alpha = Soave2019) #modifying alpha function
model = RK(["water","ethanol"],translation = RackettTranslation) #modifying translation
model = RK(["water","ethanol"],mixing = KayRule) #using another mixing rule
model = RK(["water","ethanol"],mixing = WSRule, activity = NRTL) #using advanced EoS+gᴱ mixing rule
# Passing a prebuilt model
my_alpha = SoaveAlpha(["ethane","butane"],userlocations = Dict(:acentricfactor => [0.1,0.2]))
model = RK(["ethane","butane"],alpha = my_alpha) #this is efectively now an SRK model
# User-provided parameters, passing files or folders
# Passing files or folders
model = RK(["neon","hydrogen"]; userlocations = ["path/to/my/db","cubic/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = RK(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
Mw = [20.17, 2.],
acentricfactor = [-0.03,-0.21]
k = [0. 0.18; 0.18 0.], #k,l can be ommited in single-component models.
l = [0. 0.01; 0.01 0.])
)
References
- Redlich, O., & Kwong, J. N. S. (1949). On the thermodynamics of solutions; an equation of state; fugacities of gaseous solutions. Chemical Reviews, 44(1), 233–244. doi:10.1021/cr60137a013
Clapeyron.PR
— TypePR(components;
idealmodel = BasicIdeal,
alpha = PRAlpha,
mixing = vdW1fRule,
activity = nothing,
translation = NoTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
verbose = false,
reference_state = nothing)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
k
: Pair Parameter (Float64
) (optional)l
: Pair Parameter (Float64
) (optional)
Model Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
a
: Pair Parameter (Float64
)b
: Pair Parameter (Float64
)
Input models
idealmodel
: Ideal Modelalpha
: Alpha modelmixing
: Mixing modelactivity
: Activity Model, used in the creation of the mixing model.translation
: Translation Model
Description
Peng-Robinson Equation of state.
P = RT/(V-Nb) + a•α(T)/(V-Nb₁)(V-Nb₂)
b₁ = (1 + √2)b
b₂ = (1 - √2)b
Model Construction Examples
# Using the default database
model = PR("water") #single input
model = PR(["water","ethanol"]) #multiple components
model = PR(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
model = PR(["water","ethanol"],alpha = Soave2019) #modifying alpha function
model = PR(["water","ethanol"],translation = RackettTranslation) #modifying translation
model = PR(["water","ethanol"],mixing = KayRule) #using another mixing rule
model = PR(["water","ethanol"],mixing = WSRule, activity = NRTL) #using advanced EoS+gᴱ mixing rule
# Passing a prebuilt model
my_alpha = PR78Alpha(["ethane","butane"],userlocations = Dict(:acentricfactor => [0.1,0.2]))
model = PR(["ethane","butane"],alpha = my_alpha)
# User-provided parameters, passing files or folders
model = PR(["neon","hydrogen"]; userlocations = ["path/to/my/db","cubic/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = PR(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
Mw = [20.17, 2.],
acentricfactor = [-0.03,-0.21]
k = [0. 0.18; 0.18 0.], #k,l can be ommited in single-component models.
l = [0. 0.01; 0.01 0.])
)
References
- Peng, D.Y., & Robinson, D.B. (1976). A New Two-Constant Equation of State. Industrial & Engineering Chemistry Fundamentals, 15, 59-64. doi:10.1021/I160057A011
Clapeyron.RKPR
— TypeRKPR(components;
idealmodel = BasicIdeal,
alpha = RKPRAlpha,
mixing = vdW1fRule,
activity = nothing,
translation = NoTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Critical Volume[m3/mol]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
k
: Pair Parameter (Float64
) (optional)l
: Pair Parameter (Float64
) (optional)
Model Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Critical Volume[m3/mol]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
a
: Pair Parameter (Float64
)b
: Pair Parameter (Float64
)c
: Pair Parameter (Float64
)
Input models
idealmodel
: Ideal Modelalpha
: Alpha modelmixing
: Mixing modelactivity
: Activity Model, used in the creation of the mixing model.translation
: Translation Model
Description
Redlich-Kwong-Peng-Robinson Equation of state.
P = RT/(v-b) + a•α(T)/((v + Δ₁b)*(v + Δ₂b))
Δ₁ = δ
Δ₂ = (1 - δ)/(1 + δ)
δ = ∑cᵢxᵢ
aᵢᵢ = Ωaᵢ(R²Tcᵢ²/Pcᵢ)
bᵢᵢ = Ωbᵢ(R²Tcᵢ/Pcᵢ)
Ωaᵢ = (3*yᵢ^2 + 3yᵢdᵢ + dᵢ^2 + dᵢ - 1)/(3yᵢ + dᵢ - 1)^2
Ωbᵢ = 1/(3yᵢ + dᵢ - 1)
dᵢ = (1 + cᵢ^2)/(1 + cᵢ)
yᵢ = 1 + (2(1 + cᵢ))^(1/3) + (4/(1 + cᵢ))^(1/3)
cᵢ
is fitted to match:
if Zcᵢ[exp] > 0.29
cᵢ = √2 - 1
else
Zcᵢ = 1.168Zcᵢ[exp]
f(cᵢ) == 0
f(cᵢ) = Zcᵢ - yᵢ/(3yᵢ + dᵢ - 1)
Model Construction Examples
# Using the default database
model = RKPR("water") #single input
model = RKPR(["water","ethanol"]) #multiple components
model = RKPR(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
model = RKPR(["water","ethanol"],alpha = TwuAlpha) #modifying alpha function
model = RKPR(["water","ethanol"],translation = RackettTranslation) #modifying translation
model = RKPR(["water","ethanol"],mixing = KayRule) #using another mixing rule
model = RKPR(["water","ethanol"],mixing = WSRule, activity = NRTL) #using advanced EoS+gᴱ mixing rule
# Passing a prebuilt model
my_alpha = PR78Alpha(["ethane","butane"],userlocations = Dict(:acentricfactor => [0.1,0.2]))
model = RKPR(["ethane","butane"],alpha = my_alpha)
# User-provided parameters, passing files or folders
model = RKPR(["neon","hydrogen"]; userlocations = ["path/to/my/db","cubic/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = RKPR(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
Vc = [4.25e-5, 6.43e-5],
Mw = [20.17, 2.],
acentricfactor = [-0.03,-0.21]
k = [0. 0.18; 0.18 0.], #k,l can be ommited in single-component models.
l = [0. 0.01; 0.01 0.])
)
References
- Cismondi, M., & Mollerup, J. (2005). Development and application of a three-parameter RK–PR equation of state. Fluid Phase Equilibria, 232(1–2), 74–89. doi:10.1016/j.fluid.2005.03.020
- Tassin, N. G., Mascietti, V. A., & Cismondi, M. (2019). Phase behavior of multicomponent alkane mixtures and evaluation of predictive capacity for the PR and RKPR EoS’s. Fluid Phase Equilibria, 480, 53–65. doi:10.1016/j.fluid.2018.10.005
Clapeyron.PatelTeja
— TypePatelTeja(components;
idealmodel = BasicIdeal,
alpha = NoAlpha,
mixing = vdW1fRule,
activity = nothing,
translation = PatelTejaTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Critical Volume[m3/mol]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
k
: Pair Parameter (Float64
) (optional)l
: Pair Parameter (Float64
) (optional)
Model Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Critical Volume[m3/mol]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
a
: Pair Parameter (Float64
)b
: Pair Parameter (Float64
)c
: Pair Parameter (Float64
)
Description
Patel-Teja Equation of state.
P = RT/(v-b) + a•α(T)/((v - Δ₁b)*(v - Δ₂b))
aᵢᵢ = Ωaᵢ(R²Tcᵢ²/Pcᵢ)
bᵢᵢ = Ωbᵢ(R²Tcᵢ/Pcᵢ)
cᵢ = Ωcᵢ(R²Tcᵢ/Pcᵢ)
Zcᵢ = Pcᵢ*Vcᵢ/(R*Tcᵢ)
Ωaᵢ = 3Zcᵢ² + 3(1 - 2Zcᵢ)Ωbᵢ + Ωbᵢ² + 1 - 3Zcᵢ
0 = -Zcᵢ³ + (3Zcᵢ²)*Ωbᵢ + (2 - 3Zcᵢ)*Ωbᵢ² + Ωbᵢ³
Ωcᵢ = 1 - 3Zcᵢ
γ = ∑cᵢxᵢ/∑bᵢxᵢ
δ = 1 + 6γ + γ²
ϵ = 1 + γ
Δ₁ = -(ϵ + √δ)/2
Δ₂ = -(ϵ - √δ)/2
Model Construction Examples
# Using the default database
model = PatelTeja("water") #single input
model = PatelTeja(["water","ethanol"]) #multiple components
model = PatelTeja(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
model = PatelTeja(["water","ethanol"],alpha = TwuAlpha) #modifying alpha function
model = PatelTeja(["water","ethanol"],translation = RackettTranslation) #modifying translation
model = PatelTeja(["water","ethanol"],mixing = KayRule) #using another mixing rule
model = PatelTeja(["water","ethanol"],mixing = WSRule, activity = NRTL) #using advanced EoS+gᴱ mixing rule
# Passing a prebuilt model
my_alpha = PR78Alpha(["ethane","butane"],userlocations = Dict(:acentricfactor => [0.1,0.2]))
model = PatelTeja(["ethane","butane"],alpha = my_alpha)
# User-provided parameters, passing files or folders
model = PatelTeja(["neon","hydrogen"]; userlocations = ["path/to/my/db","cubic/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = PatelTeja(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
Vc = [4.25e-5, 6.43e-5],
Mw = [20.17, 2.],
acentricfactor = [-0.03,-0.21]
k = [0. 0.18; 0.18 0.], #k,l can be ommited in single-component models.
l = [0. 0.01; 0.01 0.])
)
References
- Patel, N. C., & Teja, A. S. (1982). A new cubic equation of state for fluids and fluid mixtures. Chemical Engineering Science, 37(3), 463–473. doi:10.1016/0009-2509(82)80099-7
Clapeyron.KU
— TypeKU(components;
idealmodel = BasicIdeal,
alpha = KUAlpha,
mixing = vdW1fRule,
activity = nothing,
translation = NoTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Critical Volume[m^3]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
k
: Pair Parameter (Float64
) (optional)l
: Pair Parameter (Float64
) (optional)
Model Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Critical Volume[m^3]
omega_a
: Single Parameter (Float64
) - Critical Constant for a - No unitsomega_b
: Single Parameter (Float64
) - Critical Constant for b - No unitsMw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
a
: Pair Parameter (Float64
)b
: Pair Parameter (Float64
)
Input models
idealmodel
: Ideal Modelalpha
: Alpha modelmixing
: Mixing modelactivity
: Activity Model, used in the creation of the mixing model.translation
: Translation Model
Description
Kumar-Upadhyay Cubic Equation of state. Ωa
and Ωb
are component-dependent
P = RT/(v-b) + a•κ(T)/((v²-1.6bv - 0.8b²)
a = Ωa(R²Tcᵢ²/Pcᵢ)
b = Ωb(R²Tcᵢ²/Pcᵢ)
Ωa = Zc[(1 + 1.6α - 0.8α²)²/((1 - α²)(2 + 1.6α))]
χ = ∛[√(1458Zc³ - 1701Zc² + 540Zc -20)/32√3Zc² - (729Zc³ - 216Zc + 8)/1728Zc³]
α = [χ + (81Zc² - 72Zc + 4)/144Zc²χ + (3Zc - 2)/12Zc]
Ωa = Zc[(1 + 1.6α - 0.8α²)²/((1 - α²)(2 + 1.6α))]
Ωb = αZc
Model Construction Examples
# Using the default database
model = KU("water") #single input
model = KU(["water","ethanol"]) #multiple components
model = KU(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
model = KU(["water","ethanol"],alpha = TwuAlpha) #modifying alpha function
model = KU(["water","ethanol"],translation = RackettTranslation) #modifying translation
model = KU(["water","ethanol"],mixing = KayRule) #using another mixing rule
model = KU(["water","ethanol"],mixing = WSRule, activity = NRTL) #using advanced EoS+gᴱ mixing rule
# Passing a prebuilt model
my_alpha = PR78Alpha(["ethane","butane"],userlocations = Dict(:acentricfactor => [0.1,0.2]))
model = KU(["ethane","butane"],alpha = my_alpha)
# User-provided parameters, passing files or folders
model = KU(["neon","hydrogen"]; userlocations = ["path/to/my/db","cubic/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = KU(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
Vc = [4.25e-5, 6.43e-5],
Mw = [20.17, 2.],
acentricfactor = [-0.03,-0.21]
k = [0. 0.18; 0.18 0.], #k,l can be ommited in single-component models.
l = [0. 0.01; 0.01 0.])
)
References
- Kumar, A., & Upadhyay, R. (2021). A new two-parameters cubic equation of state with benefits of three-parameters. Chemical Engineering Science, 229(116045), 116045. doi:10.1016/j.ces.2020.116045
Variant Models
Clapeyron.Berthelot
— TypeBerthelot(components;
idealmodel = BasicIdeal,
alpha = ClausiusAlpha,
mixing = vdW1fRule,
activity = nothing,
translation = NoTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
Vc
: Single Parameter (Float64
) - Molar Volume[m^3/mol]
k
: Pair Parameter (Float64
) (optional)
Model Parameters
a
: Pair Parameter (Float64
)b
: Pair Parameter (Float64
)Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Molar Volume[m^3/mol]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
Input models
idealmodel
: Ideal Modelalpha
: Alpha modelmixing
: Mixing modelactivity
: Activity Model, used in the creation of the mixing model.translation
: Translation Model
Description
Berthelot Equation of state. it uses the Volume-Pressure Based mixing rules, that is:
a = 8*Pc*Vc^2
b = Vc/3
R = (8/3)*Pc*Vc/Tc
P = RT/(V-Nb) + a•α(T)/V²
α(T) = Tc/T
References
- Berthelot, D. (1899). Sur une méthode purement physique pour la détermination des poids moléculaires des gaz et des poids atomiques de leurs éléments. Journal de Physique Théorique et Appliquée, 8(1), 263–274. doi:10.1051/jphystap:018990080026300
Clapeyron.SRK
— FunctionSRK(components::Vector{String};
idealmodel = BasicIdeal,
alpha = SoaveAlpha,
mixing = vdW1fRule,
activity = nothing,
translation = NoTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Description
Soave-Redlich-Kwong equation of state. it uses the following models:
- Translation Model:
NoTranslation
- Alpha Model:
SoaveAlpha
- Mixing Rule Model:
vdW1fRule
Model Construction Examples
# Using the default database
model = SRK("water") #single input
model = SRK(["water","ethanol"]) #multiple components
model = SRK(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
model = SRK(["water","ethanol"],alpha = Soave2019) #modifying alpha function, using 2019 correlation
model = SRK(["water","ethanol"],translation = RackettTranslation) #modifying translation
model = SRK(["water","ethanol"],mixing = KayRule) #using another mixing rule
model = SRK(["water","ethanol"],mixing = WSRule, activity = NRTL) #using advanced EoS+gᴱ mixing rule
# Passing a prebuilt model
my_alpha = SoaveAlpha(["ethane","butane"],userlocations = Dict(:acentricfactor => [0.1,0.2]))
model = SRK(["ethane","butane"],alpha = my_alpha)
# User-provided parameters, passing files or folders
model = SRK(["neon","hydrogen"]; userlocations = ["path/to/my/db","cubic/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = SRK(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
Mw = [20.17, 2.],
acentricfactor = [-0.03,-0.21]
k = [0. 0.18; 0.18 0.], #k,l can be ommited in single-component models.
l = [0. 0.01; 0.01 0.])
)
References
- Soave, G. (1972). Equilibrium constants from a modified Redlich-Kwong equation of state. Chemical Engineering Science, 27(6), 1197–1203. doi:10.1016/0009-2509(72)80096-4
Clapeyron.PSRK
— Functionfunction PSRK(components;
idealmodel = BasicIdeal,
alpha = SoaveAlpha,
mixing = PSRKRule,
activity = PSRKUNIFAC,
translation = PenelouxTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Description
Predictive Soave-Redlich-Kwong equation of state. it uses the following models:
- Translation Model:
NoTranslation
- Alpha Model:
SoaveAlpha
- Mixing Rule Model:
PSRKRule
withPSRKUNIFAC
activity model
References
- Horstmann, S., Jabłoniec, A., Krafczyk, J., Fischer, K., & Gmehling, J. (2005). PSRK group contribution equation of state: comprehensive revision and extension IV, including critical constants and α-function parameters for 1000 components. Fluid Phase Equilibria, 227(2), 157–164. doi:10.1016/j.fluid.2004.11.002
Clapeyron.tcRK
— FunctiontcRK(components::Vector{String};
idealmodel = BasicIdeal,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false,
estimate_alpha = true,
estimate_translation = true)
translated and consistent Redlich-Kwong equation of state. it uses the following models:
- Translation Model:
ConstantTranslation
- Alpha Model:
TwuAlpha
- Mixing Rule Model:
vdW1fRule
If Twu parameters are not provided, they can be estimated from the acentric factor (acentricfactor
). If translation is not provided, it can be estimated, using Rackett compresibility Factor (ZRA
) or the acentric factor (acentricfactor
).
The use of estimates for the alpha function and volume translation can be turned off by passing estimate_alpha = false
or estimate_translation = false
.
References
- Le Guennec, Y., Privat, R., & Jaubert, J.-N. (2016). Development of the translated-consistent tc-PR and tc-RK cubic equations of state for a safe and accurate prediction of volumetric, energetic and saturation properties of pure compounds in the sub- and super-critical domains. Fluid Phase Equilibria, 429, 301–312. doi:10.1016/j.fluid.2016.09.003
- Pina-Martinez, A., Le Guennec, Y., Privat, R., Jaubert, J.-N., & Mathias, P. M. (2018). Analysis of the combinations of property data that are suitable for a safe estimation of consistent twu α-function parameters: Updated parameter values for the translated-consistent tc-PR and tc-RK cubic equations of state. Journal of Chemical and Engineering Data, 63(10), 3980–3988. doi:10.1021/acs.jced.8b00640
- Piña-Martinez, A., Privat, R., & Jaubert, J.-N. (2022). Use of 300,000 pseudo‐experimental data over 1800 pure fluids to assess the performance of four cubic equations of state: SRK , PR , tc ‐RK , and tc ‐PR. AIChE Journal. American Institute of Chemical Engineers, 68(2). doi:10.1002/aic.17518
Clapeyron.PR78
— FunctionPR78(components::Vector{String};
idealmodel = BasicIdeal,
alpha = PR78Alpha,
mixing = vdW1fRule,
activity = nothing,
translation = NoTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Peng Robinson (1978) equation of state. it uses the following models:
- Translation Model:
NoTranslation
- Alpha Model:
PR78Alpha
- Mixing Rule Model:
vdW1fRule
Model Construction Examples
# Using the default database
model = PR78("water") #single input
model = PR78(["water","ethanol"]) #multiple components
model = PR78(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
model = PR78(["water","ethanol"],alpha = Soave2019) #modifying alpha function
model = PR78(["water","ethanol"],translation = RackettTranslation) #modifying translation
model = PR78(["water","ethanol"],mixing = KayRule) #using another mixing rule
model = PR78(["water","ethanol"],mixing = WSRule, activity = NRTL) #using advanced EoS+gᴱ mixing rule
# Passing a prebuilt model
my_alpha = PRAlpha(["ethane","butane"],userlocations = Dict(:acentricfactor => [0.1,0.2]))
model = PR78(["ethane","butane"],alpha = my_alpha) #this model becomes a normal PR EoS
# User-provided parameters, passing files or folders
model = PR78(["neon","hydrogen"]; userlocations = ["path/to/my/db","cubic/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = PR78(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
Mw = [20.17, 2.],
acentricfactor = [-0.03,-0.21]
k = [0. 0.18; 0.18 0.], #k,l can be ommited in single-component models.
l = [0. 0.01; 0.01 0.])
)
References
- Robinson DB, Peng DY. The characterization of the heptanes and heavier fractions for the GPA Peng-Robinson programs. Tulsa: Gas Processors Association; 1978
Clapeyron.PTV
— TypePTV(components;
idealmodel = BasicIdeal,
alpha = NoAlpha,
mixing = vdW1fRule,
activity = nothing,
translation = PTVTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
verbose = false,
reference_state = nothing)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Critical Volume[m3/mol]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
k
: Pair Parameter (Float64
) (optional)l
: Pair Parameter (Float64
) (optional)
Model Parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) - Critical Volume[m3/mol]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
a
: Pair Parameter (Float64
)b
: Pair Parameter (Float64
)c
: Pair Parameter (Float64
)
Description
Patel-Teja-Valderrama Equation of state.
P = RT/(v-b) + a•α(T)/((v - Δ₁b)*(v - Δ₂b))
aᵢᵢ = Ωaᵢ(R²Tcᵢ²/Pcᵢ)
bᵢᵢ = Ωbᵢ(R²Tcᵢ/Pcᵢ)
cᵢ = Ωcᵢ(R²Tcᵢ/Pcᵢ)
Zcᵢ = Pcᵢ*Vcᵢ/(R*Tcᵢ)
Ωaᵢ = 0.66121 - 0.76105Zcᵢ
Ωbᵢ = 0.02207 + 0.20868Zcᵢ
Ωcᵢ = 0.57765 - 1.87080Zcᵢ
γ = ∑cᵢxᵢ/∑bᵢxᵢ
δ = 1 + 6γ + γ²
ϵ = 1 + γ
Δ₁ = -(ϵ + √δ)/2
Δ₂ = -(ϵ - √δ)/2
Model Construction Examples
# Using the default database
model = PTV("water") #single input
model = PTV(["water","ethanol"]) #multiple components
model = PTV(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
model = PTV(["water","ethanol"],alpha = TwuAlpha) #modifying alpha function
model = PTV(["water","ethanol"],translation = RackettTranslation) #modifying translation
model = PTV(["water","ethanol"],mixing = KayRule) #using another mixing rule
model = PTV(["water","ethanol"],mixing = WSRule, activity = NRTL) #using advanced EoS+gᴱ mixing rule
# Passing a prebuilt model
my_alpha = PR78Alpha(["ethane","butane"],userlocations = Dict(:acentricfactor => [0.1,0.2]))
model = PTV(["ethane","butane"],alpha = my_alpha)
# User-provided parameters, passing files or folders
model = PTV(["neon","hydrogen"]; userlocations = ["path/to/my/db","cubic/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = PTV(["neon","hydrogen"];
userlocations = (;Tc = [44.492,33.19],
Pc = [2679000, 1296400],
Vc = [4.25e-5, 6.43e-5],
Mw = [20.17, 2.],
acentricfactor = [-0.03,-0.21]
k = [0. 0.18; 0.18 0.], #k,l can be ommited in single-component models.
l = [0. 0.01; 0.01 0.])
)
References
- Valderrama, J. O. (1990). A generalized Patel-Teja equation of state for polar and nonpolar fluids and their mixtures. Journal of Chemical Engineering of Japan, 23(1), 87–91. doi:10.1252/jcej.23.87
Clapeyron.EPPR78
— FunctionEPPR78(components_or_groups;
idealmodel = BasicIdeal,
alpha = PR78Alpha,
mixing = PPR78Rule,
activity = nothing,
translation = NoTranslation,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Enhanced Predictive Peng Robinson equation of state. it uses the following models:
- Translation Model:
NoTranslation
- Alpha Model:
PR78Alpha
- Mixing Rule Model:
PPR78Rule
References
- Jaubert, J.-N., Privat, R., & Mutelet, F. (2010). Predicting the phase equilibria of synthetic petroleum fluids with the PPR78 approach. AIChE Journal. American Institute of Chemical Engineers, 56(12), 3225–3235. doi:10.1002/aic.12232
- Jaubert, J.-N., Qian, J.-W., Lasala, S., & Privat, R. (2022). The impressive impact of including enthalpy and heat capacity of mixing data when parameterising equations of state. Application to the development of the E-PPR78 (Enhanced-Predictive-Peng-Robinson-78) model. Fluid Phase Equilibria, (113456), 113456. doi:10.1016/j.fluid.2022.113456
Clapeyron.UMRPR
— FunctionUMRPR(components;
idealmodel = BasicIdeal,
userlocations = String[],
group_userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Universal Mixing Rule Peng Robinson equation of state. it uses the following models:
- Translation Model:
MTTranslation
- Alpha Model:
MTAlpha
- Mixing Rule Model:
UMRRule
withUNIFAC
activity
References
- Voutsas, E., Magoulas, K., & Tassios, D. (2004). Universal mixing rule for cubic equations of state applicable to symmetric and asymmetric systems: Results with the Peng−Robinson equation of state. Industrial & Engineering Chemistry Research, 43(19), 6238–6246. doi:10.1021/ie049580p
Clapeyron.VTPR
— FunctionVTPR(components;
idealmodel = BasicIdeal,
userlocations = String[],
group_userlocations = String[]
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Volume-translated Peng Robinson equation of state. it uses the following models:
- Translation Model:
RackettTranslation
- Alpha Model:
TwuAlpha
- Mixing Rule Model:
VTPRRule
withVTPRUNIFAC
activity
References
- Ahlers, J., & Gmehling, J. (2001). Development of an universal group contribution equation of state. Fluid Phase Equilibria, 191(1–2), 177–188. doi:10.1016/s0378-3812(01)00626-4
Clapeyron.tcPR
— FunctiontcPR(components;
idealmodel = BasicIdeal,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false,
estimate_alpha = true,
estimate_translation = true)
translated and consistent Peng Robinson equation of state. it uses the following models:
- Translation Model:
ConstantTranslation
- Alpha Model:
TwuAlpha
- Mixing Rule Model:
vdW1fRule
If Twu parameters are not provided, they can be estimated from the acentric factor (acentricfactor
). If translation is not provided, it can be estimated, using Rackett compresibility Factor (ZRA
) or the acentric factor (acentricfactor
).
References
- Le Guennec, Y., Privat, R., & Jaubert, J.-N. (2016). Development of the translated-consistent tc-PR and tc-RK cubic equations of state for a safe and accurate prediction of volumetric, energetic and saturation properties of pure compounds in the sub- and super-critical domains. Fluid Phase Equilibria, 429, 301–312. doi:10.1016/j.fluid.2016.09.003
- Pina-Martinez, A., Le Guennec, Y., Privat, R., Jaubert, J.-N., & Mathias, P. M. (2018). Analysis of the combinations of property data that are suitable for a safe estimation of consistent twu α-function parameters: Updated parameter values for the translated-consistent tc-PR and tc-RK cubic equations of state. Journal of Chemical and Engineering Data, 63(10), 3980–3988. doi:10.1021/acs.jced.8b00640
- Piña-Martinez, A., Privat, R., & Jaubert, J.-N. (2022). Use of 300,000 pseudo‐experimental data over 1800 pure fluids to assess the performance of four cubic equations of state: SRK , PR , tc ‐RK , and tc ‐PR. AIChE Journal. American Institute of Chemical Engineers, 68(2). doi:10.1002/aic.17518
Clapeyron.tcPRW
— FunctiontcPR(components;
idealmodel = BasicIdeal,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
translated and consistent Peng Robinson equation of state,with an gE mixing rule. it uses the following models:
- Translation Model:
ConstantTranslation
- Alpha Model:
TwuAlpha
- Mixing Rule Model:
gErRule
- Activity:
tcPRWilsonRes
If Twu parameters are not provided, they can be estimated from the acentric factor (acentricfactor
). If translation is not provided, it can be estimated, using Rackett compresibility Factor (ZRA
) or the acentric factor (acentricfactor
).
References
- Le Guennec, Y., Privat, R., & Jaubert, J.-N. (2016). Development of the translated-consistent tc-PR and tc-RK cubic equations of state for a safe and accurate prediction of volumetric, energetic and saturation properties of pure compounds in the sub- and super-critical domains. Fluid Phase Equilibria, 429, 301–312. doi:10.1016/j.fluid.2016.09.003
- Pina-Martinez, A., Le Guennec, Y., Privat, R., Jaubert, J.-N., & Mathias, P. M. (2018). Analysis of the combinations of property data that are suitable for a safe estimation of consistent twu α-function parameters: Updated parameter values for the translated-consistent tc-PR and tc-RK cubic equations of state. Journal of Chemical and Engineering Data, 63(10), 3980–3988. doi:10.1021/acs.jced.8b00640
- Piña-Martinez, A., Privat, R., & Jaubert, J.-N. (2022). Use of 300,000 pseudo‐experimental data over 1800 pure fluids to assess the performance of four cubic equations of state: SRK , PR , tc ‐RK , and tc ‐PR. AIChE Journal. American Institute of Chemical Engineers, 68(2). doi:10.1002/aic.17518
- Piña-Martinez, A., Privat, R., Nikolaidis, I. K., Economou, I. G., & Jaubert, J.-N. (2021). What is the optimal activity coefficient model to be combined with the translated–consistent Peng–Robinson equation of state through advanced mixing rules? Cross-comparison and grading of the Wilson, UNIQUAC, and NRTL aE models against a benchmark database involving 200 binary systems. Industrial & Engineering Chemistry Research, 60(47), 17228–17247. doi:10.1021/acs.iecr.1c03003
Clapeyron.cPR
— FunctioncPR(components::Vector{String};
idealmodel = BasicIdeal,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
consistent Peng Robinson equation of state. it uses the following models:
- Translation Model:
NoTranslation
- Alpha Model:
TwuAlpha
- Mixing Rule Model:
vdW1fRule
References
- Bell, I. H., Satyro, M., & Lemmon, E. W. (2018). Consistent twu parameters for more than 2500 pure fluids from critically evaluated experimental data. Journal of Chemical and Engineering Data, 63(7), 2402–2409. doi:10.1021/acs.jced.7b00967
Clapeyron.QCPR
— FunctionQCPR(components;
idealmodel = BasicIdeal,
userlocations = String[],
ideal_userlocations = String[],
alpha_userlocations = String[],
mixing_userlocations = String[],
activity_userlocations = String[],
translation_userlocations = String[],
reference_state = nothing,
verbose = false)
Quantum-corrected Peng Robinson equation of state. it uses the following models:
- Translation Model:
ConstantTranslation
- Alpha Model:
TwuAlpha
- Mixing Rule Model:
QCPRRule
References
- Aasen, A., Hammer, M., Lasala, S., Jaubert, J.-N., & Wilhelmsen, Ø. (2020). Accurate quantum-corrected cubic equations of state for helium, neon, hydrogen, deuterium and their mixtures. Fluid Phase Equilibria, 524(112790), 112790. doi:10.1016/j.fluid.2020.112790
Alpha (α(T))
Models
Clapeyron.α_function
— Functionα_function(model::CubicModel,V,T,z,αmodel::AlphaModel)
Interface function used in cubic models. it should return a vector of αᵢ(T).
Example:
function α_function(model::CubicModel,V,T,z,alpha_model::RKAlphaModel)
return 1 ./ sqrt.(T ./ Tc)
end
Clapeyron.NoAlpha
— TypeNoAlpha(args...) <: NoAlphaModel
Input Parameters
None
Description
Cubic alpha (α(T))
model. Default for vdW
EoS
αᵢ = 1 ∀ i
Model Construction Examples
# Because this model does not have parameters, all those constructors are equivalent:
alpha = NoAlpha()
alpha = NoAlpha("water")
alpha = NoAlpha(["water","carbon dioxide"])
Clapeyron.ClausiusAlpha
— TypeClausiusAlpha <: ClausiusAlphaModel
ClausiusAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
- none
Description
Cubic alpha (α(T))
model. Default for Clausius
and [Berthelot
]
αᵢ = 1/Trᵢ
Trᵢ = T/Tcᵢ
Model Construction Examples
# Because this model does not have parameters, all those constructors are equivalent:
alpha = ClausiusAlpha()
alpha = ClausiusAlpha("water")
alpha = ClausiusAlpha(["water","carbon dioxide"])
Clapeyron.RKAlpha
— TypeRKAlpha <: RKAlphaModel
RKAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
- none
Description
Cubic alpha (α(T))
model. Default for RK
EoS.
αᵢ = 1/√(Trᵢ)
Trᵢ = T/Tcᵢ
Model Construction Examples
# Because this model does not have parameters, all those constructors are equivalent:
alpha = RKAlpha()
alpha = RKAlpha("water")
alpha = RKAlpha(["water","carbon dioxide"])
Clapeyron.SoaveAlpha
— TypeSoaveAlpha <: SoaveAlphaModel
SoaveAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
acentricfactor
: Single Parameter (Float64
)
Description
Cubic alpha (α(T))
model. Default for SRK
EoS.
αᵢ = (1+mᵢ(1-√(Trᵢ)))^2
Trᵢ = T/Tcᵢ
mᵢ = 0.480 + 1.547ωᵢ - 0.176ωᵢ^2
To use different polynomial coefficients for mᵢ
, overload Clapeyron.α_m(::CubicModel,::SoaveAlphaModel) = (c₁,c₂,...cₙ)
Model Construction Examples
# Using the default database
alpha = SoaveAlpha("water") #single input
alpha = SoaveAlpha(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = SoaveAlpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/acentric.csv"])
# Passing parameters directly
alpha = SoaveAlpha(["neon","hydrogen"];userlocations = (;acentricfactor = [-0.03,-0.21]))
Clapeyron.Soave2019Alpha
— TypeSoave2019Alpha <: SoaveAlphaModel
Soave2019Alpha(components::Vector{String};
userlocations = String[],
verbose::Bool=false)
Input Parameters
acentricfactor
: Single Parameter (Float64
)
Description
Cubic alpha (α(T))
model. Updated m(ω) correlations for PR
and SRK
with better results for heavy molecules.
αᵢ = (1+mᵢ(1-√(Trᵢ)))^2
Trᵢ = T/Tcᵢ
mᵢ = 0.37464 + 1.54226ωᵢ - 0.26992ωᵢ^2
where, for Peng-robinson:
mᵢ = 0.3919 + 1.4996ωᵢ - 0.2721ωᵢ^2 + 0.1063ωᵢ^3
and, for Redlich-Kwong:
mᵢ = 0.4810 + 1.5963ωᵢ - 0.2963ωᵢ^2 + 0.1223ωᵢ^3
Model Construction Examples
# Using the default database
alpha = Soave2019Alpha("water") #single input
alpha = Soave2019Alpha(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = Soave2019Alpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/acentric.csv"])
# Passing parameters directly
alpha = Soave2019Alpha(["neon","hydrogen"];userlocations = (;acentricfactor = [-0.03,-0.21]))
References
- Pina-Martinez, A., Privat, R., Jaubert, J.-N., & Peng, D.-Y. (2019). Updated versions of the generalized Soave α-function suitable for the Redlich-Kwong and Peng-Robinson equations of state. Fluid Phase Equilibria, 485, 264–269. doi:10.1016/j.fluid.2018.12.007
Clapeyron.PRAlpha
— TypePRAlpha <: SoaveAlphaModel
PRAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
acentricfactor
: Single Parameter (Float64
)
Description
Cubic alpha (α(T))
model. Default for PR
EoS.
αᵢ = (1+mᵢ(1-√(Trᵢ)))^2
Trᵢ = T/Tcᵢ
mᵢ = 0.37464 + 1.54226ωᵢ - 0.26992ωᵢ^2
It is equivalent to SoaveAlpha
.
Model Construction Examples
# Using the default database
alpha = PRAlpha("water") #single input
alpha = PRAlpha(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = PRAlpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/acentric.csv"])
# Passing parameters directly
alpha = PRAlpha(["neon","hydrogen"];userlocations = (;acentricfactor = [-0.03,-0.21]))
Clapeyron.PR78Alpha
— TypePR78Alpha <: PR78AlphaModel
PR78Alpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
acentricfactor
: Single Parameter (Float64
)
Description
Cubic alpha (α(T))
model. Default for PR78
and EPPR78
EoS.
αᵢ = (1+mᵢ(1-√(Trᵢ)))^2
Trᵢ = T/Tcᵢ
if ωᵢ ≤ 0.491
mᵢ = 0.37464 + 1.54226ωᵢ - 0.26992ωᵢ^2
else
mᵢ = 0.379642 + 1.487503ωᵢ - 0.164423ωᵢ^2 - 0.016666ωᵢ^3
Model Construction Examples
# Using the default database
alpha = PR78Alpha("water") #single input
alpha = PR78Alpha(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = PR78Alpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/acentric.csv"])
# Passing parameters directly
alpha = PR78Alpha(["neon","hydrogen"];userlocations = (;acentricfactor = [-0.03,-0.21]))
Clapeyron.CPAAlpha
— TypeCPAAlpha <: CPAAlphaModel
CPAAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
c1
: Single Parameter (Float64
)
Description
Cubic alpha (α(T))
model. Default for CPA
EoS.
αᵢ = (1+c¹ᵢ(1-√(Trᵢ)))^2
Model Construction Examples
# Using the default database
alpha = CPAAlpha("water") #single input
alpha = CPAAlpha(["water","carbon dioxide"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = CPAAlpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","cpa/alpha.csv"])
# Passing parameters directly
alpha = CPAAlpha(["water","carbon dioxide"];userlocations = (;c1 = [0.67,0.76]))
Clapeyron.sCPAAlpha
— TypesCPAAlpha <: sCPAAlphaModel
sCPAAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
c1
: Single Parameter
Description
Cubic alpha (α(T))
model. Default for sCPA
EoS.
αᵢ = (1+c¹ᵢ(1-√(Trᵢ)))^2
Model Construction Examples
# Using the default database
alpha = sCPAAlpha("water") #single input
alpha = sCPAAlpha(["water","carbon dioxide"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = sCPAAlpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","scpa/alpha.csv"])
# Passing parameters directly
alpha = sCPAAlpha(["water","carbon dioxide"];userlocations = (;c1 = [0.67,0.76]))
Clapeyron.MTAlpha
— TypeMTAlpha <: MTAlphaModel
MTAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
acentricfactor
: Single Parameter (Float64
)
Description
Magoulas & Tassios Cubic alpha (α(T))
model. Default for UMRPR
EoS.
αᵢ = (1+mᵢ(1-√(Trᵢ)))^2
Trᵢ = T/Tcᵢ
mᵢ = 0.384401 + 1.52276ωᵢ - 0.213808ωᵢ^2 + 0.034616ωᵢ^3 - 0.001976ωᵢ^4
Model Construction Examples
# Using the default database
alpha = MTAlpha("water") #single input
alpha = MTAlpha(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = MTAlpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/acentric.csv"])
# Passing parameters directly
alpha = MTAlpha(["neon","hydrogen"];userlocations = (;acentricfactor = [-0.03,-0.21]))
References
- Magoulas, K., & Tassios, D. (1990). Thermophysical properties of n-Alkanes from C1 to C20 and their prediction for higher ones. Fluid Phase Equilibria, 56, 119–140. doi:10.1016/0378-3812(90)85098-u
Clapeyron.BMAlpha
— TypeBMAlpha <: BMAlphaModel
BMAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
acentricfactor
: Single Parameter (Float64
)
Description
Boston Mathias Cubic alpha (α(T))
model.
if Trᵢ > 1
αᵢ = (exp((1-2/(2+mᵢ))*(1-Trᵢ^(1+mᵢ/2))))^2
else
αᵢ = (1+mᵢ*(1-√(Trᵢ)))^2
Trᵢ = T/Tcᵢ
for PR models:
mᵢ = 0.37464 + 1.54226ωᵢ - 0.26992ωᵢ^2
for RK models:
mᵢ = 0.480 + 1.547ωᵢ - 0.176ωᵢ^2
Model Construction Examples
# Using the default database
alpha = BMAlpha("water") #single input
alpha = BMAlpha(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = BMAlpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/acentric.csv"])
# Passing parameters directly
alpha = BMAlpha(["neon","hydrogen"];userlocations = (;acentricfactor = [-0.03,-0.21]))
References
- .M. Boston, P.M. Mathias, Proceedings of the 2nd International Conference on Phase Equilibria and Fluid Properties in the Chemical Process Industries, West Berlin, March, 1980, pp. 823–849
Clapeyron.TwuAlpha
— TypeTwuAlpha <: TwuAlphaModel
Twu91Alpha = TwuAlpha
TwuAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
M
: Single ParameterN
: Single ParameterL
: Single Parameter
Description
Cubic alpha (α(T))
model. Default for VTPR
EoS. Also known as Twu-91 alpha
αᵢ = Trᵢ^(N*(M-1))*exp(L*(1-Trᵢ^(N*M))
Trᵢ = T/Tcᵢ
Model Construction Examples
# Using the default database
alpha = TwuAlpha("water") #single input
alpha = Twu91Alpha("water") #same function
alpha = TwuAlpha(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = TwuAlpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","twu.csv"])
# Passing parameters directly
alpha = TwuAlpha(["neon","hydrogen"];
userlocations = (;L = [0.40453, 156.21],
M = [0.95861, -0.0062072],
N = [0.8396, 5.047])
)
References
- Twu, C. H., Lee, L. L., & Starling, K. E. (1980). Improved analytical representation of argon thermodynamic behavior. Fluid Phase Equilibria, 4(1–2), 35–44. doi:10.1016/0378-3812(80)80003-3
Clapeyron.Twu88Alpha
— FunctionTwu88Alpha::TwuAlpha
Twu88Alpha(components::Vector{String};
userlocations = String[],
verbose::Bool=false)
Input Parameters
M
: Single ParameterN
: Single Parameter (optional)L
: Single Parameter
Model Parameters
M
: Single ParameterN
: Single ParameterL
: Single Parameter
Description
Cubic alpha (α(T))
model. Also known as Twu-88 alpha.
αᵢ = Trᵢ^(N*(M-1))*exp(L*(1-Trᵢ^(N*M))
N = 2
Trᵢ = T/Tcᵢ
if N
is specified, it will be used instead of the default value of 2.
Model Construction Examples
# Using the default database
alpha = Twu88Alpha("water") #single input
alpha = Twu88Alpha(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = Twu88Alpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","twu88.csv"])
# Passing parameters directly
alpha = Twu88Alpha(["neon","hydrogen"];
userlocations = (;L = [0.40453, 156.21],
M = [0.95861, -0.0062072],
N = [0.8396, 5.047]) #if we don't pass N, then is assumed N = 2
)
References
- Twu, C. H., Lee, L. L., & Starling, K. E. (1980). Improved analytical representation of argon thermodynamic behavior. Fluid Phase Equilibria, 4(1–2), 35–44. doi:10.1016/0378-3812(80)80003-3
Clapeyron.PatelTejaAlpha
— TypePatelTejaAlpha <: SoaveAlphaModel
PatelTejaAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
acentricfactor
: Single Parameter (Float64
)
Description
Cubic alpha (α(T))
model. Default for PatelTeja
EoS.
αᵢ = (1+mᵢ(1-√(Trᵢ)))^2
Trᵢ = T/Tcᵢ
mᵢ = 0.452413 + 1.30982ωᵢ - 0.295937ωᵢ^2
Model Construction Examples
# Using the default database
alpha = PatelTejaAlpha("water") #single input
alpha = PatelTejaAlpha(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = PatelTejaAlpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/acentric.csv"])
# Passing parameters directly
alpha = PatelTejaAlpha(["neon","hydrogen"];userlocations = (;acentricfactor = [-0.03,-0.21]))
Clapeyron.KUAlpha
— TypeKUAlpha <: AlphaModel
KUAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
acentricfactor
: Single Parameter (Float64
)
Description
Cubic alpha (α(T))
model. Default for KU
EoS.
For Tr < 1
αᵢ = (1+mᵢ(1-√(Trᵢ))^nᵢ)^2
Trᵢ = T/Tcᵢ
mᵢ = 0.37790 + 1.51959ωᵢ - 0.46904ωᵢ^2 + 0.015679ωᵢ^3
nᵢ = 0.97016 + 0.05495ωᵢ - 0.1293ωᵢ^2 + 0.0172028ωᵢ^3
For Tr > 1
is a 6th order taylor expansion around T = Tc
.
Model Construction Examples
# Using the default database
alpha = KUAlpha("water") #single input
alpha = KUAlpha(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = KUAlpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/acentric.csv"])
# Passing parameters directly
alpha = KUAlpha(["neon","hydrogen"];userlocations = (;acentricfactor = [-0.03,-0.21]))
References
- Kumar, A., & Upadhyay, R. (2021). A new two-parameters cubic equation of state with benefits of three-parameters. Chemical Engineering Science, 229(116045), 116045. doi:10.1016/j.ces.2020.116045
Clapeyron.RKPRAlpha
— TypeRKPRAlpha <: RKPRAlphaModel
RKPRAlpha(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
acentricfactor
: Single Parameter (Float64
)
Description
Cubic alpha (α(T))
model. Default for RKPR
EoS.
αᵢ = (3/(2 + Trᵢ))^kᵢ
kᵢ = (12.504Zc -2.7238) + (7.4513*Zc + 1.9681)ωᵢ + (-2.4407*Zc + 0.0017)ωᵢ^2
Trᵢ = T/Tcᵢ
Model Construction Examples
# Using the default database
alpha = RKPRAlpha("water") #single input
alpha = RKPRAlpha(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
alpha = RKPRAlpha(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/acentric.csv"])
# Passing parameters directly
alpha = RKPRAlpha(["neon","hydrogen"];userlocations = (;acentricfactor = [-0.03,-0.21]))
Volume Translation Models
Clapeyron.translation
— Functiontranslation(model::CubicModel,V,T,z,translation_model::TranslationModel)
Interface function used in cubic models. it should return a vector of cᵢ. such as Ṽ = V + mixing(c,z)
Example:
function translation(model::CubicModel,V,T,z,translation_model::RackettTranslation)
Tc = model.params.Tc.values
Pc = model.params.Pc.values
Vc = translation_model.params.Vc.values
R = Clapeyron.R̄
Zc = Pc .* Vc ./ (R .* Tc)
c = 0.40768 .* (0.29441 .- Zc) .* R .* Tc ./ Pc
return c
end
Clapeyron.NoTranslation
— TypeNoTranslation(args...) <: TranslationModel
Input Parameters
None
Description
Default volume translation model for cubic models. it performs no translation:
V = V₀ + mixing_rule(cᵢ)
cᵢ = 0 ∀ i
Model Construction Examples
# Because this model does not have parameters, all those constructors are equivalent:
translation = NoTranslation()
translation = NoTranslation("water")
translation = NoTranslation(["water","carbon dioxide"])
Clapeyron.ConstantTranslation
— TypeConstantTranslation <: ConstantTranslationModel
ConstantTranslation(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
v_shift
: Single Parameter (Float64
) - Volume shift[m³/mol]
Description
Constant Translation model for cubics:
V = V₀ + mixing_rule(cᵢ)
where cᵢ
is constant. It does not have parameters by default, the volume shifts must be user-supplied.
Model Construction Examples
# Using user-provided parameters
# Passing files or folders
translation = ConstantTranslation(["neon","hydrogen"]; userlocations = ["path/to/my/db","properties/critical.csv"])
# Passing parameters directly
translation = ConstantTranslation(["neon","hydrogen"];userlocations = (;Vc = [4.25e-5, 6.43e-5]))
Clapeyron.RackettTranslation
— TypeRackettTranslation <: RackettTranslationModel
RackettTranslation(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
Vc
: Single Parameter (Float64
) - Critical Volume[m³/mol]
Model Parameters
Vc
: Single Parameter (Float64
) - Critical Volume[m³/mol]
v_shift
: Single Parameter (Float64
) - Volume shift[m³/mol]
Description
Rackett Translation model for cubics:
V = V₀ + mixing_rule(cᵢ)
cᵢ = 0.40768*RTcᵢ/Pcᵢ*(0.29441-Zcᵢ)
Zcᵢ = Pcᵢ*Vcᵢ/(RTcᵢ)
Model Construction Examples
# Using the default database
translation = RackettTranslation("water") #single input
translation = RackettTranslation(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
translation = RackettTranslation(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/Vc.csv"])
# Passing parameters directly
translation = RackettTranslation(["neon","hydrogen"];userlocations = (;Vc = [4.25e-5, 6.43e-5]))
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.PenelouxTranslation
— TypePenelouxTranslation <: PenelouxTranslationModel
PenelouxTranslation(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
Vc
: Single Parameter (Float64
) - Critical Volume[m³/mol]
Model Parameters
Vc
: Single Parameter (Float64
) - Critical Volume[m³/mol]
v_shift
: Single Parameter (Float64
) - Volume shift[m³/mol]
Description
Peneloux Translation model for cubics:
V = V₀ + mixing_rule(cᵢ)
cᵢ = -0.252*RTcᵢ/Pcᵢ*(1.5448Zcᵢ - 0.4024)
Zcᵢ = Pcᵢ*Vcᵢ/(RTcᵢ)
Model Construction Examples
# Using the default database
translation = PenelouxTranslation("water") #single input
translation = PenelouxTranslation(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
translation = PenelouxTranslation(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/Vc.csv"])
# Passing parameters directly
translation = PenelouxTranslation(["neon","hydrogen"];userlocations = (;Vc = [4.25e-5, 6.43e-5]))
References
- Péneloux A, Rauzy E, Fréze R. (1982) A consistent correction for Redlich‐Kwong‐Soave volumes. Fluid Phase Equilibria 1, 8(1), 7–23. doi:10.1016/0378-3812(82)80002-2
Clapeyron.MTTranslation
— TypeMTTranslation <: MTTranslationModel
MTTranslation(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
acentricfactor
: Single Parameter (Float64
)
Description
Magoulas Tassios Translation model for cubics:
V = V₀ + mixing_rule(cᵢ)
cᵢ = T₀ᵢ+(T̄cᵢ-T̄₀ᵢ)*exp(β*abs(1-Trᵢ))
Trᵢ = T/T̄cᵢ
T̄cᵢ = (RTcᵢ/Pcᵢ)*(0.3074-Zcᵢ)
T̄₀ᵢ = (RTcᵢ/Pcᵢ)*(-0.014471 + 0.067498ωᵢ - 0.084852ωᵢ^2 + 0.067298ωᵢ^3 - 0.017366ωᵢ^4)
Zcᵢ = 0.289 - 0.0701ωᵢ - 0.0207ωᵢ^2
βᵢ = -10.2447 - 28.6312ωᵢ
Model Construction Examples
# Using the default database
translation = MTTranslation("water") #single input
translation = MTTranslation(["water","ethanol"]) #multiple components
# Using user-provided parameters
# Passing files or folders
translation = MTTranslation(["neon","hydrogen"]; userlocations = ["path/to/my/db","critical/acentric.csv"])
# Passing parameters directly
translation = MTTranslation(["neon","hydrogen"];userlocations = (;acentricfactor = [-0.03,-0.21]))
References
- Magoulas, K., & Tassios, D. (1990). Thermophysical properties of n-Alkanes from C1 to C20 and their prediction for higher ones. Fluid Phase Equilibria, 56, 119–140. doi:10.1016/0378-3812(90)85098-u
Mixing Rule Models
Clapeyron.vdW1fRule
— TypevdW1fRule <: vdW1fRuleModel
vdW1fRule(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Description
van der Wals One-Fluid mixing rule for cubic parameters:
aᵢⱼ = √(aᵢaⱼ)(1 - kᵢⱼ)
bᵢⱼ = (1 - lᵢⱼ)(bᵢ + bⱼ)/2
ā = ∑aᵢⱼxᵢxⱼ√(αᵢ(T)αⱼ(T))
b̄ = ∑bᵢⱼxᵢxⱼ
c̄ = ∑cᵢxᵢ
Model Construction Examples
# Because this model does not have parameters, all those constructors are equivalent:
mixing = vdW1fRule()
mixing = vdW1fRule("water")
mixing = vdW1fRule(["water","carbon dioxide"])
Clapeyron.KayRule
— TypeKayRule <: KayRuleModel
KayRule(components;
userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Description
Kay mixing rule for cubic parameters:
aᵢⱼ = √(aᵢaⱼ)(1 - kᵢⱼ)
bᵢⱼ = (1 - lᵢⱼ)(bᵢ + bⱼ)/2
ā = b̄*(∑[aᵢⱼxᵢxⱼ√(αᵢ(T)αⱼ(T))/bᵢⱼ])^2
b̄ = (∑∛(bᵢⱼ)xᵢxⱼ)^3
c̄ = ∑cᵢxᵢ
Model Construction Examples
# Because this model does not have parameters, all those constructors are equivalent:
mixing = KayRule()
mixing = KayRule("water")
mixing = KayRule(["water","carbon dioxide"])
Clapeyron.HVRule
— TypeHVRule{γ} <: HVRuleModel
HVRule(components;
activity = Wilson,
userlocations = String[],
activity_userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Input models
activity
: Activity Model
Description
Huron-Vidal Mixing Rule
b̄ = ∑bᵢⱼxᵢxⱼ
c̄ = ∑cᵢxᵢ
ā = b̄(∑[xᵢaᵢᵢαᵢ/(bᵢᵢ)] - gᴱ/λ)
for Redlich-Kwong:
λ = log(2) (0.6931471805599453)
for Peng-Robinson:
λ = 1/(2√(2))log((2+√(2))/(2-√(2))) (0.6232252401402305)
λ
is a coefficient indicating the relation between gᴱ
and gᴱ(cubic)
at infinite pressure. see [1] for more information. it can be customized by defining HV_λ(::HVRuleModel,::CubicModel,z)
References
- Huron, M.-J., & Vidal, J. (1979). New mixing rules in simple equations of state for representing vapour-liquid equilibria of strongly non-ideal mixtures. Fluid Phase Equilibria, 3(4), 255–271. doi:10.1016/0378-3812(79)80001-1
Model Construction Examples
# Using the default database
mixing = HVRule(["water","carbon dioxide"]) #default: Wilson Activity Coefficient.
mixing = HVRule(["water","carbon dioxide"],activity = NRTL) #passing another Activity Coefficient Model.
mixing = HVRule([("ethane",["CH3" => 2]),("butane",["CH2" => 2,"CH3" => 2])],activity = UNIFAC) #passing a GC Activity Coefficient Model.
# Passing a prebuilt model
act_model = NRTL(["water","ethanol"],userlocations = (a = [0.0 3.458; -0.801 0.0],b = [0.0 -586.1; 246.2 0.0], c = [0.0 0.3; 0.3 0.0]))
mixing = HVRule(["water","ethanol"],activity = act_model)
# Using user-provided parameters
# Passing files or folders
mixing = HVRule(["water","ethanol"]; activity = NRTL, activity_userlocations = ["path/to/my/db","nrtl_ge.csv"])
# Passing parameters directly
mixing = HVRule(["water","ethanol"];
activity = NRTL,
userlocations = (a = [0.0 3.458; -0.801 0.0],
b = [0.0 -586.1; 246.2 0.0],
c = [0.0 0.3; 0.3 0.0])
)
Clapeyron.MHV1Rule
— TypeMHV1Rule{γ} <: MHV1RuleModel
MHV1Rule(components;
activity = Wilson,
userlocations = String[],
activity_userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Input models
activity
: Activity Model
Description
Modified Huron-Vidal Mixing Rule, First Order
bᵢⱼ = (1 - lᵢⱼ)(bᵢ + bⱼ)/2
b̄ = ∑bᵢⱼxᵢxⱼ
c̄ = ∑cᵢxᵢ
ā = b̄RT(∑[xᵢaᵢᵢαᵢ/(RTbᵢᵢ)] - [gᴱ/RT + ∑log(bᵢᵢ/b̄)]/q)
if the model is Peng-Robinson:
q = 0.53
if the model is Redlich-Kwong:
q = 0.593
to use different values for q
, overload Clapeyron.MHV1q(::CubicModel,::MHV1Model) = q
Model Construction Examples
# Using the default database
mixing = MHV1Rule(["water","carbon dioxide"]) #default: Wilson Activity Coefficient.
mixing = MHV1Rule(["water","carbon dioxide"],activity = NRTL) #passing another Activity Coefficient Model.
mixing = MHV1Rule([("ethane",["CH3" => 2]),("butane",["CH2" => 2,"CH3" => 2])],activity = UNIFAC) #passing a GC Activity Coefficient Model.
# Passing a prebuilt model
act_model = NRTL(["water","ethanol"],userlocations = (a = [0.0 3.458; -0.801 0.0],b = [0.0 -586.1; 246.2 0.0], c = [0.0 0.3; 0.3 0.0]))
mixing = MHV1Rule(["water","ethanol"],activity = act_model)
# Using user-provided parameters
# Passing files or folders
mixing = MHV1Rule(["water","ethanol"]; activity = NRTL, activity_userlocations = ["path/to/my/db","nrtl_ge.csv"])
# Passing parameters directly
mixing = MHV1Rule(["water","ethanol"];
activity = NRTL,
userlocations = (a = [0.0 3.458; -0.801 0.0],
b = [0.0 -586.1; 246.2 0.0],
c = [0.0 0.3; 0.3 0.0])
)
References
- Michelsen, M. L. (1990). A modified Huron-Vidal mixing rule for cubic equations of state. Fluid Phase Equilibria, 60(1–2), 213–219. doi:10.1016/0378-3812(90)85053-d
Clapeyron.MHV2Rule
— TypeMHV2Rule{γ} <: MHV2RuleModel
MHV2Rule(components;
activity = Wilson,
userlocations = String[],
activity_userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Input models
activity
: Activity Model
Description
Modified Huron-Vidal Mixing Rule, Second Order.
bᵢⱼ = (1 - lᵢⱼ)(bᵢ + bⱼ)/2
b̄ = ∑bᵢⱼxᵢxⱼ
c̄ = ∑cᵢxᵢ
ᾱᵢ = aᵢαᵢ/bᵢRT
ċ = -q₁*Σᾱᵢxᵢ - q₂*Σᾱᵢxᵢ^2 - gᴱ/RT - ∑log(bᵢᵢ/b̄)
ā = (-q₁ - √(q₁^2 - 4q₂ċ))/(2q₂)
if the model is Peng-Robinson:
q₁ = -0.4347, q₂ = -0.003654
if the model is Redlich-Kwong:
q₁ = -0.4783, q₂ = -0.0047
(-0.4783,-0.0047)
to use different values for q₁
and q₂
, overload Clapeyron.MHV1q(::CubicModel,::MHV2Model) = (q₁,q₂)
Model Construction Examples
# Using the default database
mixing = MHV2Rule(["water","carbon dioxide"]) #default: Wilson Activity Coefficient.
mixing = MHV2Rule(["water","carbon dioxide"],activity = NRTL) #passing another Activity Coefficient Model.
mixing = MHV2Rule([("ethane",["CH3" => 2]),("butane",["CH2" => 2,"CH3" => 2])],activity = UNIFAC) #passing a GC Activity Coefficient Model.
# Passing a prebuilt model
act_model = NRTL(["water","ethanol"],userlocations = (a = [0.0 3.458; -0.801 0.0],b = [0.0 -586.1; 246.2 0.0], c = [0.0 0.3; 0.3 0.0]))
mixing = MHV2Rule(["water","ethanol"],activity = act_model)
# Using user-provided parameters
# Passing files or folders
mixing = MHV2Rule(["water","ethanol"]; activity = NRTL, activity_userlocations = ["path/to/my/db","nrtl_ge.csv"])
# Passing parameters directly
mixing = MHV1Rule(["water","ethanol"];
activity = NRTL,
userlocations = (a = [0.0 3.458; -0.801 0.0],
b = [0.0 -586.1; 246.2 0.0],
c = [0.0 0.3; 0.3 0.0])
)
References
- Michelsen, M. L. (1990). A modified Huron-Vidal mixing rule for cubic equations of state. Fluid Phase Equilibria, 60(1–2), 213–219. doi:10.1016/0378-3812(90)85053-d
Clapeyron.LCVMRule
— TypeLCVMRule{γ} <: LCVMRuleModel
LCVMRule(components;
activity = Wilson,
userlocations = String[],
activity_userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Input models
activity
: Activity Model
Description
Linear Combination of Vidal and Michaelsen (LCVM) Mixing Rule
bᵢⱼ = (1 - lᵢⱼ)(bᵢ + bⱼ)/2
b̄ = ∑bᵢⱼxᵢxⱼ
c̄ = ∑cᵢxᵢ
ᾱᵢ = aᵢαᵢ/bᵢRT
ċ = -q₁*Σᾱᵢxᵢ - q₂*Σᾱᵢxᵢ^2 - gᴱ/RT - ∑log(bᵢᵢ/b̄)
ā = b̄RT(-1.827[gᴱ/RT - 0.3∑log(bᵢᵢ/b̄)] + Σᾱᵢxᵢ)
Model Construction Examples
# Using the default database
mixing = LCVMRule(["water","carbon dioxide"]) #default: Wilson Activity Coefficient.
mixing = LCVMRule(["water","carbon dioxide"],activity = NRTL) #passing another Activity Coefficient Model.
mixing = LCVMRule([("ethane",["CH3" => 2]),("butane",["CH2" => 2,"CH3" => 2])],activity = UNIFAC) #passing a GC Activity Coefficient Model.
# Passing a prebuilt model
act_model = NRTL(["water","ethanol"],userlocations = (a = [0.0 3.458; -0.801 0.0],b = [0.0 -586.1; 246.2 0.0], c = [0.0 0.3; 0.3 0.0]))
mixing = LCVMRule(["water","ethanol"],activity = act_model)
# Using user-provided parameters
# Passing files or folders
mixing = LCVMRule(["water","ethanol"]; activity = NRTL, activity_userlocations = ["path/to/my/db","nrtl_ge.csv"])
# Passing parameters directly
mixing = LCVMRule(["water","ethanol"];
activity = NRTL,
userlocations = (a = [0.0 3.458; -0.801 0.0],
b = [0.0 -586.1; 246.2 0.0],
c = [0.0 0.3; 0.3 0.0])
)
References
- Boukouvalas, C., Spiliotis, N., Coutsikos, P., Tzouvaras, N., & Tassios, D. (1994). Prediction of vapor-liquid equilibrium with the LCVM model: a linear combination of the Vidal and Michelsen mixing rules coupled with the original UNIFAC. Fluid Phase Equilibria, 92, 75–106. doi:10.1016/0378-3812(94)80043-x
Clapeyron.WSRule
— TypeWSRule{γ} <: WSRuleModel
WSRule(components;
activity = Wilson,
userlocations = String[],
activity_userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Input models
activity
: Activity Model
Description
Wong-Sandler Mixing Rule.
aᵢⱼ = √(aᵢaⱼ)(1 - kᵢⱼ)
bᵢⱼ = (bᵢ + bⱼ)/2
c̄ = ∑cᵢxᵢ
B̄ = ΣxᵢxⱼB̄ᵢⱼ
B̄ᵢⱼ = (1 - kᵢⱼ)((bᵢ - aᵢ/RT) + (bⱼ - aⱼ/RT))/2
b̄ = B̄/(1 - gᴱ/λRT - Σxᵢaᵢαᵢ/bᵢRT)
ā = RT(b̄ - B̄)
for Redlich-Kwong:
λ = log(2) (0.6931471805599453)
for Peng-Robinson:
λ = 1/(2√(2))log((2+√(2))/(2-√(2))) (0.6232252401402305)
λ
is a coefficient indicating the relation between gᴱ
and gᴱ(cubic)
at infinite pressure. see [1] for more information. it can be customized by defining WS_λ(::WSRuleModel,::CubicModel)
Model Construction Examples
# Using the default database
mixing = WSRule(["water","carbon dioxide"]) #default: Wilson Activity Coefficient.
mixing = WSRule(["water","carbon dioxide"],activity = NRTL) #passing another Activity Coefficient Model.
mixing = WSRule([("ethane",["CH3" => 2]),("butane",["CH2" => 2,"CH3" => 2])],activity = UNIFAC) #passing a GC Activity Coefficient Model.
# Passing a prebuilt model
act_model = NRTL(["water","ethanol"],userlocations = (a = [0.0 3.458; -0.801 0.0],b = [0.0 -586.1; 246.2 0.0], c = [0.0 0.3; 0.3 0.0]))
mixing = WSRule(["water","ethanol"],activity = act_model)
# Using user-provided parameters
# Passing files or folders
mixing = WSRule(["water","ethanol"]; activity = NRTL, activity_userlocations = ["path/to/my/db","nrtl_ge.csv"])
# Passing parameters directly
mixing = WSRule(["water","ethanol"];
activity = NRTL,
userlocations = (a = [0.0 3.458; -0.801 0.0],
b = [0.0 -586.1; 246.2 0.0],
c = [0.0 0.3; 0.3 0.0])
)
References
- Wong, D. S. H., & Sandler, S. I. (1992). A theoretically correct mixing rule for cubic equations of state. AIChE journal. American Institute of Chemical Engineers, 38(5), 671–680. doi:10.1002/aic.690380505
Clapeyron.modWSRule
— TypemodWSRule{γ} <: WSRuleModel
modWSRule(components;
activity = Wilson,
userlocations = String[],
activity_userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Input models
activity
: Activity Model
Description
modified Wong-Sandler Mixing Rule.
aᵢⱼ = √(aᵢaⱼ)(1 - kᵢⱼ)
bᵢⱼ = (bᵢ + bⱼ)/2
c̄ = ∑cᵢxᵢ
B̄ = Σxᵢxⱼ(bᵢⱼ - aᵢⱼ√(αᵢαⱼ)/RT)
b̄ = B̄/(1 - gᴱ/λRT - Σxᵢaᵢαᵢ/bᵢRT)
ā = RT(b̄ - B̄)
for Redlich-Kwong:
λ = log(2) (0.6931471805599453)
for Peng-Robinson:
λ = 1/(2√(2))log((2+√(2))/(2-√(2))) (0.6232252401402305)
λ
is a coefficient indicating the relation between gᴱ
and gᴱ(cubic)
at infinite pressure. see [1] for more information. it can be customized by defining WS_λ(::WSRuleModel,::CubicModel,z)
Model Construction Examples
# Using the default database
mixing = modWSRule(["water","carbon dioxide"]) #default: Wilson Activity Coefficient.
mixing = modWSRule(["water","carbon dioxide"],activity = NRTL) #passing another Activity Coefficient Model.
mixing = modWSRule([("ethane",["CH3" => 2]),("butane",["CH2" => 2,"CH3" => 2])],activity = UNIFAC) #passing a GC Activity Coefficient Model.
# Passing a prebuilt model
act_model = NRTL(["water","ethanol"],userlocations = (a = [0.0 3.458; -0.801 0.0],b = [0.0 -586.1; 246.2 0.0], c = [0.0 0.3; 0.3 0.0]))
mixing = modWSRule(["water","ethanol"],activity = act_model)
# Using user-provided parameters
# Passing files or folders
mixing = modWSRule(["water","ethanol"]; activity = NRTL, activity_userlocations = ["path/to/my/db","nrtl_ge.csv"])
# Passing parameters directly
mixing = modWSRule(["water","ethanol"];
activity = NRTL,
userlocations = (a = [0.0 3.458; -0.801 0.0],
b = [0.0 -586.1; 246.2 0.0],
c = [0.0 0.3; 0.3 0.0])
)
References
- Wong, D. S. H., & Sandler, S. I. (1992). A theoretically correct mixing rule for cubic equations of state. AIChE journal. American Institute of Chemical Engineers, 38(5), 671–680. doi:10.1002/aic.690380505
- Orbey, H., & Sandler, S. I. (1995). Reformulation of Wong-Sandler mixing rule for cubic equations of state. AIChE journal. American Institute of Chemical Engineers, 41(3), 683–690. doi:10.1002/aic.690410325
Clapeyron.VTPRRule
— TypeVTPRRule{γ} <: VTPRRuleModel
VTPRRule(components;
activity = UNIFAC,
userlocations = String[],
activity_userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Input models
activity
: Activity Model
Description
Mixing Rule used by the Volume-translated Peng-Robinson VTPR
equation of state. only works with activity models that define an excess residual gibbs energy function Clapeyron.excess_g_res(model,P,T,z)
function (like UNIQUAC
and UNIFAC
models)
aᵢⱼ = √(aᵢaⱼ)(1-kᵢⱼ)
bᵢⱼ = ((bᵢ^(3/4) + bⱼ^(3/4))/2)^(4/3)
log(γʳ)ᵢ = lnγ_res(model.activity,V,T,z)
gᴱᵣₑₛ = ∑RTlog(γʳ)ᵢxᵢ
b̄ = ∑bᵢⱼxᵢxⱼ
c̄ = ∑cᵢxᵢ
ā = b̄RT(∑[xᵢaᵢᵢαᵢ/(RTbᵢᵢ)] - gᴱᵣₑₛ/(0.53087RT))
Model Construction Examples
# Note: this model was meant to be used exclusively with the VTPRUNIFAC activity model.
# Using the default database
mixing = VTPRRule(["water","carbon dioxide"]) #default: VTPRUNIFAC Activity Coefficient.
mixing = VTPRRule(["water","carbon dioxide"],activity = NRTL) #passing another Activity Coefficient Model
mixing = VTPRRule([("ethane",["CH3" => 2]),("butane",["CH2" => 2,"CH3" => 2])],activity = UNIFAC) #passing a GC Activity Coefficient Model.
# Passing a prebuilt model
act_model = NRTL(["water","ethanol"],userlocations = (a = [0.0 3.458; -0.801 0.0],b = [0.0 -586.1; 246.2 0.0], c = [0.0 0.3; 0.3 0.0]))
mixing = VTPRRule(["water","ethanol"],activity = act_model)
# Using user-provided parameters
# Passing files or folders
mixing = VTPRRule(["water","ethanol"]; activity = NRTL, activity_userlocations = ["path/to/my/db","nrtl_ge.csv"])
# Passing parameters directly
mixing = VTPRRule(["water","ethanol"];
activity = NRTL,
userlocations = (a = [0.0 3.458; -0.801 0.0],
b = [0.0 -586.1; 246.2 0.0],
c = [0.0 0.3; 0.3 0.0])
)
References
- Ahlers, J., & Gmehling, J. (2001). Development of an universal group contribution equation of state. Fluid Phase Equilibria, 191(1–2), 177–188. doi:10.1016/s0378-3812(01)00626-4
Clapeyron.PSRKRule
— TypePSRKRule{γ} <: MHV1RuleModel
PSRKRule(components;
activity = Wilson,
userlocations = String[],
activity_userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Input models
activity
: Activity Model
Description
Mixing Rule used by the Predictive Soave-Redlich-Kwong PSRK
EoS, derived from the First Order modified Huron-Vidal Mixing Rule.
Model Construction Examples
# Note: this model was meant to be used exclusively with the PSRKUNIFAC activity model.
# Using the default database
mixing = PSRKRule(["water","carbon dioxide"]) #default: PSRKUNIFAC Activity Coefficient.
mixing = PSRKRule(["water","carbon dioxide"],activity = NRTL) #passing another Activity Coefficient Model
mixing = PSRKRule([("ethane",["CH3" => 2]),("butane",["CH2" => 2,"CH3" => 2])],activity = UNIFAC) #passing a GC Activity Coefficient Model.
# Passing a prebuilt model
act_model = NRTL(["water","ethanol"],userlocations = (a = [0.0 3.458; -0.801 0.0],b = [0.0 -586.1; 246.2 0.0], c = [0.0 0.3; 0.3 0.0]))
mixing = PSRKRule(["water","ethanol"],activity = act_model)
# Using user-provided parameters
# Passing files or folders
mixing = PSRKRule(["water","ethanol"]; activity = NRTL, activity_userlocations = ["path/to/my/db","nrtl_ge.csv"])
# Passing parameters directly
mixing = PSRKRule(["water","ethanol"];
activity = NRTL,
userlocations = (a = [0.0 3.458; -0.801 0.0],
b = [0.0 -586.1; 246.2 0.0],
c = [0.0 0.3; 0.3 0.0])
)
Clapeyron.UMRRule
— TypeUMRRule{γ} <: UMRRuleModel
UMRRule(components;
activity = UNIFAC,
userlocations = String[],
activity_userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Input models
activity
: Activity Model
Description
Mixing Rule used by the Universal Mixing Rule Peng-Robinson UMRPR
equation of state.
aᵢⱼ = √(aᵢaⱼ)(1 - kᵢⱼ)
bᵢⱼ = (1 - lᵢⱼ)((√bᵢ +√bⱼ)/2)^2
b̄ = ∑bᵢⱼxᵢxⱼ
c̄ = ∑cᵢxᵢ
ā = b̄RT(∑[xᵢaᵢᵢαᵢ/(RTbᵢᵢ)] - [gᴱ/RT]/0.53)
## Model Construction Examples
Note: this model was meant to be used exclusively with the UNIFAC activity model.
Using the default database
mixing = VTPRRule(["water","carbon dioxide"]) #default: UNIFAC Activity Coefficient. mixing = VTPRRule(["water","carbon dioxide"],activity = NRTL) #passing another Activity Coefficient Model mixing = VTPRRule([("ethane",["CH3" => 2]),("butane",["CH2" => 2,"CH3" => 2])],activity = ogUNIFAC) #passing a GC Activity Coefficient Model.
Passing a prebuilt model
actmodel = NRTL(["water","ethanol"],userlocations = (a = [0.0 3.458; -0.801 0.0],b = [0.0 -586.1; 246.2 0.0], c = [0.0 0.3; 0.3 0.0])) mixing = VTPRRule(["water","ethanol"],activity = actmodel)
Using user-provided parameters
Passing files or folders
mixing = VTPRRule(["water","ethanol"]; activity = NRTL, activityuserlocations = ["path/to/my/db","nrtlge.csv"])
Passing parameters directly
mixing = VTPRRule(["water","ethanol"]; activity = NRTL, userlocations = (a = [0.0 3.458; -0.801 0.0], b = [0.0 -586.1; 246.2 0.0], c = [0.0 0.3; 0.3 0.0]) )
Clapeyron.QCPRRule
— TypeQCPRRule <: MHV2RuleModel
QCPRRule(components;
activity = Wilson,
userlocations = String[],
activity_userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Input models
activity
: Activity Model
Description
Quantum-Corrected Mixing Rule, used by QCPR
EoS:
aᵢⱼ = √(aᵢaⱼ)(1 - kᵢⱼ)
bᵢⱼ = (1 - lᵢⱼ)(bqᵢ + bqⱼ)/2
bqᵢ = bᵢβᵢ(T)
βᵢ(T) = (1 + Aᵢ/(T + Bᵢ))^3 / (1 + Aᵢ/(Tcᵢ + Bᵢ))^3
ā = ∑aᵢⱼxᵢxⱼ√(αᵢ(T)αⱼ(T))
b̄ = ∑bᵢⱼxᵢxⱼ
c̄ = ∑cᵢxᵢ
References
- Aasen, A., Hammer, M., Lasala, S., Jaubert, J.-N., & Wilhelmsen, Ø. (2020). Accurate quantum-corrected cubic equations of state for helium, neon, hydrogen, deuterium and their mixtures. Fluid Phase Equilibria, 524(112790), 112790. doi:10.1016/j.fluid.2020.112790
Clapeyron.PPR78Rule
— TypePPR78Rule <: PPR78RuleModel
PPR78Rule(components;
userlocations = String[],
group_userlocations = String[]
verbose::Bool=false)
Input Parameters
A
: Pair Parameter (Float64
) - Fitted Parameter[K]
B
: Pair Parameter (Float64
) - Fitted Parameter[K]
Description
PPR78 Mixing Rule, Uses E-PPR78 Group params. Default for EPPR78
EoS.
aᵢⱼ = √(aᵢaⱼ)
bᵢⱼ = (bᵢ +bⱼ)/2
b̄ = ∑bᵢⱼxᵢxⱼ
c̄ = ∑cᵢxᵢ
ā = b̄(∑[xᵢaᵢαᵢ/(bᵢᵢ)] - ∑xᵢxⱼbᵢbⱼEᵢⱼ/2b̄)
Eᵢⱼ = ∑(z̄ᵢₖ - z̄ⱼₖ)(z̄ᵢₗ - z̄ⱼₗ) × Aₖₗ × (298.15/T)^(Aₖₗ/Bₖₗ - 1)
References
- Jaubert, J.-N., Privat, R., & Mutelet, F. (2010). Predicting the phase equilibria of synthetic petroleum fluids with the PPR78 approach. AIChE Journal. American Institute of Chemical Engineers, 56(12), 3225–3235. doi:10.1002/aic.12232
- Jaubert, J.-N., Qian, J.-W., Lasala, S., & Privat, R. (2022). The impressive impact of including enthalpy and heat capacity of mixing data when parameterising equations of state. Application to the development of the E-PPR78 (Enhanced-Predictive-Peng-Robinson-78) model. Fluid Phase Equilibria, (113456), 113456. doi:10.1016/j.fluid.2022.113456
Clapeyron.gErRule
— TypegErRule{γ} <: gErRuleModel
gErRule(components;
activity = NRTL,
userlocations = String[],
activity_userlocations = String[],
verbose::Bool=false)
Input Parameters
None
Input models
activity
: Activity Model
Description
Mixing rule that uses the residual part of the activity coefficient model:
bᵢⱼ = ( (bᵢ^(2/3) + bⱼ^(2/3)) / 2 )^(3/2)
b̄ = ∑bᵢⱼxᵢxⱼ
c̄ = ∑cᵢxᵢ
ā/b̄ = b̄RT*(∑xᵢaᵢ/bᵢ + gᴱᵣ/Λ
Λ = 1/(r₂ - r₁) * log((1 - r₂)/(1 - r₁))
Model Construction Examples
# Using the default database
mixing = gErRule(["water","carbon dioxide"]) #default: NRTL Activity Coefficient.
mixing = gErRule(["water","carbon dioxide"],activity = Wilson) #passing another Activity Coefficient Model.
mixing = gErRule([("ethane",["CH3" => 2]),("butane",["CH2" => 2,"CH3" => 2])],activity = UNIFAC) #passing a GC Activity Coefficient Model.
# Passing a prebuilt model
act_model = NRTL(["water","ethanol"],userlocations = (a = [0.0 3.458; -0.801 0.0],b = [0.0 -586.1; 246.2 0.0], c = [0.0 0.3; 0.3 0.0]))
mixing = gErRule(["water","ethanol"],activity = act_model)
# Using user-provided parameters
# Passing files or folders
mixing = gErRule(["water","ethanol"]; activity = NRTL, activity_userlocations = ["path/to/my/db","nrtl_ge.csv"])
# Passing parameters directly
mixing = gErRule(["water","ethanol"];
activity = NRTL,
userlocations = (a = [0.0 3.458; -0.801 0.0],
b = [0.0 -586.1; 246.2 0.0],
c = [0.0 0.3; 0.3 0.0])
)
References
- Piña-Martinez, A., Privat, R., Nikolaidis, I. K., Economou, I. G., & Jaubert, J.-N. (2021). What is the optimal activity coefficient model to be combined with the translated–consistent Peng–Robinson equation of state through advanced mixing rules? Cross-comparison and grading of the Wilson, UNIQUAC, and NRTL aE models against a benchmark database involving 200 binary systems. Industrial & Engineering Chemistry Research, 60(47), 17228–17247. doi:10.1021/acs.iecr.1c03003