Index
Clapeyron.ECS
Clapeyron.LKP
Clapeyron.PeTS
Clapeyron.SLKRule
Clapeyron.SLk0k1lMixingRule
Clapeyron.SanchezLacombe
Clapeyron.SPUNG
Clapeyron.mix_vε
Clapeyron.shape_factors
Extended Corresponding states Model
Clapeyron.ECS
— TypeExtendedCorrespondingStates <: EoSModel
function ECS(components,
refmodel=PropaneRef(),
shapemodel=SRK(components),
shaperef = SRK(refmodel.components))
Input Models
shape_model
: shape modelshape_ref
: shape reference. is the same type of EoS thatshape_model
model_ref
: Reference model
Description
A Extended Corresponding states method.
The idea is to use a "shape model" that provides a corresponding states parameters and a "reference model" that implements a helmholtz energy function, so that:
eos(shape_model,v,T,x)/RT = eos(model_ref,v₀,T₀)/RT₀
where:
T₀ = T/f
v₀ = v/h
f,h = shape_factors(model::ECS,shape_ref::EoSModel,V,T,z)
shape_factors
can be used to create custom Extended Corresponding state models.
References
.1 Mollerup, J. (1998). Unification of the two-parameter equation of state and the principle of corresponding states. Fluid Phase Equilibria, 148(1–2), 1–19. doi:10.1016/s0378-3812(98)00230-1
Clapeyron.shape_factors
— Functionshape_factors(model::ECS,V,T,z=SA[1.0])
shape_factors(model::ECS,shape_ref::ABCubicModel,V,T,z=SA[1.0])
shape_factors(model::ECS,shape_ref::EoSModel,V,T,z=SA[1.0])
Returns f
and h
scaling factors, used by the ECS
Equation of state.
eos(shape_model,v,T,x)/RT = eos(model_ref,v₀,T₀)/RT₀
where:
T₀ = T/f
v₀ = v/h
For cubics, a general procedure is defined in [1]:
h = b/b₀
fh = a(T)/a₀(T₀)
For general EoS, there is no existent publications on how to obtain shape factors. However, we can "map" any EoS to a cubic with:
b ≈ lb_volume(model,z)
a ≈ RT*(b - B)
B = second_virial_coefficient(model,T)
This is not tested extensively and it is considered an Experimental feature, subject to future changes.
References
- Mollerup, J. (1998). Unification of the two-parameter equation of state and the principle of corresponding states. Fluid Phase Equilibria, 148(1–2), 1–19. doi:10.1016/s0378-3812(98)00230-1
Clapeyron.SPUNG
— Functionfunction function SPUNG(components,
refmodel=PropaneRef(),
shapemodel=SRK(components),
shaperef = SRK(refmodel.components))
Description
SPUNG: State Research Program for Utilization of Natural Gas
ECS
method. It uses SRK
as the shape model and PropaneRef
as the reference model.
References
- Wilhelmsen, Ø., Skaugen, G., Jørstad, O., & Li, H. (2012). Evaluation of SPUNG* and other equations of state for use in carbon capture and storage modelling. Energy Procedia, 23, 236–245. doi:10.1016/j.egypro.2012.06.024
Clapeyron.LKP
— TypeLKP <: EmpiricHelmholtzModel
LKP(components;
idealmodel=BasicIdeal,
verbose=false)
Input parameters
Tc
: Single Parameter (Float64
) - Critical Temperature[K]
Pc
: Single Parameter (Float64
) - Critical Pressure[Pa]
Vc
: Single Parameter (Float64
) (optional) - Critical Volume[m^3]
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
acentricfactor
: Single Parameter (Float64
) - Acentric Factor (no units)k
: Pair Parameter (Float64
) (optional) - binary interaction parameter (no units)
Input models
idealmodel
: Ideal Model
Description
Lee-Kesler-Plöker equation of state. corresponding states using interpolation between a simple, spherical fluid (methane, ∅
) and a reference fluid (n-octane, ref
):
αᵣ = (1 - ωᵣ)*αᵣ(δr,τ,params(∅)) + ωᵣ*αᵣ(δr,τ,params(ref))
τ = Tr/T
δr = Vr/V/Zr
Zr = Pr*Vr/(R*Tr)
Pr = (0.2905 - 0.085*ω̄)*R*Tr/Vr
ωᵣ = (ω̄ - ω(∅))/(ω(ref) - ω(∅))
ω̄ = ∑xᵢωᵢ
Tr = ∑xᵢ*xⱼ*Tcᵢⱼ*Vcᵢⱼ^η * (1-kᵢⱼ)
Vr = ∑xᵢ*xⱼ*Tcᵢⱼ*Vcᵢⱼ
Tcᵢⱼ = √Tcᵢ*Tcⱼ
Vcᵢⱼ = 0.125*(∛Vcᵢ + ∛Vcⱼ)^3
η = 0.25
Model Construction Examples
# Using the default database
model = LKP("water") #single input
model = LKP(["water","ethanol"]) #multiple components
model = LKP(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
# Passing a prebuilt model
my_idealmodel = MonomerIdeal(["neon","hydrogen"];userlocations = (;Mw = [20.17, 2.]))
model = LKP(["neon","hydrogen"],idealmodel = my_idealmodel)
# User-provided parameters, passing files or folders
model = LKP(["neon","hydrogen"]; userlocations = ["path/to/my/db","lkp/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = LKP(["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.
)
References
- Plöcker, U., Knapp, H., & Prausnitz, J. (1978). Calculation of high-pressure vapor-liquid equilibria from a corresponding-states correlation with emphasis on asymmetric mixtures. Industrial & Engineering Chemistry Process Design and Development, 17(3), 324–332. doi:10.1021/i260067a020
Sanchez–Lacombe Model
Clapeyron.SanchezLacombe
— TypeSanchezLacombe(components;
idealmodel = BasicIdeal,
mixing = SLk0k1lMixingRule,
userlocations = String[],
ideal_userlocations = String[],
mixing_userlocations = String[],
reference_state = false,
verbose = false)
Input parameters
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
segment
: Single Parameter (Float64
) - Number of segments (no units)epsilon
: Single Parameter (Float64
) - Nonbonded interaction energy per monomer[J/mol]
vol
: Single Parameter (Float64
) - Closed Packed Specific volume[m^3/mol]
Model Parameters
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
segment
: Single Parameter (Float64
) - Number of segments (no units)epsilon
: Pair Parameter (Float64
) - Nonbonded interaction energy per monomer[J/mol]
vol
: Pair Parameter (Float64
) - Closed Packed Specific volume[m^3/mol]
Input models
idealmodel
: Ideal Modelmixing
: Mixing model
Description
Sanchez-Lacombe Lattice Fluid Equation of State.
xᵢ = zᵢ/∑zᵢ
r̄ = ∑xᵢrᵢ
vᵣ,εᵣ = mix_vε(model,V,T,z,model.mixing,r̄,∑zᵢ)
ρ̃ = r̄*vᵣ/v
T̃ = R̄*T/εᵣ
aᵣ = r̄*(- ρ̃ /T̃ + (1/ρ̃ - 1)*log(1 - ρ̃ ) + 1)
References
- Neau, E. (2002). A consistent method for phase equilibrium calculation using the Sanchez–Lacombe lattice–fluid equation-of-state. Fluid Phase Equilibria, 203(1–2), 133–140. doi:10.1016/s0378-3812(02)00176-0
Clapeyron.mix_vε
— Functionmix_vε(model::SanchezLacombeModel,V,T,z,mix::SLMixingRule,r̄ = @f(rmix),∑z = sum(z))
Function used to dispatch on the different mixing rules available for Sanchez-Lacombe.
Example:
function mix_vε(model::SanchezLacombe,V,T,z,mix::SLKRule,r̄,Σz = sum(z))
v = model.params.vol.values
ε = model.params.epsilon.values
r = model.params.segment.values
k = mix.k.values
x = z ./ Σz
ϕ = @. r * x / r̄
εᵣ = sum(ε[i,j]*(1-k[i,j])*ϕ[i]*ϕ[j] for i ∈ @comps for j ∈ @comps)
vᵣ = sum(v[i,j]*ϕ[i]*ϕ[j] for i ∈ @comps for j ∈ @comps)
return vᵣ,εᵣ
Clapeyron.SLKRule
— TypeSLKRule(components; userlocations = String[], verbose = false)
Input parameters
k
: Pair Parameter (Float64
) (optional) - Binary Interaction Parameter (no units)
Constant Kᵢⱼ mixing rule for Sanchez-Lacombe:
εᵢⱼ = √εᵢεⱼ*(1-kᵢⱼ)
vᵢⱼ = (vᵢ + vⱼ)/2
ϕᵢ = rᵢ*xᵢ/r̄
εᵣ = ΣΣϕᵢϕⱼεᵢⱼ
vᵣ = ΣΣϕᵢϕⱼvᵢⱼ
Clapeyron.SLk0k1lMixingRule
— TypeSLKRule(components; userlocations = String[], verbose = false)
Input parameters
k0
,k
: Pair Parameter (Float64
, optional) - Binary Interaction Parameter (no units)k1
: Pair Parameter (Float64
, optional) - Binary Interaction Parameter (no units)l
: Pair Parameter (Float64
,optional) - Binary Interaction Parameter (no units)
Neau's Consistent k₀,k₁,l mixing rule for Sanchez-Lacombe:
εᵢⱼ = √εᵢεⱼ
vᵢⱼ = (1 - lᵢⱼ)(vᵢ + vⱼ)/2
ϕᵢ = rᵢ*xᵢ/r̄
εᵣ = ΣΣϕᵢϕⱼεᵢⱼ*(1 - k₀ᵢⱼ + (1 - δᵢⱼ)(Σϕₖk₁ᵢₖ + Σϕₖk₁ₖⱼ))
vᵣ = ΣΣϕᵢϕⱼvᵢⱼ
Where δᵢⱼ
is i == j ? 1 : 0
References
- Neau, E. (2002). A consistent method for phase equilibrium calculation using the Sanchez–Lacombe lattice–fluid equation-of-state. Fluid Phase Equilibria, 203(1–2), 133–140. doi:10.1016/s0378-3812(02)00176-0
Other molecular Models
Clapeyron.PeTS
— TypePeTSModel <: EoSModel
PeTS(components;
idealmodel = BasicIdeal,
userlocations = String[],
ideal_userlocations = String[],
reference_state = nothing,
verbose = false)
Input parameters
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
segment
: Single Parameter (Float64
) - Number of segments (no units)sigma
: Single Parameter (Float64
) - Segment Diameter [A°
]epsilon
: Single Parameter (Float64
) - Reduced dispersion energy[K]
k
: Pair Parameter (Float64
) (optional) - Binary Interaction Paramater (no units)
Model Parameters
Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
segment
: Single Parameter (Float64
) - Number of segments (no units)sigma
: Pair Parameter (Float64
) - Mixed segment Diameter[m]
epsilon
: Pair Parameter (Float64
) - Mixed reduced dispersion energy[K]
Input models
idealmodel
: Ideal Model
Description
Perturbed, Truncated and Shifted (PeTS) Equation of State.
References
- Heier, M., Stephan, S., Liu, J., Chapman, W. G., Hasse, H., & Langenbach, K. (2018). Equation of state for the Lennard-Jones truncated and shifted fluid with a cut-off radius of 2.5 σ based on perturbation theory and its applications to interfacial thermodynamics. Molecular Physics, 116(15–16), 2083–2094. doi:10.1080/00268976.2018.1447153