Index
Clapeyron.CompositeModelClapeyron.ECSClapeyron.FluidCorrelationClapeyron.GammaPhiClapeyron.LKPClapeyron.LKPSJTClapeyron.PeTSClapeyron.SLClapeyron.SLKRuleClapeyron.SLk0k1lMixingRuleClapeyron.SanchezLacombeClapeyron.ZeroResidualClapeyron.enhancedLKPClapeyron.SPUNGClapeyron.mix_vεClapeyron.shape_factors
Extended Corresponding states Model
Clapeyron.ECS — Type
ExtendedCorrespondingStates <: 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 — Function
shape_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/hFor 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,T,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 — Function
function 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 — Type
LKP <: 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³]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.25Model 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
Clapeyron.LKPSJT — Type
LKPSJT <: LKPModel
LKPSJT(components;
idealmodel=BasicIdeal,
userlocations = String[],
ideal_userlocations = String[],
verbose=false,
reference_state = nothing)
enhancedLKP(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³]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, Sabozin-Jäger-Thol enhancement. Corresponding states using interpolation between a simple, spherical fluid (methane, ∅) and a reference fluid (n-octane, ref):
Instead of using the original BWR formulation, the reference equations of state of methane and octane are used.
αᵣ = (1 - ωᵣ)*αᵣ(δ,τ,params(∅)) + ωᵣ*αᵣ(δ,τ,params(ref))
τ = Tr/T
δ = Vr/V
ωᵣ = (ω̄ - ω(∅))/(ω(ref) - ω(∅))
ω̄ = ∑xᵢωᵢ
Tr = ∑xᵢ*xⱼ*Tcᵢⱼ*Vcᵢⱼ^η * (1-kᵢⱼ)
Vr = ∑xᵢ*xⱼ*Tcᵢⱼ*Vcᵢⱼ
Tcᵢⱼ = √Tcᵢ*Tcⱼ
Vcᵢⱼ = 0.125*(∛Vcᵢ + ∛Vcⱼ)^3
η = 0.25Model Construction Examples
# Using the default database
model = LKPSJT("water") #single input
model = LKPSJT(["water","ethanol"]) #multiple components
model = LKPSJT(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
model = enhancedLKP(["water","ethanol"], idealmodel = ReidIdeal) #enhancedLKP is just an alias to LKPSJT
# Passing a prebuilt model
my_idealmodel = MonomerIdeal(["neon","hydrogen"];userlocations = (;Mw = [20.17, 2.]))
model = LKPSJT(["neon","hydrogen"],idealmodel = my_idealmodel)
# User-provided parameters, passing files or folders
model = LKPSJT(["neon","hydrogen"]; userlocations = ["path/to/my/db","lkp/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = LKPSJT(["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
- Sabozin, F., Jäger, A., & Thol, M. (2024). Enhancement of the Lee–Kesler–Plöcker equation of state for calculating thermodynamic properties of long-chain alkanes. International Journal of Thermophysics, 45(5). doi:10.1007/s10765-024-03360-0
Clapeyron.enhancedLKP — Type
LKPSJT <: LKPModel
LKPSJT(components;
idealmodel=BasicIdeal,
userlocations = String[],
ideal_userlocations = String[],
verbose=false,
reference_state = nothing)
enhancedLKP(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³]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, Sabozin-Jäger-Thol enhancement. Corresponding states using interpolation between a simple, spherical fluid (methane, ∅) and a reference fluid (n-octane, ref):
Instead of using the original BWR formulation, the reference equations of state of methane and octane are used.
αᵣ = (1 - ωᵣ)*αᵣ(δ,τ,params(∅)) + ωᵣ*αᵣ(δ,τ,params(ref))
τ = Tr/T
δ = Vr/V
ωᵣ = (ω̄ - ω(∅))/(ω(ref) - ω(∅))
ω̄ = ∑xᵢωᵢ
Tr = ∑xᵢ*xⱼ*Tcᵢⱼ*Vcᵢⱼ^η * (1-kᵢⱼ)
Vr = ∑xᵢ*xⱼ*Tcᵢⱼ*Vcᵢⱼ
Tcᵢⱼ = √Tcᵢ*Tcⱼ
Vcᵢⱼ = 0.125*(∛Vcᵢ + ∛Vcⱼ)^3
η = 0.25Model Construction Examples
# Using the default database
model = LKPSJT("water") #single input
model = LKPSJT(["water","ethanol"]) #multiple components
model = LKPSJT(["water","ethanol"], idealmodel = ReidIdeal) #modifying ideal model
model = enhancedLKP(["water","ethanol"], idealmodel = ReidIdeal) #enhancedLKP is just an alias to LKPSJT
# Passing a prebuilt model
my_idealmodel = MonomerIdeal(["neon","hydrogen"];userlocations = (;Mw = [20.17, 2.]))
model = LKPSJT(["neon","hydrogen"],idealmodel = my_idealmodel)
# User-provided parameters, passing files or folders
model = LKPSJT(["neon","hydrogen"]; userlocations = ["path/to/my/db","lkp/my_k_values.csv"])
# User-provided parameters, passing parameters directly
model = LKPSJT(["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
- Sabozin, F., Jäger, A., & Thol, M. (2024). Enhancement of the Lee–Kesler–Plöcker equation of state for calculating thermodynamic properties of long-chain alkanes. International Journal of Thermophysics, 45(5). doi:10.1007/s10765-024-03360-0
Sanchez–Lacombe Model
Clapeyron.SanchezLacombe — Type
SanchezLacombe(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³·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³·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.SL — Type
SanchezLacombe(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³·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³·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ε — Function
mix_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 — Type
SLKRule(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 — Type
SLKRule(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 — Type
PeTSModel <: 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[Å]epsilon: Single Parameter (Float64) - Reduced dispersion energy[K]k: Pair Parameter (Float64) (optional) - Binary Interaction Parameter (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
Composite and Utility Models
Clapeyron.CompositeModel — Type
function CompositeModel(components ; mapping = nothing, liquid = nothing, gas = nothing, fluid = nothing, solid = nothing, saturation = nothing, melting = nothing, sublimation = nothing, gasuserlocations = String[], liquiduserlocations = String[], fluiduserlocations = String[], soliduserlocations = String[], saturationuserlocations = String[], meltinguserlocations = String[], sublimationuserlocations = String[], verbose = false, referencestate = nothing)
Model that holds representations of fluid (and/or solid) that aren't evaluated using the Helmholtz energy-based approach used in the rest of the library.
It contains a fluid model, a solid model (optional), and a mapping between the solid and liquid components (if necessary).
There are three available representations for the fluid model:
- A Helmholtz-based EoS.
- Fluid Correlations, consisting in a gas model, a correlation for obtaining the saturation pressure, and a liquid model. Both gas and liquid models can optionally be Helmholtz models too, but correlations for saturated liquid and vapour are also allowed.
- Activity models, consisting of a liquid activity and a model for the fluid. The fluid model can be a Helmholtz-based model, or another
CompositeModelcontaining correlations.
When the solid field is specified, some properties (like volume) start taking in account the solid phase in their calculations. Optionally, there are other models that provide specific correlations for SLE equilibria (like SolidHfus).
The solid model is optional and does not impact VLE (and LLE) calculations. There are two available representations for the solid model:
- a Helmholtz-based EoS, can be used to calculate both melting/sublimation and solubilities.
- Chemical Potential models, can be used for solubilities.
Examples:
- Saturation pressure calculated using Correlations:
#Rackett correlation for liquids, DIPPR 101 correlation for the saturation pressure, ideal gas for the vapour volume.
julia> model = CompositeModel(["water"],liquid = RackettLiquid,saturation = DIPPR101Sat,gas = BasicIdeal)
Composite Model (Correlation-Based) with 1 component:
Gas Model: BasicIdeal()
Liquid Model: RackettLiquid("water")
Saturation Model: DIPPR101Sat("water")
julia> saturation_pressure(model,373.15)
(101260.56298096628, 1.8234039042643886e-5, 0.030639190960720403)- Bubble Pressure, calculated using fluid correlations and a Raoult solver:
julia> model = CompositeModel(["octane","heptane"],liquid = RackettLiquid,saturation = DIPPR101Sat,gas = BasicIdeal)
Composite Model (Correlation-Based) with 2 components:
Gas Model: BasicIdeal()
Liquid Model: RackettLiquid("octane", "heptane")
Saturation Model: DIPPR101Sat("octane", "heptane")
julia> bubble_pressure(model,300.15,[0.9,0.1])
(2552.3661540464022, 0.00015684158726046333, 0.9777538974501402, [0.7376170278676232, 0.2623829721323768])- Bubble Pressure, using an Activity Model along with another model for fluid properties:
#using a Helmholtz-based fluid
julia> model = CompositeModel(["octane","heptane"],liquid = UNIFAC,fluid = PR)
Composite Model (γ-ϕ) with 2 components:
Activity Model: UNIFAC{PR{BasicIdeal, PRAlpha, NoTranslation, vdW1fRule}}("octane", "heptane")
Fluid Model: PR{BasicIdeal, PRAlpha, NoTranslation, vdW1fRule}("octane", "heptane")
julia> bubble_pressure(model,300.15,[0.9,0.1])
(2694.150594740186, 0.00016898441224336215, 0.9239727973658585, [0.7407077952279438, 0.2592922047720562])
#using a correlation-based fluid
julia> fluidmodel = CompositeModel(["octane","heptane"],liquid = RackettLiquid,saturation = DIPPR101Sat,gas = BasicIdeal);
model2 = CompositeModel(["octane","heptane"],liquid = UNIFAC, fluid = fluidmodel)
Composite Model (γ-ϕ) with 2 components:
Activity Model: UNIFAC{PR{BasicIdeal, PRAlpha, NoTranslation, vdW1fRule}}("octane", "heptane")
Fluid Model: FluidCorrelation{BasicIdeal, RackettLiquid, DIPPR101Sat}("octane", "heptane")
julia> bubble_pressure(model2,300.15,[0.9,0.1])
(2551.6008524130893, 0.00015684158726046333, 0.9780471551726359, [0.7378273929683233, 0.2621726070316766])Clapeyron.FluidCorrelation — Type
FluidCorrelation{V,L,Sat,Cp} <: RestrictedEquilibriaModelWrapper struct to signal that a CompositeModel uses correlations for calculation of saturation points, vapour and liquid phase volumes.
Clapeyron.GammaPhi — Type
GammaPhi{γ,Φ} <: RestrictedEquilibriaModelWrapper struct to signal that a CompositeModel uses an activity model in conjunction with a fluid.
Clapeyron.ZeroResidual — Type
ZeroResidual <: EoSModel
ZeroResidual(components;
idealmodel = BasicIdeal,
ideal_userlocations = String[],
verbose = false)Input parameters
None
Description
Zero residual model.
a_res(model,V,T,z) = 0