Contents
Index
Clapeyron.ActivityBubblePressure
Clapeyron.ActivityBubbleTemperature
Clapeyron.ActivityDewPressure
Clapeyron.ActivityDewTemperature
Clapeyron.ChemPotBubblePressure
Clapeyron.ChemPotBubbleTemperature
Clapeyron.ChemPotDewPressure
Clapeyron.ChemPotDewTemperature
Clapeyron.DETPFlash
Clapeyron.FugBubblePressure
Clapeyron.FugBubbleTemperature
Clapeyron.FugDewPressure
Clapeyron.FugDewTemperature
Clapeyron.MCFlashJL
Clapeyron.MichelsenTPFlash
Clapeyron.MultiPhaseTPFlash
Clapeyron.RRTPFlash
Clapeyron.LLE_pressure
Clapeyron.LLE_temperature
Clapeyron.UCEP_mix
Clapeyron.UCST_mix
Clapeyron.VLLE_pressure
Clapeyron.VLLE_temperature
Clapeyron.VT_chemical_stability
Clapeyron.VT_diffusive_stability
Clapeyron.VT_mechanical_stability
Clapeyron.azeotrope_pressure
Clapeyron.azeotrope_temperature
Clapeyron.bubble_pressure
Clapeyron.bubble_temperature
Clapeyron.crit_mix
Clapeyron.cross_second_virial
Clapeyron.dew_pressure
Clapeyron.dew_temperature
Clapeyron.equivol_cross_second_virial
Clapeyron.eutectic_point
Clapeyron.gibbs_duhem
Clapeyron.gibbs_solvation
Clapeyron.isstable
Clapeyron.numphases
Clapeyron.sle_solubility
Clapeyron.slle_solubility
Clapeyron.supports_reduction
Clapeyron.tp_flash
Clapeyron.tpd
Multi component properties
Clapeyron.bubble_pressure
— Functionbubble_pressure(model::EoSModel, T, x, method = ChemPotBubblePressure())
Calculates the bubble pressure and properties at a given temperature. Returns a tuple, containing:
- Bubble Pressure
[Pa]
- liquid volume at Bubble Point [
m³
] - vapour volume at Bubble Point [
m³
] - Gas composition at Bubble Point
By default, uses equality of chemical potentials, via ChemPotBubblePressure
Clapeyron.bubble_temperature
— Functionbubble_temperature(model::EoSModel, p, x,method::BubblePointMethod = ChemPotBubbleTemperature())
Calculates the bubble temperature and properties at a given pressure. Returns a tuple, containing:
- Bubble Temperature
[K]
- liquid volume at Bubble Point [
m³
] - vapour volume at Bubble Point [
m³
] - Gas composition at Bubble Point
By default, uses equality of chemical potentials, via ChemPotBubbleTemperature
Clapeyron.dew_pressure
— Functiondew_pressure(model::EoSModel, T, y,method = ChemPotDewPressure())
Calculates the dew pressure and properties at a given temperature. Returns a tuple, containing:
- Dew Pressure
[Pa]
- liquid volume at Dew Point [
m³
] - vapour volume at Dew Point [
m³
] - Liquid composition at Dew Point
By default, uses equality of chemical potentials, via ChemPotDewPressure
Clapeyron.dew_temperature
— Functiondew_temperature(model::EoSModel, p, y, method = ChemPotDewTemperature())
calculates the dew temperature and properties at a given pressure. Returns a tuple, containing:
- Dew Temperature
[K]
- liquid volume at Dew Point [
m³
] - vapour volume at Dew Point [
m³
] - Liquid composition at Dew Point
By default, uses equality of chemical potentials, via ChemPotDewTemperature
Clapeyron.azeotrope_pressure
— Functionazeotrope_pressure(model::EoSModel, T; v0 = x0_azeotrope_pressure(model,T))
calculates the azeotrope pressure and properties at a given temperature. Returns a tuple, containing:
- Azeotrope Pressure
[Pa]
- liquid volume at Azeotrope Point [
m³
] - vapour volume at Azeotrope Point [
m³
] - Azeotrope composition
Clapeyron.azeotrope_temperature
— Functionazeotrope_temperature(model::EoSModel, T; v0 = x0_bubble_pressure(model,T,[0.5,0.5]))
Calculates the azeotrope temperature and properties at a given pressure. Returns a tuple, containing:
- Azeotrope Temperature
[K]
- liquid volume at Azeotrope Point [
m³
] - vapour volume at Azeotrope Point [
m³
] - Azeotrope composition
Clapeyron.LLE_pressure
— FunctionLLE_pressure(model::EoSModel, T, x; v0 = x0_LLE_pressure(model,T,x))
calculates the Liquid-Liquid equilibrium pressure and properties at a given temperature.
Returns a tuple, containing:
- LLE Pressure
[Pa]
- liquid volume of composition
x₁ = x
at LLE Point [m³
] - liquid volume of composition
x₂
at LLE Point [m³
] - Liquid composition
x₂
Clapeyron.LLE_temperature
— FunctionLLE_temperature(model::EoSModel, p, x; T0 = x0_LLE_temperature(model,p,x))
calculates the Liquid-Liquid equilibrium temperature and properties at a given pressure.
Returns a tuple, containing:
- LLE Pressure
[Pa]
- liquid volume of composition
x₁ = x
at LLE Point [m³
] - liquid volume of composition
x₂
at LLE Point [m³
] - Liquid composition
x₂
Clapeyron.VLLE_pressure
— FunctionVLLE_pressure(model::EoSModel, T; v0 = x0_LLE_pressure(model,T))
calculates the Vapor-Liquid-Liquid equilibrium pressure and properties of a binary mixture at a given temperature.
Returns a tuple, containing:
- VLLE Pressure
[Pa]
- Liquid volume of composition
x₁
at VLLE Point [m³
] - Liquid volume of composition
x₂
at VLLE Point [m³
] - Vapour volume of composition
y
at VLLE Point [m³
] - Liquid composition
x₁
- Liquid composition
x₂
- Liquid composition
y
Clapeyron.VLLE_temperature
— FunctionVLLE_temperature(model::EoSModel, p; T0 = x0_LLE_temperature(model,p))
calculates the Vapor-Liquid-Liquid equilibrium temperature and properties of a binary mixture at a given temperature.
Returns a tuple, containing:
- VLLE temperature
[K]
- Liquid volume of composition
x₁
at VLLE Point [m³
] - Liquid volume of composition
x₂
at VLLE Point [m³
] - Vapour volume of composition
y
at VLLE Point [m³
] - Liquid composition
x₁
- Liquid composition
x₂
- Liquid composition
y
Clapeyron.crit_mix
— Functioncrit_mix(model::EoSModel,z;v0=x=x0_crit_mix(model,z))
Returns the critical mixture point at a ginven composition.
Returns a tuple, containing:
- Critical Mixture Temperature
[K]
- Critical Mixture Pressure
[Pa]
- Critical Mixture Volume
[m³]
Clapeyron.UCEP_mix
— FunctionUCEP_mix(model::EoSModel;v0=x0_UCEP_mix(model))
Calculates the Upper Critical End Point of a binary mixture.
returns:
- UCEP Temperature [
K
] - UCEP Pressure [
Pa
] - liquid volume at UCEP Point [
m³
] - vapour volume at UCEP Point [
m³
] - liquid molar composition at UCEP Point
- vapour molar composition at UCEP Point
Clapeyron.UCST_mix
— FunctionUCST_mix(model::EoSModel,T;v0=x0_UCST_mix(model,T))
Calculates the Upper critical solution point of a mixture at a given Temperature.
returns:
- UCST Pressure [
Pa
] - volume at UCST Point [
m³
] - molar composition at UCST Point
Clapeyron.gibbs_solvation
— Functiongibbs_solvation(model::EoSModel, T; threaded=true, vol0=(nothing,nothing))
Calculates the solvation free energy as:
g_solv = -R̄*T*log(K)
where the first component is the solvent and second is the solute.
Clapeyron.cross_second_virial
— Functioncross_second_virial(model,T,z)
Default units: [m^3]
Calculates the second cross virial coefficient (B₁₂) of a binary mixture, using the definition:
B̄ = x₁^2*B₁₁ + 2x₁x₂B₁₂ + x₂^2*B₂₂
B₁₂ = (B̄ - x₁^2*B₁₁ - x₂^2*B₂₂)/2x₁x₂
!!! info composition-dependent The second cross virial coefficient calculated from a equation of state can present a dependency on composition [1], but normally, experiments for obtaining the second virial coefficient are made by mixing the same volume of two gases. you can calculate B12 in this way by using (Clapeyron.equivolcrosssecond_virial)[@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.equivol_cross_second_virial
— Functionequivol_cross_second_virial(model::EoSModel,T,p_exp = 200000.0)
calculates the second cross virial coefficient, by simulating the mixing of equal volumes of pure gas, at T,P conditions. The equal volume of each pure gas sets an specific molar amount for each component. Details of the experiment can be found at [1].
Example
model = SAFTVRQMie(["helium","neon"])
B12 = equivol_cross_second_virial(model,)
References
- Brewer, J., & Vaughn, G. W. (1969). Measurement and correlation of some interaction second virial coefficients from − 125° to 50°C. I. The Journal of Chemical Physics, 50(7), 2960–2968. doi:10.1063/1.1671491
Clapeyron.sle_solubility
— Functionsle_solubility(model::CompositeModel, p, T, z; solute)
Calculates the solubility of each component within a solution of the other components, at a given temperature and composition. Returns a matrix containing the composition of the SLE phase boundary for each component. If solute
is specified, returns only the solubility of the specified component.
Can only function when solid and fluid models are specified within a CompositeModel.
sle_solubility(model::CompositeModel, p, T, z; solute)
Calculates the solubility of each component within a solution of the other components, at a given temperature and composition. Returns a matrix containing the composition of the SLE phase boundary for each component. If solute
is specified, returns only the solubility of the specified component.
Can only function when solid and fluid models are specified within a CompositeModel.
Clapeyron.slle_solubility
— Functionslle_solubility(model::CompositeModel, p, T)
Calculates the phase boundary for solid-liquid-liquid equilibriumm of a ternary mixture, at a given temperature and pressure. Returns a matrix containing the composition of the two liquids phases.
Can only function when solid and liquid models are specified within a CompositeModel and when the third component is the solute.
Clapeyron.eutectic_point
— Functioneutectic_point(model::CompositeModel, p)
Calculates the eutectic point of a binary mixture (at a given pressure). Returns a tuple containing the eutectic temperature and the composition at the eutectic point.
Can only function when solid and liquid models are specified within a CompositeModel.
Bubble/Dew methods
Clapeyron.ChemPotBubblePressure
— TypeChemPotBubblePressure(kwargs...)
Function to compute bubble_pressure
via chemical potentials. It directly solves the equality of chemical potentials system of equations.
Inputs:
y0 = nothing
: optional, initial guess for the vapor phase compositionp0 = nothing
: optional, initial guess for the bubble pressure [Pa
]vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesatol = 1e-8
: optional, absolute tolerance of the non linear system of equationsrtol = 1e-12
: optional, relative tolerance of the non linear system of equationsmax_iters = 10000
: optional, maximum number of iterationsnonvolatiles = nothing
: optional, Vector of strings containing non volatile compounds. those will be set to zero on the vapour phase.
Clapeyron.FugBubblePressure
— TypeFugBubblePressure(kwargs...)
Function to compute bubble_pressure
via fugacity coefficients. First it uses successive substitution to update the phase composition and a outer newtown loop to update the pressure. If no convergence is reached after itmax_newton
iterations, the system is solved using a multidimensional non-linear system of equations.
Inputs:
y0 = nothing
: optional, initial guess for the vapor phase compositionp0 = nothing
: optional, initial guess for the bubble pressure [Pa
]vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesitmax_newton = 10
: optional, number of iterations to update the pressure using newton's methoditmax_ss = 5
: optional, number of iterations to update the liquid phase composition using successive substitutiontol_x = 1e-8
: optional, tolerance to stop successive substitution cycletol_p = 1e-8
: optional, tolerance to stop newton cycletol_of = 1e-8
: optional, tolerance to check if the objective function is zero.nonvolatiles = nothing
: optional, Vector of strings containing non volatile compounds. those will be set to zero on the vapour phase.
Clapeyron.ActivityBubblePressure
— TypeActivityBubblePressure(kwargs...)
Function to compute bubble_pressure
using Activity Coefficients. On activity coefficient models it solves the problem via succesive substitucion. On helmholtz-based models, it uses the Chapman approximation for activity coefficients.
Inputs:
gas_fug = true
: if the solver uses gas fugacity coefficients. onActivityModel
is set by default tofalse
poynting = true
: if the solver use the poynting correction on the liquid fugacity coefficients. onActivityModel
is set by default tofalse
y0 = nothing
: optional, initial guess for the vapor phase compositionp0 = nothing
: optional, initial guess for the bubble pressure [Pa
]vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesatol = 1e-8
: optional, absolute tolerance of the non linear system of equationsrtol = 1e-12
: optional, relative tolerance of the non linear system of equationsitmax_ss = 40
: optional, maximum number of sucesive substitution iterations
Clapeyron.ChemPotBubbleTemperature
— TypeChemPotBubbleTemperature(kwargs...)
Function to compute bubble_temperature
via chemical potentials. It directly solves the equality of chemical potentials system of equations.
Inputs:
y = nothing
: optional, initial guess for the vapor phase composition.T0 = nothing
: optional, initial guess for the bubble temperature [K
].vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesatol = 1e-8
: optional, absolute tolerance of the non linear system of equationsrtol = 1e-12
: optional, relative tolerance of the non linear system of equationsmax_iters = 10000
: optional, maximum number of iterationsnonvolatiles = nothing
: optional, Vector of strings containing non volatile compounds. those will be set to zero on the vapour phase.
Clapeyron.FugBubbleTemperature
— TypeFugBubbleTemperature(kwargs...)
Method to compute bubble_temperature
via fugacity coefficients. First it uses successive substitution to update the phase composition and a outer newtown loop to update the temperature. If no convergence is reached after itmax_newton
iterations, the system is solved using a multidimensional non-linear system of equations.
Inputs:
y = nothing
: optional, initial guess for the vapor phase composition.T0 = nothing
: optional, initial guess for the bubble temperature [K
].vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesitmax_newton = 10
: optional, number of iterations to update the temperature using newton's methoditmax_ss = 5
: optional, number of iterations to update the liquid phase composition using successive substitutiontol_x = 1e-8
: optional, tolerance to stop successive substitution cycletol_T = 1e-8
: optional, tolerance to stop newton cycletol_of = 1e-8
: optional, tolerance to check if the objective function is zero.nonvolatiles
: optional, Vector of strings containing non volatile compounds. those will be set to zero on the vapour phase.
Clapeyron.ActivityBubbleTemperature
— TypeActivityBubbleTemperature(kwargs...)
Function to compute bubble_temperature
using Activity Coefficients. On activity coefficient models it solves the problem via succesive substitucion. On helmholtz-based models, it uses the Chapman approximation for activity coefficients.
Inputs:
gas_fug = true
: if the solver uses gas fugacity coefficients. onActivityModel
is set by default tofalse
poynting = true
: if the solver use the poynting correction on the liquid fugacity coefficients. onActivityModel
is set by default tofalse
y0 = nothing
: optional, initial guess for the vapor phase compositionT0 = nothing
: optional, initial guess for the bubble temperature [K
]vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesatol = 1e-8
: optional, absolute tolerance of the non linear system of equationsrtol = 1e-12
: optional, relative tolerance of the non linear system of equationsitmax_ss = 40
: optional, maximum number of sucesive substitution iterations
Clapeyron.ChemPotDewPressure
— TypeChemPotDewPressure(kwargs...)
Function to compute dew_pressure
via chemical potentials. It directly solves the equality of chemical potentials system of equations.
Inputs:
x0 = nothing
: optional, initial guess for the liquid phase compositionp0 = nothing
: optional, initial guess for the dew pressure [Pa
]vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesatol = 1e-8
: optional, absolute tolerance of the non linear system of equationsrtol = 1e-12
: optional, relative tolerance of the non linear system of equationsmax_iters = 10000
: optional, maximum number of iterationsnoncondensables = nothing
: optional, Vector of strings containing non condensable compounds. those will be set to zero on the liquid phase.
Clapeyron.FugDewPressure
— TypeFugDewPressure(kwargs...)
Method to compute dew_pressure
via fugacity coefficients. First it uses successive substitution to update the phase composition and a outer newtown loop to update the pressure. If no convergence is reached after itmax_newton
iterations, the system is solved using a multidimensional non-linear system of equations.
Inputs:
x0 = nothing
: optional, initial guess for the liquid phase compositionp0 = nothing
: optional, initial guess for the dew pressure [Pa
]vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesitmax_newton = 10
: optional, number of iterations to update the pressure using newton's methoditmax_ss = 5
: optional, number of iterations to update the liquid phase composition using successive substitutiontol_x = 1e-8
: optional, tolerance to stop successive substitution cycletol_p = 1e-8
: optional, tolerance to stop newton cycletol_of = 1e-8
: optional, tolerance to check if the objective function is zero.noncondensables = nothing
: optional, Vector of strings containing non condensable compounds. those will be set to zero on the liquid phase.
Clapeyron.ActivityDewPressure
— TypeActivityDewPressure(kwargs...)
Function to compute dew_pressure
using Activity Coefficients. On activity coefficient models it solves the problem via succesive substitucion. On helmholtz-based models, it uses the Chapman approximation for activity coefficients.
Inputs:
gas_fug = true
: if the solver uses gas fugacity coefficients. onActivityModel
is set by default tofalse
poynting = true
: if the solver use the poynting correction on the liquid fugacity coefficients. onActivityModel
is set by default tofalse
x0 = nothing
: optional, initial guess for the liquid phase compositionp0 = nothing
: optional, initial guess for the dew pressure [Pa
]vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesatol = 1e-8
: optional, absolute tolerance of the non linear system of equationsrtol = 1e-12
: optional, relative tolerance of the non linear system of equationsitmax_ss = 40
: optional, maximum number of sucesive substitution iterations
Clapeyron.ChemPotDewTemperature
— TypeChemPotDewTemperature(kwargs...)
Function to compute dew_temperature
via chemical potentials. It directly solves the equality of chemical potentials system of equations.
Inputs:
x0 = nothing
: optional, initial guess for the liquid phase compositionT0 =nothing
: optional, initial guess for the dew temperature [K
]vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesatol = 1e-8
: optional, absolute tolerance of the non linear system of equationsrtol = 1e-12
: optional, relative tolerance of the non linear system of equationsmax_iters = 10000
: optional, maximum number of iterationsnoncondensables = nothing
: optional, Vector of strings containing non condensable compounds. those will be set to zero on the liquid phase.
Clapeyron.FugDewTemperature
— TypeFugDewTemperature(kwargs...)
Method to compute dew_temperature
via fugacity coefficients. First it uses successive substitution to update the phase composition and a outer newtown loop to update the temperature. If no convergence is reached after itmax_newton
iterations, the system is solved using a multidimensional non-linear system of equations.
Inputs:
x0 = nothing
: optional, initial guess for the liquid phase compositionT0 = nothing
: optional, initial guess for the dew temperature [K
]vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesitmax_newton = 10
: optional, number of iterations to update the temperature using newton's methoditmax_ss = 5
: optional, number of iterations to update the liquid phase composition using successive substitutiontol_x = 1e-8
: optional, tolerance to stop successive substitution cycletol_T = 1e-8
: optional, tolerance to stop newton cycletol_of = 1e-8
: optional, tolerance to check if the objective function is zero.noncondensables = nothing
: optional, Vector of strings containing non condensable compounds. those will be set to zero on the liquid phase.
Clapeyron.ActivityDewTemperature
— TypeActivityDewTemperature(kwargs...)
Function to compute dew_temperature
using Activity Coefficients. On activity coefficient models it solves the problem via succesive substitucion. On helmholtz-based models, it uses the Chapman approximation for activity coefficients.
Inputs:
gas_fug = true
: if the solver uses gas fugacity coefficients. onActivityModel
is set by default tofalse
poynting = true
: if the solver use the poynting correction on the liquid fugacity coefficients. onActivityModel
is set by default tofalse
x0 = nothing
: optional, initial guess for the liquid phase compositionT0 = nothing
: optional, initial guess for the dew temperature [K
]vol0 = nothing
: optional, initial guesses for the liquid and vapor phase volumesatol = 1e-8
: optional, absolute tolerance of the non linear system of equationsrtol = 1e-12
: optional, relative tolerance of the non linear system of equationsitmax_ss = 40
: optional, maximum number of sucesive substitution iterations
Consistency and Stability
Clapeyron.gibbs_duhem
— Functiongibbs_duhem(model,V,T,z=[1.0])
performs a Gibbs-Duhem check on the input conditions:
∑zᵢμᵢ - G ≈ 0
Where G
is the total gibbs energy. it can help diagnose if a user-defined eos is consistent.
return |∑zᵢμᵢ - G|, ∑zᵢμᵢ and G at the specified conditions.
Clapeyron.isstable
— Functionisstable(model,V,T,z)::Bool
Performs stability tests for a (V,T,z) pair, and warn if any tests fail. returns true/false
.
Checks:
- mechanical stability: isothermal compressibility is not negative.
- diffusive stability: all eigenvalues of
∂²A/∂n²
are positive. - chemical stability: there isn't any other combinations of compositions at p(V,T),T that are more stable than the input composition.
Clapeyron.VT_mechanical_stability
— FunctionVT_mechanical_stability(model,V,T,z = SA[1.0])::Bool
Performs a mechanical stability for a (V,T,z) pair, returns true/false
. Checks if isothermal compressibility is not negative.
This function does not have a p
,T
counterpart, because if we calculate the volume via volume(model,p,T,z)
, it will be, by definition, a mechanically stable phase.
Clapeyron.VT_diffusive_stability
— FunctionVT_diffusive_stability(model,V,T,z)::Bool
Performs a diffusive stability for a (V,T,z) pair, returns true/false
. Checks if all eigenvalues of ∂²A/∂n²
are positive.
Clapeyron.VT_chemical_stability
— FunctionVT_chemical_stability(model,V,T,z)::Bool
Performs a chemical stability check using the tangent plane distance criterion, using the tpd function
Clapeyron.tpd
— Functiontpd(model,p,T,z;break_first = false,lle = false,tol_trivial = 1e-5, di = nothing)
Calculates the Tangent plane distance function (tpd
). It returns:
- a vector with trial phase compositions where
tpd < 0
- a vector with the
tpd
values - a vector with symbols indicating the phase of the input composition
- a vector with symbols indicating the phase of the trial composition
It iterates over each two-phase combination, starting from pure trial compositions, it does succesive substitution, then Gibbs optimization.
If the vectors are empty, then the procedure couldn't find a negative tpd
. That is an indication that the phase is (almost) surely stable.
TP Flash
Clapeyron.tp_flash
— Functiontp_flash(model, p, T, n, method::TPFlashMethod = DETPFlash())
Routine to solve non-reactive multicomponent flash problem. The default method uses Global Optimization. see DETPFlash
Inputs:
- T, Temperature
- p, Pressure
- n, vector of number of moles of each species
Outputs - Tuple containing:
- xᵢⱼ, Array of mole fractions of species j in phase i
- nᵢⱼ, Array of mole numbers of species j in phase i, [mol]
- G, Gibbs Free Energy of Equilibrium Mixture [J]
Clapeyron.DETPFlash
— TypeDETPFlash(; numphases = 2,
max_steps = 1e4*(numphases-1),
population_size =20,
time_limit = Inf,
verbose = false,
logspace = false,
equilibrium = :auto)
Method to solve non-reactive multicomponent flash problem by finding global minimum of Gibbs Free Energy via Differential Evolution.
User must assume a number of phases, numphases
. If true number of phases is smaller than numphases, model should predict either (a) identical composition in two or more phases, or (b) one phase with negligible total number of moles. If true number of phases is larger than numphases, a thermodynamically unstable solution will be predicted.
The optimizer will stop at max_steps
evaluations or at time_limit
seconds
The equilibrium
keyword allows to restrict the search of phases to just liquid-liquid equilibria (equilibrium = :lle
). the default searches for liquid and gas phases.
Clapeyron.RRTPFlash
— TypeRRTPFlash{T}(;kwargs...)
Method to solve non-reactive multicomponent flash problem by Rachford-Rice equation.
Only two phases are supported. if K0
is nothing
, it will be calculated via the Wilson correlation.
Keyword Arguments:
- equilibrium = equilibrium type ":vle" for liquid vapor equilibria, ":lle" for liquid liquid equilibria
K0
(optional), initial guess for the constants Kx0
(optional), initial guess for the composition of phase xy0
= optional, initial guess for the composition of phase yvol0
= optional, initial guesses for phase x and phase y volumesK_tol
= tolerance to stop the calculationmax_iters
= number of Successive Substitution iterations to performnacc
= accelerate successive substitution method every nacc steps. Should be a integer bigger than 3. Set to 0 for no acceleration.noncondensables
= arrays with names (strings) of components non allowed on the liquid phase. In the case of LLE equilibria, corresponds to thex
phasenonvolatiles
= arrays with names (strings) of components non allowed on the vapour phase. In the case of LLE equilibria, corresponds to they
phase
Clapeyron.MichelsenTPFlash
— TypeMichelsenTPFlash{T}(;kwargs...)
Method to solve non-reactive multicomponent flash problem by Michelsen's method.
Only two phases are supported. if K0
is nothing
, it will be calculated via the Wilson correlation.
Keyword Arguments:
- equilibrium = equilibrium type ":vle" for liquid vapor equilibria, ":lle" for liquid liquid equilibria
K0
(optional), initial guess for the constants Kx0
(optional), initial guess for the composition of phase xy0
= optional, initial guess for the composition of phase yvol0
= optional, initial guesses for phase x and phase y volumesK_tol
= tolerance to stop the calculationss_iters
= number of Successive Substitution iterations to performnacc
= accelerate successive substitution method every nacc steps. Should be a integer bigger than 3. Set to 0 for no acceleration.second_order
= wheter to solve the gibbs energy minimization using the analytical hessian or notnoncondensables
= arrays with names (strings) of components non allowed on the liquid phase. In the case of LLE equilibria, corresponds to thex
phasenonvolatiles
= arrays with names (strings) of components non allowed on the vapour phase. In the case of LLE equilibria, corresponds to they
phase
Clapeyron.MultiPhaseTPFlash
— TypeMultiPhaseTPFlash(;kwargs...)
Method to solve non-reactive multiphase (np
phases), multicomponent (nc
components) flash problem.
The flash algorithm uses successive stability tests to find new phases [1], and then tries to solve the system via rachford-rice and succesive substitution for nc * np * ss_iters
iterations.
If the Rachford-Rice SS fails to converge, it proceeds to solve the system via gibbs minimization in VT-space using lnK-β-ρ as variables [3].
The algorithm finishes when SS or the gibbs minimization converges and all resulting phases are stable.
If the result of the phase equilibria is not stable, then it proceeds to add/remove phases again, for a maximum of phase_iters
iterations.
Keyword Arguments:
K0
(optional), initial guess for the constants Kx0
(optional), initial guess for the composition of phase xy0
(optional), initial guess for the composition of phase yn0
(optional), initial guess for all compositions. it can be a matrix or a vector of vectors.K_tol = sqrt(eps(Float64))
, tolerance to stop the calculation (norm(lnK,1) < K_tol
)ss_iters = 4
, number of Successive Substitution iterations to performnacc = 3
, accelerate successive substitution method every nacc steps. Should be a integer bigger than 3. Set to 0 for no acceleration.second_order = true
, whether to solve the gibbs energy minimization using the analytical hessian or not. If set tofalse
, the gibbs minimization will be done using L-BFGS.full_tpd
= false, whether to start with a simple K-split or using an intensive TPD search first.max_phases = typemax(Int)
, the algorithm stops if there are more thanmin(max_phases,nc)
phasesphase_iters = 20
, the maximum number of solve-add/remove-phase iterations
References
- Thermopack - Thermodynamic equilibrium algorithms reimplemented in a new framework. (2020, September 08). https://github.com/thermotools/thermopack. Retrieved May 4, 2024, from https://github.com/thermotools/thermopack/blob/main/docs/memo/flash/flash.pdf
- Okuno, R., Johns, R. T. T., & Sepehrnoori, K. (2010). A new algorithm for Rachford-Rice for multiphase compositional simulation. SPE Journal, 15(02), 313–325. doi:10.2118/117752-pa
- Adhithya, T. B., & Venkatarathnam, G. (2021). New pressure and density based methods for isothermal-isobaric flash calculations. Fluid Phase Equilibria, 537(112980), 112980. doi:10.1016/j.fluid.2021.112980
Clapeyron.MCFlashJL
— TypeMCFlashJL(; method = MultiComponentFlash.SSIFlash(), storage = nothing, V = NaN, K = nothing, kwargs.... )
Uses MultiComponentFlash.jl
two-phase flash solver. allows passing storage to minimize allocations. That storage can be created by calling MultiComponentFlash.flash_storage(model,p,T,z,method::MCFlashJL)
This method requires MultiComponentFlash
to be loaded in the current session (using MultiComponentFlash
) and julia >= v1.9
Clapeyron.numphases
— Functionnumphases(method::TPFlashMethod)
Return the number of phases supported by the TP flash method. by default its set to 2. it the method allows it, you can set the number of phases by doing method(;numphases = n)
.
Clapeyron.supports_reduction
— Functionsupports_reduction(method::TPFlashMethod)::Bool
Checks if a TP Flash method supports index reduction (the ability to prune model components with compositions lower than a threshold). All current Clapeyron methods support index reduction, but some methods that alllow passing a cache could have problems