Index
Clapeyron.AsymmetricMixing
Clapeyron.EmpiricDeparture
Clapeyron.EmpiricIdeal
Clapeyron.GEDeparture
Clapeyron.LJRef
Clapeyron.LinearMixing
Clapeyron.MultiFluid
Clapeyron.QuadraticDeparture
Clapeyron.SingleFluid
Clapeyron.SingleFluidIdeal
Clapeyron.Ammonia2023
Clapeyron.EOS_CG
Clapeyron.EOS_LNG
Clapeyron.GERG2008
Clapeyron.HelmAct
Clapeyron.IAPWS95
Clapeyron.LorentzBerthelotMixing
Clapeyron.PropaneRef
Clapeyron.TholLJ
Clapeyron.XiangDeiters
Clapeyron.create_departure
Clapeyron.departure_functions
Empiric Helmholtz Models
Empiric (or MultiParameter) models in Clapeyron are composed of three different, but interacting parts:
- Pure Fluid parameters
- Mixing volume and temperature
- Departure model
Pure Fluids are instantiated from CoolProp JSON files, via the SingleFluid
struct. In theory, any pure fluid should be supported. Furthermore, there is support for using directly the fluids defined in the CoolProp library:
julia> SingleFluid("Ethanol")
ERROR: cannot found component file R113. Try loading the CoolProp library by loading it.
Stacktrace:
....
julia> using CoolProp # loads the CoolProp library and allows access to their JSON.
julia> SingleFluid("Ethanol")
MultiParameter Equation of state for Ethanol:
Polynomial power terms: 6
Exponential terms: 10
Gaussian bell-shaped terms: 9
MultiComponent models are a collection of SingleFluid
models + a mixing model + a departure model:
julia> model = GERG2008(["water","carbon dioxide"])
MultiFluid{EmpiricAncillary, AsymmetricMixing, EmpiricDeparture} with 2 components:
"water"
"carbon dioxide"
Contains parameters: Mw, Tc, Pc, Vc, acentricfactor, lb_volume
julia> model.pures
2-element Vector{SingleFluid{EmpiricAncillary}}:
SingleFluid{EmpiricAncillary}("water")
SingleFluid{EmpiricAncillary}("carbon dioxide")
julia> model.mixing
AsymmetricMixing with 2 components:
"water"
"carbon dioxide"
Contains parameters: gamma_T, gamma_v, beta_T, beta_v
julia> model.departure
EmpiricDeparture with 2 components:
"water"
"carbon dioxide"
Contains parameters: F, parameters
Generic Models
Clapeyron.SingleFluid
— TypeSingleFluid(components;
userlocations = String[],
ancillaries = nothing,
ancillaries_userlocations = String[],
estimate_pure = false,
coolprop_userlocations = true,
Rgas = nothing,
verbose = false)
Input parameters
- JSON data (CoolProp and teqp format)
Input models
ancillaries
: a model that provides initial guesses for saturation calculations. ifnothing
, then they will be parsed from the input JSON.
Description
Instantiates a single-component Empiric EoS model. Rgas
can be used to set the value of the gas constant that is used during property calculations.
If coolprop_userlocations
is true, then Clapeyron will try to look if the fluid is present in the CoolProp library.
The properties, ideal and residual terms can be accessed via the properties
, ideal
and residual
fields respectively:
julia> model = SingleFluid("water")
MultiParameter Equation of state for water:
Polynomial power terms: 7
Exponential terms: 44
Gaussian bell-shaped terms: 3
Non Analytic terms: 2
julia> model.ideal
Ideal MultiParameter coefficients:
Lead terms: -8.3204464837497 + 6.6832105275932*τ + 3.00632*log(τ)
Plank-Einstein terms: 5
julia> model.residual
Residual MultiParameter coefficients:
Polynomial power terms: 7
Exponential terms: 44
Gaussian bell-shaped terms: 3
Non Analytic terms: 2
Clapeyron.SingleFluidIdeal
— TypeSingleFluidIdeal(components;
userlocations = String[],
Rgas = nothing,
verbose = false,
coolprop_userlocations = true)
Input parameters
- JSON data (CoolProp and teqp format)
Input models
ancillaries
: a model that provides initial guesses for saturation calculations. ifnothing
, then they will be parsed from the input JSON.
Description
Instantiates the ideal part of a single-component Empiric EoS model. Rgas
can be used to set the value of the gas constant that is used during property calculations.
If coolprop_userlocations
is true, then Clapeyron will try to look if the fluid is present in the CoolProp library.
The properties and ideal terms can be accessed via the properties
and ideal
fields respectively:
julia> model = SingleFluidIdeal("water")
Ideal MultiParameter Equation of state for water:
Lead terms: -8.3204464837497 + 6.6832105275932*τ + 3.00632*log(τ)
Plank-Einstein terms: 5
julia> model.ideal
Ideal MultiParameter coefficients:
Lead terms: -8.3204464837497 + 6.6832105275932*τ + 3.00632*log(τ)
Plank-Einstein terms: 5
Clapeyron.MultiFluid
— TypeMultiFluid(components;
idealmodel = nothing,
ideal_userlocations = String[],
pure_userlocations = String[],
mixing = AsymmetricMixing,
departure = EmpiricDeparture,
mixing_userlocations = String[],
departure_userlocations = String[],
estimate_pure = false,
estimate_mixing = :off,
coolprop_userlocations = true,
Rgas = nothing,
reference_state = nothing,
verbose = false)
Input parameters
- JSON data (CoolProp and teqp format)
Input models
idealmodel
: Ideal Model. if it isnothing
, then it will parse the ideal model from the input JSON.mixing
: mixing model for temperature and volume.departure
: departure model
Description
Instantiates a multi-component Empiric EoS model. Rgas
can be used to set the value of the gas constant that is used during property calculations.
If coolprop_userlocations
is true, then Clapeyron will try to look if the fluid is present in the CoolProp library.
If estimate_pure
is true, then, if a JSON is not found, the pure model will be estimated, using the XiangDeiters
model
estimate_mixing
is used to fill missing mixing values in the case of using AsymmetricMixing
. on other mixing models it has no effect.
estimate_mixing = :off
will perform no calculation of mixing parameter, throwing an error if missing values are found.estimate_mixing = :lb
will perform Lorentz-Berthelot estimation of missing mixing parameters. (γT = βT = γv = βv = 1.0). additionally, you can passLorentzBerthelotMixing
to usek
andl
BIP instead.estimate_mixing = :linear
will perform averaging of γT and γv so thatT(x) = ∑xᵢTᵢ
andV(x) = ∑xᵢVᵢ
on missing mixing parameters. Additionally, you can useLinearMixing
to perform this directly.
Rgas
sets the value of the gas constant to be used by the multifluid. The default is the following:
- If
Rgas
is not specified and the input is a single component model, then the value ofRgas
will be taken from the fluid json file. - If
Rgas
is not specified and the input is a multi-component model, then the value ofRgas
will be set toClapeyron.R̄ = Rgas() = 8.31446261815324
(2019 defined constant value)
Clapeyron.EmpiricIdeal
— TypeEmpiricIdeal(components;
pure_userlocations = String[],
estimate_pure = false,
coolprop_userlocations = true,
Rgas = R̄,
verbose = false)
Input parameters
- JSON data (CoolProp and teqp format)
Description
Instantiates the ideal part of a multi-component Empiric EoS model. Rgas
can be used to set the value of the gas constant that is used during property calculations.
If coolprop_userlocations
is true, then Clapeyron will try to look if the fluid is present in the CoolProp library.
If estimate_pure
is true, then, if a JSON is not found, the pure model will be estimated, using the XiangDeiters
model
SingleFluid Models
Clapeyron.XiangDeiters
— FunctionXiangDeiters::SingleFluid
XiangDeiters(component;
idealmodel = BasicIdeal,
userlocations = String[],
ideal_userlocations = String[],
Rgas = 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]
acentricfactor
: Single Parameter (Float64
)
Input models
idealmodel
: Ideal Model
Description
Xiang-Deiters model. estimates a single component Empiric EoS Model from critical values and the acentric factor.
Zc = PcVc/RTc
θ = (Zc - 0.29)^2
aᵣ = a₀(δ,τ) + ω*a₁(δ,τ) + θ*a₂(δ,τ)
Rgas
can be used to set the value of the gas constant that is used during property calculations.
Model Construction Examples
# Using the default database
model = XiangDeiters("water") #single input
model = XiangDeiters(["water"]) #single input, as a vector
model = XiangDeiters(["water"], idealmodel = ReidIdeal) #modifying ideal model
# Passing a prebuilt model
my_idealmodel = MonomerIdeal(["ethane"])
model = XiangDeiters(["ethane"],idealmodel = my_idealmodel)
# User-provided parameters, passing files or folders
model = XiangDeiters(["hydrogen"]; userlocations = ["path/to/my/db","critical.csv"])
# User-provided parameters, passing parameters directly
model = XiangDeiters(["hydrogen"];
userlocations = (;Tc = [44.492],
Pc = [2679000],
Vc = [4.25e-5],
Mw = [2.0],
acentricfactor = [-0.21])
)
references
- Xiang, H. W., & Deiters, U. K. (2008). A new generalized corresponding-states equation of state for the extension of the Lee–Kesler equation to fluids consisting of polar and larger nonpolar molecules. Chemical Engineering Science, 63(6), 1490–1496. doi:10.1016/j.ces.2007.11.029
Clapeyron.IAPWS95
— FunctionIAPWS95 <: EmpiricHelmholtzModel
IAPWS95()
Input parameters
None
Description
IAPWS95 (International Association for the Properties of Water and Steam) Pure water Model, 2018 update.
δ = ρ/ρc
τ = T/Tc
a⁰(δ,τ) = log(δ) + n⁰₁ + n⁰₂τ + n⁰₃log(τ) + ∑n⁰ᵢ(1-exp(-γ⁰ᵢτ)), i ∈ 4:8
aʳ(δ,τ) = aʳ₁+ aʳ₂ + aʳ₃ + aʳ₄
aʳ₁(δ,τ) = ∑nᵢδ^(dᵢ)τ^(tᵢ), i ∈ 1:7
aʳ₂(δ,τ) = ∑nᵢexp(-δ^cᵢ)δ^(dᵢ)τ^(tᵢ), i ∈ 8:51
aʳ₃(δ,τ) = ∑nᵢexp(-αᵢ(δ - εᵢ)^2 - βᵢ(τ - γᵢ)^2)δ^(dᵢ)τ^(tᵢ), i ∈ 52:54
aʳ₄(δ,τ) = ∑nᵢδΨΔ^(bᵢ), i ∈ 55:56
Δ = θ^2 + Bᵢ[(δ - 1)^2]^aᵢ
θ = (1 - τ) + Aᵢ[(δ - 1)^2]^(1/2βᵢ)
Ψ = exp(-Cᵢ(δ - 1)^2 - Dᵢ(τ - 1)^2)
parameters n⁰
,γ⁰
,n
,t
,d
,c
,α
,β
,γ
,ε
,A
,B
,C
,D
where obtained via fitting.
References
- Wagner, W., & Pruß, A. (2002). The IAPWS formulation 1995 for the thermodynamic properties of ordinary water substance for general and scientific use. Journal of physical and chemical reference data, 31(2), 387–535. doi:10.1063/1.1461829
- IAPWS R6-95 (2018). Revised Release on the IAPWS Formulation 1995 for the Thermodynamic Properties of Ordinary Water Substance for General and Scientific Use
Clapeyron.PropaneRef
— FunctionPropaneRef <: EmpiricHelmholtzModel
PropaneRef()
Input parameters
None
Description
Propane Reference Equation of State
δ = ρ/ρc
τ = T/Tc
a⁰(δ,τ) = log(δ) + n⁰₁ + n⁰₂τ + n⁰₃log(τ) + ∑n⁰ᵢ(1-exp(-γ⁰ᵢτ)), i ∈ 4:7
aʳ(δ,τ) = aʳ₁+ aʳ₂ + aʳ₃
aʳ₁(δ,τ) = ∑nᵢδ^(dᵢ)τ^(tᵢ), i ∈ 1:5
aʳ₂(δ,τ) = ∑nᵢexp(-δ^cᵢ)δ^(dᵢ)τ^(tᵢ), i ∈ 6:11
aʳ₃(δ,τ) = ∑nᵢexp(-ηᵢ(δ - εᵢ)^2 - βᵢ(τ - γᵢ)^2)δ^(dᵢ)τ^(tᵢ), i ∈ 12:18
parameters n⁰
,γ⁰
,n
,t
,d
,c
,η
,β
,γ
,ε
where obtained via fitting.
References
- Lemmon, E. W., McLinden, M. O., & Wagner, W. (2009). Thermodynamic properties of propane. III. A reference equation of state for temperatures from the melting line to 650 K and pressures up to 1000 MPa. Journal of Chemical and Engineering Data, 54(12), 3141–3180. doi:10.1021/je900217v
Clapeyron.TholLJ
— FunctionTholLJ()
Lennard-Jones Reference equation of state. valid from 0.5 < T/Tc < 7 and pressures up to p/pc = 500. ``` τᵢ = 1.32/T δᵢ = n/0.31V a⁰ᵢ(δ,τ) = log(δᵢ) + 1.5log(τᵢ) - 1.515151515τᵢ + 6.262265814 a⁰(δ,τ,z) = ∑xᵢ(a⁰ᵢ + log(xᵢ)) aʳ(δ,τ) = aʳ₁+ aʳ₂ + aʳ₃ + aʳ₄ aʳ₁(δ,τ) = ∑nᵢδ^(dᵢ)τ^(tᵢ), i ∈ 1:6 aʳ₂(δ,τ) = ∑nᵢexp(-δ^cᵢ)δ^(dᵢ)τ^(tᵢ), i ∈ 7:12 aʳ₃(δ,τ) = ∑nᵢexp(-ηᵢ(δ - εᵢ)^2 - βᵢ(τ - γᵢ)^2)δ^(dᵢ)τ^(tᵢ), i ∈ 13:23
Clapeyron.Ammonia2023
— FunctionAmmonia2023 <: EmpiricHelmholtzModel
Ammonia2023()
Input parameters
None
Description
Ammonia Reference Equation of State (2023)
δ = ρ/ρc
τ = T/Tc
a⁰(δ,τ) = log(δ) + n⁰₁ + n⁰₂τ + n⁰₃log(τ) + ∑n⁰ᵢ(1-exp(-γ⁰ᵢτ)), i ∈ 4:7
aʳ(δ,τ) = aʳ₁+ aʳ₂ + aʳ₃
aʳ₁(δ,τ) = ∑nᵢδ^(dᵢ)τ^(tᵢ), i ∈ 1:5
aʳ₂(δ,τ) = ∑nᵢexp(-δ^cᵢ)δ^(dᵢ)τ^(tᵢ), i ∈ 6:8
aʳ₃(δ,τ) = ∑nᵢexp(-ηᵢ(δ - εᵢ)^2 - βᵢ(τ - γᵢ)^2)δ^(dᵢ)τ^(tᵢ), i ∈ 9:18
aʳ₃(δ,τ) = ∑nᵢexp(-ηᵢ(δ - εᵢ)^2 - 1/(βᵢ*(τ -γᵢ)^2 + bᵢ))δ^(dᵢ)τ^(tᵢ), i ∈ 19:20
parameters n⁰
,γ⁰
,n
,t
,d
,c
,η
,β
,γ
,ε
where obtained via fitting.
References
- Gao, K., Wu, J., Bell, I. H., Harvey, A. H., & Lemmon, E. W. (2023). A reference equation of state with an associating term for the thermodynamic properties of ammonia. Journal of Physical and Chemical Reference Data, 52(1), 013102. doi:10.1063/5.0128269
MultiComponent Models
Clapeyron.LJRef
— TypeLJRef <: EmpiricHelmholtzModel
LJRef(components;
userlocations = String[],
verbose = false)
Input parameters
sigma
: Single Parameter (Float64
) - particle size [Å]epsilon
: Single Parameter (Float64
) - dispersion energy [K
]Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
k
: Pair Parameter (Float64
) (optional) -sigma
mixing coefficient
Model Parameters
sigma
: Pair Parameter (Float64
) - particle size [m]epsilon
: Pair Parameter (Float64
) - dispersion energy [K
]Mw
: Single Parameter (Float64
) - Molecular Weight[g/mol]
Description
Lennard-Jones Reference equation of state. valid from 0.5 < T/Tc < 7 and pressures up to p/pc = 500.
σᵢⱼ = (σᵢ + σⱼ)/2
ϵᵢⱼ = (1-kᵢⱼ)√(ϵⱼϵⱼ)
σ^3 = Σxᵢxⱼσᵢⱼ^3
ϵ = Σxᵢxⱼϵᵢⱼσᵢⱼ^3/σ^3
τᵢ = 1.32ϵᵢ/T
δᵢ = n(Nₐσᵢ^3)/0.31V
a⁰ᵢ(δ,τ) = log(δᵢ) + 1.5log(τᵢ) - 1.515151515τᵢ + 6.262265814
a⁰(δ,τ,z) = ∑xᵢ(a⁰ᵢ + log(xᵢ))
τ = 1.32ϵ/T
δ = n(Nₐσ^3)/0.31V
aʳ(δ,τ) = aʳ₁+ aʳ₂ + aʳ₃ + aʳ₄
aʳ₁(δ,τ) = ∑nᵢδ^(dᵢ)τ^(tᵢ), i ∈ 1:6
aʳ₂(δ,τ) = ∑nᵢexp(-δ^cᵢ)δ^(dᵢ)τ^(tᵢ), i ∈ 7:12
aʳ₃(δ,τ) = ∑nᵢexp(-ηᵢ(δ - εᵢ)^2 - βᵢ(τ - γᵢ)^2)δ^(dᵢ)τ^(tᵢ), i ∈ 13:23
parameters n
,t
,d
,c
,η
,β
,γ
,ε
where obtained via fitting.
The original model was done with only one component in mind. to support multiple components, a VDW 1-fluid mixing rule (shown above) is implemented, but it is not tested.
References
- Thol, M., Rutkai, G., Köster, A., Lustig, R., Span, R., & Vrabec, J. (2016). Equation of state for the Lennard-Jones fluid. Journal of physical and chemical reference data, 45(2), 023101. doi:10.1063/1.4945000
Clapeyron.GERG2008
— FunctionGERG2008::MultiFluid
GERG2008(components;
Rgas = 8.314472,
reference_state = nothing,
verbose = false)
input Parameters
None
Description
The GERG-2008 Wide-Range Equation of State for Natural Gases and Other Mixtures. valid for 21 compounds (Clapeyron.GERG2008_names
).
a = a⁰ + aʳ
a⁰ = ∑xᵢ(a⁰ᵢ(τᵢ,δᵢ) + ln(xᵢ))
δᵢ = ρ/ρcᵢ
τᵢ = Tcᵢ/T
a⁰ᵢ = ln(δᵢ) + R∗/R[n⁰ᵢ₋₁ + n⁰ᵢ₋₂τᵢ + n⁰ᵢ₋₃ln(τᵢ) + ∑n⁰ᵢ₋ₖln(abs(sinh(ϑ₀ᵢ₋ₖτᵢ))) + ∑n⁰ᵢ₋ₖln(cosh(ϑ₀ᵢ₋ₖτᵢ))]
R∗ = 8.314510
R = 8.314472
τ = Tᵣ/T
δ = ρ/ρᵣ
(1/ρᵣ) = ∑∑xᵢxⱼβᵥ₋ᵢⱼγᵥ₋ᵢⱼ[(xᵢ+xⱼ)/(xᵢβᵥ₋ᵢⱼ^2 + xⱼ)]•1/8(1/∛ρcᵢ + 1/∛ρcⱼ)^2
Tᵣ = ∑∑xᵢxⱼβₜ₋ᵢⱼγₜ₋ᵢⱼ[(xᵢ+xⱼ)/(xᵢβₜ₋ᵢⱼ^2 + xⱼ)]•√(TcᵢTcⱼ)
aʳ = ∑xᵢaᵣᵢ(τ,δ) + ∑∑xᵢxⱼFᵢⱼaʳᵢⱼ(τ,δ)
aʳᵢ = ∑nᵢ₋ₖδ^(dᵢ₋ₖ)τ^(tᵢ₋ₖ) + ∑nᵢ₋ₖδ^(dᵢ₋ₖ)τ^(tᵢ₋ₖ)exp(-δ^cᵢ₋ₖ)
aʳᵢⱼ = ∑nᵢⱼ₋ₖδ^(dᵢⱼ₋ₖ)τ^(tᵢⱼ₋ₖ) + ∑nᵢⱼ₋ₖδ^(dᵢⱼ₋ₖ)τ^(tᵢⱼ₋ₖ)exp(ηᵢⱼ₋ₖ(δ-εᵢⱼ₋ₖ)^2 + βᵢⱼ₋ₖ(δ-γᵢⱼ₋ₖ))
References
- Kunz, O., & Wagner, W. (2012). The GERG-2008 wide-range equation of state for natural gases and other mixtures: An expansion of GERG-2004. Journal of Chemical and Engineering Data, 57(11), 3032–3091. doi:10.1021/je300655b
Clapeyron.EOS_LNG
— FunctionEOS_LNG::MultiFluid
EOS_LNG(components::Vector{String};
Rgas = 8.314472,
reference_state = nothing,
verbose = false)
input Parameters
None
Description
EOS-LNG: A Fundamental Equation of State for the Calculation of Thermodynamic Properties of Liquefied Natural Gases. valid for 21 compounds (Clapeyron.GERG2008_names
). the EoS has new binary-specific parameters for methane + n-butane, methane + isobutane, methane + n-pentane, and methane + isopentane.
It uses the same functional form as GERG2008
.
References
- Thol, M., Richter, M., May, E. F., Lemmon, E. W., & Span, R. (2019). EOS-LNG: A fundamental equation of state for the calculation of thermodynamic properties of liquefied natural gases. Journal of Physical and Chemical Reference Data, 48(3), 033102. doi:10.1063/1.5093800
- Kunz, O., & Wagner, W. (2012). The GERG-2008 wide-range equation of state for natural gases and other mixtures: An expansion of GERG-2004. Journal of Chemical and Engineering Data, 57(11), 3032–3091. doi:10.1021/je300655b
Clapeyron.EOS_CG
— FunctionEOS_LNG::MultiFluid
EOS_LNG(components::Vector{String};
Rgas = R̄,
reference_state = nothing,
verbose = false)
input Parameters
None
Description
EOS-LNG: A Fundamental Equation of State for the Calculation of Thermodynamic Properties of Liquefied Natural Gases. valid for 21 compounds (Clapeyron.GERG2008_names
). the EoS has new binary-specific parameters for methane + n-butane, methane + isobutane, methane + n-pentane, and methane + isopentane.
It uses the same functional form as GERG2008
.
References
EOS-CG: : A Mixture Model for the Calculation of Thermodynamic Properties of CCS Mixtures
It uses the same functional form as GERG2008
.
References
- Neumann, T., Herrig, S., Bell, I. H., Beckmüller, R., Lemmon, E. W., Thol, M., & Span, R. (2023). EOS-CG-2021: A mixture model for the calculation of thermodynamic properties of CCS mixtures. International Journal of Thermophysics, 44(12). doi:10.1007/s10765-023-03263-6
Clapeyron.HelmAct
— FunctionHelmAct::MultiFluid
HelmAct(components;
pure_userlocations = String[],
activity = PSRKUNIFAC,
activity_userlocations = String[],
estimate_pure = false,
coolprop_userlocations = false,
Rgas = R̄,
reference_state = nothing,
verbose = verbose)
input Parameters
- CoolProp JSON pure fluid files
Input Models
activity
: activity model.
Description
Creates a Helmholtz + Activity (HelmAct) model:
aᵣ = ∑xᵢaᵣᵢ(δ,τ) + Δa
Δa = gᴱᵣ/RT - log(1+bρ)/log(1+bρref) * ∑xᵢ(aᵣᵢ(δref,τ) - aᵣᵢ(δrefᵢ,τᵢ))
τᵢ = Tcᵢ/T
δref = ρref/ρr
δrefᵢ = ρrefᵢ/ρcᵢ
b = 1/1.17ρref
where gᴱᵣ
is the residual part of the excess gibbs free energy obtained from an activity model.
References
- Jäger, A., Breitkopf, C., & Richter, M. (2021). The representation of cross second virial coefficients by multifluid mixture models and other equations of state. Industrial & Engineering Chemistry Research, 60(25), 9286–9295. doi:10.1021/acs.iecr.1c01186
Mixing Models
Clapeyron.LinearMixing
— TypeLinearMixing <: MultiFluidDepartureModel
LinearMixing(components;
userlocations = String[],
verbose = false)
Input parameters
none
Description
Linear mixing rule for MultiParameter EoS models:
τ = T̄/T
δ = V̄/V
V̄ = ∑xᵢVcⱼ
T̄ = ∑xᵢTcᵢ
Model Construction Examples
# Because this model does not have parameters, all those constructors are equivalent:
mixing = LinearMixing()
mixing = LinearMixing("water")
mixing = LinearMixing(["water","carbon dioxide"])
Clapeyron.AsymmetricMixing
— TypeAsymmetricMixing <: MultiFluidDepartureModel
AsymmetricMixing(components;
userlocations = String[],
verbose = false)
Input parameters
beta_v
: Pair Parameter (Float64
) - binary interaction parameter (no units)gamma_v
: Pair Parameter (Float64
) - binary interaction parameter (no units)beta_T
: Pair Parameter (Float64
) - binary interaction parameter (no units)gamma_T
: Pair Parameter (Float64
) - binary interaction parameter (no units)
Description
Asymmetric mixing rule for MultiParameter EoS models:
τ = T̄/T
δ = V̄/V
V̄ = ∑xᵢxⱼ * βᵛᵢⱼ * γᵛᵢⱼ * (xᵢ + xⱼ)/(xᵢ*βᵛᵢⱼ^2 + xⱼ) * Vᵣᵢⱼ
T̄ = ∑xᵢxⱼ * βᵛᵢⱼ * γᵀᵢⱼ * (xᵢ + xⱼ)/(xᵢ*βᵀᵢⱼ^2 + xⱼ) * Tᵣᵢⱼ
Vᵣᵢⱼ = 0.125*(∛Vcᵢ + ∛Vcⱼ)^3
Tᵣᵢⱼ = √(Tcᵢ*Tcⱼ)
With the asymmetry present in the β parameters:
βᵛᵢⱼ = 1/βᵛⱼᵢ
βᵀᵢⱼ = 1/βᵀⱼᵢ
If there is no data present, the parameters can be estimated:
- Linear estimation:
βᵛᵢⱼ = βᵛᵢⱼ = 1
γᵛᵢⱼ = 4*(Vcᵢ + Vcⱼ)/(∛Vcᵢ + ∛Vcⱼ)^3
γᵀᵢⱼ = 0.5*(Tcᵢ + Tcⱼ)/√(Tcᵢ*Tcⱼ)
- Lorentz-Berthelot Estimation:
βᵛᵢⱼ = βᵛᵢⱼ = γᵛᵢⱼ = γᵀᵢⱼ = 1
References
- R. Klimeck, Ph.D. dissertation, Ruhr-Universit¨at Bochum, 2000
Clapeyron.LorentzBerthelotMixing
— FunctionLorentzBerthelotMixing::AsymmetricMixing
LorentzBerthelotMixing(components;
userlocations = String[],
verbose = false)
Input parameters
k
: Pair Parameter (Float64
) - binary interaction parameter for temperature (no units)l
: Pair Parameter (Float64
) - binary interaction parameter for volume (no units)
Description
Lorentz-Berthelot Mixing for MultiParameter EoS models:
τ = T̄/T
δ = V̄/V
V̄ = ∑xᵢxⱼ * Vᵣᵢⱼ * (1 - lᵢⱼ)
T̄ = ∑xᵢxⱼ * Tᵣᵢⱼ * (1 - kᵢⱼ)
Vᵣᵢⱼ = 0.125*(∛Vcᵢ + ∛Vcⱼ)^3
Tᵣᵢⱼ = √(Tcᵢ*Tcⱼ)
missing parameters will be assumed kᵢⱼ = lᵢⱼ = 0
Departure Models
Clapeyron.EmpiricDeparture
— TypeEmpiricDeparture <: MultiFluidDepartureModel EmpiricDeparture(components; userlocations = String[], verbose = false)
Input parameters
none
F
: Pair Parameter (Float64
) - binary interaction parameter (no units)parameters
: Pair Parameter (String
) - JSON data containing the departure terms for the binary pair
Description
Departure that uses empiric departure functions:
aᵣ = ∑xᵢaᵣᵢ(δ,τ) + Δa
Δa = ∑xᵢxⱼFᵢⱼaᵣᵢⱼ(δ,τ)
aᵣᵢⱼ = ∑nᵢⱼ₋ₖδ^(dᵢⱼ₋ₖ)*τ^(tᵢⱼ₋ₖ) +
∑nᵢⱼ₋ₖδ^(dᵢⱼ₋ₖ)τ^(tᵢⱼ₋ₖ)*exp(-gᵢⱼ₋ₖδ^lᵢⱼ₋ₖ) +
∑nᵢⱼ₋ₖδ^(dᵢⱼ₋ₖ)τ^(tᵢⱼ₋ₖ)*exp(ηᵢⱼ₋ₖ(δ-εᵢⱼ₋ₖ)^2 + βᵢⱼ₋ₖ(τ-γᵢⱼ₋ₖ)^2)
Clapeyron.departure_functions
— Functiondeparture_functions(model::MultiFluid)
if the model is using a EmpiricDeparture
departure model, return the matrix of departure functions. you can set a departure in the following way:
using CoolProp #load CoolProp models
model = MultiFluid(["helium","methanol"],mixing = LorentzBerthelotMixing)
dep_mat = departure_functions(model)
dep = Dict(
#reduced CoolProp Departure format, you only need the type and parameters.
#ResidualHelmholtzPower would work too.
type => "Exponential",
n => [1,1,1,1],
t => [1,1,1,1],
d => [1,1,1,1],
l => [1,1,1,1],
)
dep_mat[1,2] = create_departure(dep,F)
#if you want to delete a departure model:
dep_mat[1,2] = nothing
using SparseArrays
Clapeyron.create_departure
— Functioncreate_departure(data,F = 1.0;verbose = false)
Creates a departure model for use in a MultiFluid
model with EmpiricDeparture
.
If data
is a String
and starts with {
or [
, it will be recognized as JSON text. the text will be parsed as a file location otherwise. You can pass a Dict
or NamedTuple
if you want to skip the JSON parsing.
Examples
d1 = create_departure("/data/EthanePropane.json",0.9) # reading from a file
dep = Dict(
#reduced CoolProp Departure format, you only need the type and parameters.
#ResidualHelmholtzPower would work too.
:type => "Exponential",
:n => [1,1,1,1],
:t => [1,1,1,1],
:d => [1,1,1,1],
:l => [1,1,1,1],
)
d2 = create_departure(dep) #F is set to 1.0
Clapeyron.GEDeparture
— TypeGEDeparture <: MultiFluidDepartureModel GEDeparture(components; activity = UNIFAC, userlocations = String[], verbose = false)
Input parameters
none
k1
: Pair Parameter (Float64
) - binary, T-dependent interaction parameter[K^-1]
Model parameters
vref
: Single Parameter (Float64
, calculated) - Reference pure molar volume[m3/mol]
Input models
activity
: activity model
Description
Departure that uses the residual excess gibbs energy from an activity model:
aᵣ = ∑xᵢaᵣᵢ(δ,τ) + Δa
Δa = gᴱᵣ/RT - log(1+bρ)/log(1+bρref) * ∑xᵢ(aᵣᵢ(δref,τ) - aᵣᵢ(δrefᵢ,τᵢ))
τᵢ = Tcᵢ/T
δref = ρref/ρr
δrefᵢ = ρrefᵢ/ρcᵢ
b = 1/1.17ρref
References
- Jäger, A., Breitkopf, C., & Richter, M. (2021). The representation of cross second virial coefficients by multifluid mixture models and other equations of state. Industrial & Engineering Chemistry Research, 60(25), 9286–9295. doi:10.1021/acs.iecr.1c01186
Clapeyron.QuadraticDeparture
— TypeQuadraticDeparture <: MultiFluidDepartureModel
QuadraticDeparture(components;
userlocations = String[],
verbose = false)
Input parameters
k0
: Pair Parameter (Float64
) - binary interaction parameter (no units)k1
: Pair Parameter (Float64
) - binary, T-dependent interaction parameter[K^-1]
Description
Departure that uses a quadratic mixing rule:
aᵣ = ∑xᵢxⱼaᵣᵢⱼ
aᵣᵢⱼ = 0.5*(aᵣᵢ + aᵣⱼ)*(1 - (k₀ + k₁T))
References
- Jäger, A., Breitkopf, C., & Richter, M. (2021). The representation of cross second virial coefficients by multifluid mixture models and other equations of state. Industrial & Engineering Chemistry Research, 60(25), 9286–9295. doi:10.1021/acs.iecr.1c01186