Contents

Index

Primitive functions

Almost all models in Clapeyron based on helmholtz free energy have at least one of the following functions defined:

Clapeyron.eosFunction
eos(model::EoSModel, V, T, z=SA[1.0])

Returns the total Helmholtz free energy.

Inputs:

  • model::EoSModel Thermodynamic model to evaluate
  • V Total volume, in [m³]
  • T Temperature, in [K]
  • z mole amounts, in [mol], by default is @SVector [1.0]

Outputs:

  • Total Helmholtz free 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 free energy
Clapeyron.eos_resFunction
eos_res(model::EoSModel, V, T, z=SA[1.0])

Returns the residual Helmholtz free energy.

Inputs:

  • model::EoSModel Thermodynamic model to evaluate
  • V Total volume, in [m³]
  • T Temperature, in [K]
  • z mole amounts, in [mol], by default is @SVector [1.0]

Outputs:

  • Residual Helmholtz free 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.idealmodelFunction
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
true
Clapeyron.a_resFunction
a_res(model::EoSModel, V, T, z,args...)

Reduced residual Helmholtz free energy.

Inputs:

  • model::EoSModel Thermodynamic model to evaluate
  • V Total volume, in [m³]
  • T Temperature, in [K]
  • z mole amounts, in [mol], by default is @SVector [1.0]

Outputs:

  • Residual Helmholtz free energy, no units

You can define your own EoS by adding a method to a_res that accepts your custom model.

Automatic Differenciation functions

All bulk properties in Clapeyron are calculated via a combination of these Automatic Differenciation Primitives over eos or eos_res

Clapeyron.∂f∂TFunction
∂f∂T(model,V,T,z=SA[1.0])

returns f and ∂f/∂T at constant total volume and composition, where f is the total helmholtz energy, given by eos(model,V,T,z)

Clapeyron.∂f∂VFunction
∂f∂V(model,V,T,z=SA[1.0])

returns f and ∂f/∂V at constant temperature and composition, where f is the total helmholtz energy, given by eos(model,V,T,z), and V is the total volume

Clapeyron.∂fFunction
∂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∂VFunction
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.∂2fFunction
∂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.∂2pFunction
∂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_hessFunction
f_hess(model,V,T,z)

returns the second order volume (V) and temperature (T) derivatives of the total helmholtz energy (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.∂²³fFunction
∂²³f(model,V,T,z=SA[1.0])

returns ∂²A/∂V² and ∂³A/∂V³, in a single ForwardDiff pass. used mainly in crit_pure objective function

Thermodynamic Method Dispatch types

Clapeyron.ThermodynamicMethodType
ThermodynamicMethod

Abstract 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 superancillary

Bubble 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 points
Clapeyron.SaturationMethodType
SaturationMethod <: ThermodynamicMethod

Abstract type for saturation_temperature and saturation_pressure routines. Should at least support passing the crit keyword, containing the critical point, if available.

Clapeyron.BubblePointMethodType
BubblePointMethod <: ThermodynamicMethod

Abstract type for bubble_pressure and bubble_temperature routines.

Should at least support passing the y0 keyword, containing an initial vapour phase, if available.

Clapeyron.DewPointMethodType
DewPointMethod <: ThermodynamicMethod

Abstract type for dew_pressure and dew_temperature routines.

Should at least support passing the x0 keyword, containing an initial vapour phase, if available.

Clapeyron.TPFlashMethodType
TPFlashMethod <: ThermodynamicMethod

Abstract type for tp_flash routines. it requires defining numphases(method) and tp_flash_impl(model,p,T,n,method). If the method accept component-dependent inputs, it also should define index_reduction(method,nonzero_indices)

Reference States

Clapeyron.ReferenceStateType
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.
  • :ashrae: h = s = 0 at -40C saturated liquid
  • :iir: h = 200.0 kJ/kg, s=1.0 kJ/kg/K at 0C saturated liquid
  • :nbp: h = s = 0 at 1 atm saturated liquid

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)
  • :saturation_pressure h = H0, s = S0, at T = T0, saturated phase (specified by the phase argument)
  • :saturation_temperature h = H0, s = S0, at p = P0, saturated phase (specified by the phase argument)

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-9
Clapeyron.reference_stateFunction
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_stateFunction
has_reference_state(model)::Bool

Checks 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))
true

Note that the default idealmodel (BasicIdeal) does not allow for setting reference states.