Contents
Index
Clapeyron.BubblePointMethodClapeyron.DewPointMethodClapeyron.EoSModelClapeyron.ReferenceStateClapeyron.SaturationMethodClapeyron.TPFlashMethodClapeyron.ThermodynamicMethodClapeyron.RgasClapeyron.a_resClapeyron.eosClapeyron.eos_resClapeyron.eosshowClapeyron.f_hessClapeyron.has_groupsClapeyron.has_reference_stateClapeyron.has_sitesClapeyron.idealmodelClapeyron.p∂p∂2pClapeyron.p∂p∂VClapeyron.reference_stateClapeyron.set_reference_state!Clapeyron.∂2fClapeyron.∂2pClapeyron.∂fClapeyron.∂f∂TClapeyron.∂f∂V
Primitive functions
Almost all models in Clapeyron based on Helmholtz energy have at least one of the following functions defined:
Clapeyron.eos — Function
eos(model::EoSModel, V, T, z=SA[1.0])Returns the total Helmholtz energy.
Inputs:
model::EoSModelThermodynamic model to evaluateVTotal volume, in[m³]TTemperature, in[K]zmole amounts, in[mol], by default is@SVector [1.0]
Outputs:
- Total Helmholtz energy, in
[J].
By default, it calls R̄*T*∑(z)*(a_ideal(ideal_model,V,T,z) + a_res(model,V,T,z)) where ideal_model == idealmodel(model), where a_res is the reduced residual Helmholtz energy and a_ideal is the reduced ideal Helmholtz energy. You can mix and match ideal models if you provide:
[idealmodel](@ref)(model): extracts the ideal model from your Thermodynamic model.[a_res](@ref)(model,V,T,z): residual reduced Helmholtz energy.
Clapeyron.eos_res — Function
eos_res(model::EoSModel, V, T, z=SA[1.0])Returns the residual Helmholtz energy.
Inputs:
model::EoSModelThermodynamic model to evaluateVTotal volume, in[m³]TTemperature, in[K]zmole amounts, in[mol], by default is@SVector [1.0]
Outputs:
- Residual Helmholtz energy, in
[J].
By default, it calls R̄*T*∑(z)*(a_res(model,V,T,z)) where a_res is the reduced residual Helmholtz energy.
Clapeyron.idealmodel — Function
idealmodel(model::EoSModel)Retrieves the ideal model from the input's model. If the model is already an idealmodel, return nothing
Examples:
julia> pr = PR(["water"],idealmodel = MonomerIdeal)
PR{MonomerIdeal, PRAlpha, NoTranslation, vdW1fRule} with 1 component:
"water"
Contains parameters: a, b, Tc, Pc, Mw
julia> ideal = idealmodel(pr)
MonomerIdeal with 1 component:
"water"
Contains parameters: Mw
julia> idealmodel(ideal) == nothing
trueClapeyron.a_res — Function
a_res(model::EoSModel, V, T, z,args...)Returns reduced residual Helmholtz energy.
Inputs:
model::EoSModelThermodynamic model to evaluateVTotal volume, in[m³]TTemperature, in[K]zmole amounts, in[mol], by default is@SVector [1.0]
Outputs:
- Residual Helmholtz energy, no units.
You can define your own EoS by adding a method to a_res that accepts your custom model.
Clapeyron.eosshow — Function
eosshow(io::IO, model::EoSModel)
eosshow(io::IO, ::MIME"text/plain", model::EoSModel)Custom pretty-printer for EoSModel instances.
This is the backend used by Base.show for models that opt into the custom display. The text/plain variant prints components, parameters, reference state, and (optionally) citations when enabled via ENV["CLAPEYRON_SHOW_REFERENCES"].
Core types and utilities
Clapeyron.EoSModel — Type
EoSModelRoot type for all equation of state (EOS) models in Clapeyron.jl.
Subtypes of EoSModel represent specific equations of state (e.g., cubic, SAFT, activity coefficient models) and must implement their core interface methods.
By default, EoSModel subtypes are considered to be helmholtz-based EoS, and as such, they must implement their core interface methods:
a_res: An implementation of the reduced residual Helmholtz energylb_volume: Lower bound volume for the equation of stateT_scale: A temperature scaling factor.
Other EoSModel subtypes may have other interfaces (like activity models or volume-based models).
Clapeyron.Rgas — Function
Rgas(model)
Rgas()Returns the gas constant used by an EoSModel.
By default, uses the current 2019 definition: R̄ = 8.31446261815324 [J⋅K⁻¹⋅mol⁻¹]. you can call Rgas() to obtain this value.
Clapeyron.has_groups — Function
has_groups(model::EoSModel)
has_groups(::Type{<:EoSModel})Returns true when a model type carries group parameters.
This checks whether the model type has a groups field with GroupParam storage and is used to enable group-contribution logic.
Clapeyron.has_sites — Function
has_sites(model::EoSModel)
has_sites(::Type{<:EoSModel})Returns true when a model type carries site parameters.
This checks whether the model type has a sites field with SiteParam storage and is used to enable site-based logic (e.g. association models).
Automatic Differentiation functions
All bulk properties in Clapeyron are calculated via a combination of these Automatic Differentiation Primitives over eos or eos_res:
Clapeyron.∂f∂T — Function
∂f∂T(model,V,T,z=SA[1.0])Returns ∂f/∂T at constant total volume and composition, where f is the total Helmholtz energy, given by eos(model,V,T,z).
Clapeyron.∂f∂V — Function
∂f∂V(model,V,T,z)Returns ∂f/∂V at constant temperature T and composition z, where f is the total Helmholtz energy, given by eos(model,V,T,z), V is the total volume.
Clapeyron.∂f — Function
∂f(model,V,T,z)Returns zeroth order (value) and first order derivative information of the total Helmholtz energy (given by eos(model,V,T,z)). The result is given in two values:
grad_f,fval = ∂2f(model,V,T,z)where:
fval = f(V,T) = eos(model,V,T,z)
grad_f = [ ∂f/∂V; ∂f/∂T]Where V is the total volume, T is the temperature and f is the total Helmholtz energy.
Clapeyron.p∂p∂V — Function
p∂p∂V(model,V,T,z=SA[1.0])Returns p and ∂p/∂V at constant temperature, where p is the pressure = pressure(model,V,T,z) and V is the total volume.
Clapeyron.∂2f — Function
∂2f(model,V,T,z)Returns zeroth order (value), first order and second order derivative information of the total Helmholtz energy (given by eos(model,V,T,z)). The result is given in three values:
hess_f,grad_f,fval = ∂2f(model,V,T,z)where:
fval = f(V,T) = eos(model,V,T,z)
grad_f = [ ∂f/∂V; ∂f/∂T]
hess_f = [ ∂²f/∂V²; ∂²f/∂V∂T
∂²f/∂V∂T; ∂²f/∂V²]Where V is the total volume, T is the temperature and f is the total Helmholtz energy.
Clapeyron.∂2p — Function
∂2p(model,V,T,z)Returns zeroth order (value), first order and second order derivative information of the pressure. the result is given in three values:
hess_p,grad_p,pval = ∂2p(model,V,T,z)where:
pval = p(V,T) = pressure(model,V,T,z)
grad_p = [ ∂p/∂V; ∂p/∂T]
hess_p = [ ∂²p/∂V²; ∂²p/∂V∂T
∂²p/∂V∂T; ∂²p/∂V²]Where V is the total volume, T is the temperature and p is the pressure.
Clapeyron.f_hess — Function
f_hess(model,V,T,z)Returns the second order volume V and temperature T derivatives of the total Helmholtz energy f (given by eos(model,V,T,z)). The result is given in a 2x2 SMatrix, in the form:
[ ∂²f/∂V² ∂²f/∂V∂T
∂²f/∂V∂T ∂²f/∂T²]Use this instead of the ∂2f if you only need second order information. ∂2f also gives zeroth and first order derivative information, but due to a bug in the used AD, it allocates more than necessary.
Clapeyron.p∂p∂2p — Function
p∂p∂2p(model,V,T,z=SA[1.0])Returns the pressure p and their first and second volume derivatives ∂p/∂V and ∂²p/∂V², in a single ForwardDiff pass.
Thermodynamic Method Dispatch types
Clapeyron.ThermodynamicMethod — Type
ThermodynamicMethodAbstract type for all thermodynamic methods.
Normally, a thermodynamic method has the form: property(model,state..,method::ThermodynamicMethod). All methods used in this way subtype ThermodynamicMethod.
Examples
Saturation pressure:
model = PR(["water"])
Tsat = 373.15
saturation_pressure(model,Tsat) #using default method (chemical potential with volume base)
saturation_pressure(model,Tsat,SuperAncSaturation()) #solve using cubic superancillaryBubble point pressure
model = PCSAFT(["methanol","cyclohexane"])
T = 313.15
z = [0.5,0.5]
bubble_pressure(model,T,z) #using default method (chemical potential equality)
bubble_pressure(model,T,z,FugBubblePressure(y0 = = [0.6,0.4], p0 = 5e4)) #using isofugacity criteria with starting pointsClapeyron.SaturationMethod — Type
SaturationMethod <: ThermodynamicMethodAbstract type for saturation_temperature and saturation_pressure routines. Should at least support passing the crit keyword, containing the critical point, if available.
Clapeyron.BubblePointMethod — Type
BubblePointMethod <: ThermodynamicMethodAbstract type for bubble_pressure and bubble_temperature routines.
Should at least support passing the y0 keyword, containing an initial guess value for vapour phase, if available.
Clapeyron.DewPointMethod — Type
DewPointMethod <: ThermodynamicMethodAbstract type for dew_pressure and dew_temperature routines.
Should at least support passing the x0 keyword, containing an initial vapour phase, if available.
Clapeyron.TPFlashMethod — Type
TPFlashMethod <: ThermodynamicMethodAbstract type for tp_flash routines.
Reference States
Clapeyron.ReferenceState — Type
ReferenceState(type::Symbol = :no_set;T0 = NaN;P0 = NaN,H0 = NaN,S0 = NaN,phase = :unknown,z0 = Float64[])Parameter used to define a reference state for enthalpy and entropy, normally stored in the ideal model. When set, it calculates a set of a0 and a1 values such as the entropy and enthalpy at a specified point are fixed.
the type argument accepts the following standalone options:
:no_set: it returns the current defaults stablished by the equation of state. Leaves theReferenceStatestruct uninitialized.:zero: also returns the current defaults, but initializes the reference state struct, for later modification:ashrae: h = s = 0 at -40 °C saturated liquid:iir: h = 200.0 kJ·kg⁻¹, s=1.0 kJ·kg⁻¹·K⁻¹ at 0 °C saturated liquid:nbp: h = s = 0 at 1 atm saturated liquid:stp: h = s = 0 at 1 bar, 0 °C fluid of the most stable phase:stp_old: h = s = 0 at 1 atm, 0 °C fluid of the most stable phase:ntp: h = s = 0 at 1 atm, 20 °C fluid of the most stable phase
it also accepts the following options, that require additional specifications:
:volume: h = H0, s = S0, at T = T0, v =volume(model,P0,T0,z0,phase = phase):ideal_gas: h = H0, s = S0, at T = T0, v =volume(Clapeyron.idealmodel(model),P0,T0,z0):saturation_pressure: h = H0, s = S0, at T = T0, saturated phase (specified by thephaseargument):saturation_temperature: h = H0, s = S0, at p = P0, saturated phase (specified by thephaseargument)
If z0 is not specified, the reference state calculation will be done for each component separately.
Examples
julia> model = PCSAFT(["water","pentane"],idealmodel = ReidIdeal,reference_state = ReferenceState(:nbp))
PCSAFT{ReidIdeal, Float64} with 2 components:
"water"
"pentane"
Contains parameters: Mw, segment, sigma, epsilon, epsilon_assoc, bondvol
julia> model2 = PCSAFT(["water","pentane"],idealmodel = ReidIdeal,reference_state = :nbp) #equivalent
PCSAFT{ReidIdeal, Float64} with 2 components:
"water"
"pentane"
Contains parameters: Mw, segment, sigma, epsilon, epsilon_assoc, bondvol
julia> pure = split_model(model)
2-element Vector{PCSAFT{ReidIdeal, Float64}}:
PCSAFT{ReidIdeal, Float64}("water")
PCSAFT{ReidIdeal, Float64}("pentane")
julia> T,vl,_ = saturation_temperature(pure[1],101325.0) #saturated liquid at 1 atm
(373.2706553019503, 2.0512186595412677e-5, 0.03006573003253086)
julia> enthalpy(pure[1],101325.0,T)
-5.477897970382323e-6
julia> entropy(pure[1],101325.0,T)
5.009221069190994e-9Clapeyron.reference_state — Function
reference_state(model)::Union{ReferenceState,Nothing}Returns the reference state of the input model, if available. Returns nothing otherwise.
Examples
julia> reference_state(PCSAFT("water"))
false
julia> has_reference_state(PCSAFT("water",idealmodel = ReidIdeal))
true
julia> reference_state(PCSAFT("water",idealmodel = MonomerIdeal)) #has reference state, it is not initialized.
ReferenceState(String[], Float64[], Float64[], NaN, NaN, Float64[], Float64[], Float64[], :unknown, :no_set)
julia> reference_state(PCSAFT("water",idealmodel = MonomerIdeal, reference_state = ReferenceState(:nbp))) #has an initialized reference state
ReferenceState(["water"], [33107.133379491206], [17.225988924236503], NaN, NaN, [0.0], [0.0], [0.0], :unknown, :nbp)Clapeyron.has_reference_state — Function
has_reference_state(model)::BoolChecks if the input EoSModel has a reference state. Returns true/false.
Examples
julia> has_reference_state(PCSAFT("water"))
false
julia> has_reference_state(PCSAFT("water",idealmodel = ReidIdeal))
trueNote that the default idealmodel (BasicIdeal) does not allow for setting reference states.
Clapeyron.set_reference_state! — Function
set_reference_state!(model::EoSModel; verbose=false)
set_reference_state!(model, new_ref; verbose=false)Initialize and apply a reference state for model.
When a ReferenceState is provided, it is initialized (if needed) and then applied to the model. For mixtures, this may trigger per-component reference state initialization via split_pure_model.