Contents

Index

Parsing Parameters from Files

Clapeyron.ParamOptionsType
ParamOptions(;kwargs...)

Struct containing all the options related to parameter parsing:

  • userlocations = String[]: List of used-defined locations to search.
  • group_userlocations = String[]: List of used-defined group locations to search.
  • asymmetricparams::Vector{String} = String[]: List of pair parameters that follow that param[i,j] ≠ param[j,i]. if not set on asymmetric pairs, the asymmetric values will be overwritten!
  • ignore_headers::Vector{String} = ["dipprnumber", "smiles"]: List of ignored headers.
  • ignore_missing_singleparams::Vector{String} = String[]: List of parameters where checking for missing parameter values (in SingleParam) or the diagonal (on PairParam) are ignored.
  • verbose::Bool = false: If true, show all operations done by getparams displayed in the terminal. this includes the warnings emmited by CSV.jl
  • species_columnreference::String ="species": column name to check for components. in pair and association params, it will check for #species#1 and #species#2, where #species# is the value of this option.
  • site_columnreference::String ="site": column name to check for sites in association params, it will check for #site#1 and #site#2, where #site# is the value of this option.
  • group_columnreference::String = "groups": column name to check for groups in group data.
  • normalisecomponents::Bool = true: If true, performs normalization of strings, on the CSV and input components
  • n_sites_columns::Dict{String,String} = Dict( "e" => "n_e","e1" => "n_e1","e2" => "n_e2","H" => "n_H"): dictionary to look number of sites. the number of sites is stored as columns in a single parameter csv file. for example, the number of sites of name e will be looked on the column n_e
  • return_sites::Bool = true: If set to false, association params will be ignored and sites will not be created, even if they exist in the list of locations.
  • component_delimiter::String = "~|~": When there are multiple component names to match, seperate them by this delimiter.
Clapeyron.getparamsFunction
params, sites = getparams(components,locations;kwargs...)

returns a Dict{String,ClapeyronParam} containing all the parameters found for the list of components in the available CSVs. locations are the locations relative to Clapeyron database. the available keywords are the ones used ∈ ParamOptions if return_sites is set to false, getparams will only return the found params.

Single to Pair promotion

When reading multiple CSVs, if a parameter name appears in a single paramter file and in a pair parameter file, the single parameter values will be promoted to be the diagonal values of the pair interaction matrix:

my_parameter_single.csv

Clapeyron Database File
like parameters
species,a
sp1,1000
sp2,700
sp3,850

my_parameter_pair.csv

Clapeyron Database File
pair parameters
species1,species2,a
sp1,sp2,875
sp2,sp3,792
sp3,sp1,960

julia> res = getparams(["sp1","sp2"],userlocations = [my_parameter_single.csv,my_parameter_pair.csv])
Dict{String, Clapeyron.ClapeyronParam} with 1 entry:
  "a" => PairParam{Int64}("a")["sp1", "sp2"]

julia> res["a"].values
2×2 Matrix{Int64}:
 1000  875
  875  700

This promotion fails only happens in Single-Pair combinations. it fails otherwise.

In-memory CSV parsing

if you pass any string starting with Clapeyron Database File, it will be parsed as a CSV instead of being used as a filepath:

julia> x = """Clapeyron Database File,
       in memory like parameters
       species,a,b
       sp1,1000,0.05
       sp2,700,0.41
       """
"Clapeyron Database File,
in memory parameters [csvtype = like,grouptype = in_memory_read]
species,a,b
sp1,1000,0.05
sp2,700,0.41
"
julia> Clapeyron.getparams(["sp1","sp2"],userlocations = [x])
Dict{String, Clapeyron.ClapeyronParam} with 2 entries:
  "b" => SingleParam{Float64}("b")["sp1", "sp2"]
  "a" => SingleParam{Int64}("a")["sp1", "sp2"]

Special prefixes

There are some special prefixes that are used by the parser to signal some specific behaviour to be done at parsing time, for one CSV or a group of them:

  • @DB: replaces the path by the current Clapeyron default database. When doing getparams(components,["location"]), the paths are lowered to getparams(components,userlocations = ["@DB/location"]).

In a way, is a path shortcut used internally by Clapeyron to parse it's own database. you can change the path where @DB points to (or add other path shortcuts), via adding a corresponding entry to the Clapeyron.SHORT_PATHS Dict.

  • @REPLACE: Any filepath starting with @REPLACE will clear all previous appearances of the parameter names found in the CSV that contains the prefix.
  • @REMOVEDEFAULTS: it is used alone, and needs to be passed at the first position of the vector of userlocations. it will skip parsing of the default parameters:

The effect of the the parser can be summarized by the following examples:

model = PCSAFT(["water"],userlocations = ["@REMOVEDEFAULTS"]) #fails, no parameters found, no CSV parsed
model = PCSAFT(["water"],userlocations = ["@REPLACE/empty_params.csv"]) #fails, no parameters found, default parameters parsed and then removed
model = PCSAFT(["water"],userlocations = ["@REPLACE/my_pcsaft_kij.csv"]) #success, default kij parameters replaced by the ones on `my_pcsaft_kij.csv`
model = PCSAFT(["water"],userlocations = ["@REMOVEDEFAULTS","@DB/SAFT/PCSAFT","@DB/properties/molarmass.csv"]) #sucess. default parameters csv removed, and parsed again, using the @DB prefix to point to the default database.

You can use the @REPLACE keyword in a in-memory CSV by adding it at the start of the string, followed by an space:

#This will replace all previous parsed occurences of `a` and `b`
x_replace = """@REPLACE Clapeyron Database File,
in memory like parameters
species,a,b
sp1,1000,0.05
sp2,700,0.41
"""

CSV type detection and group type

The second line of the csv is used for comments and to identify the type of CSV used. for example:

x = """Clapeyron Database File
       in memory like parameters
       species,a,b
       sp1,1000,0.05
       sp2,700,0.41
       """

Will be parsed as a table with single parameter data. if you want more flexibility, you can instead pass the csvtype between brackets:

x = """Clapeyron Database File
       i can write anything here, unlike, association [csvtype = like] but the csv type is already specified.
       species,a,b
       sp1,1000,0.05
       sp2,700,0.41
       """

additionaly, there are some cases when you want to absolutely sure that your types don't clash with the default values. this is the case with different group parametrizations of UNIFAC (Dormund, VTPR, PSRK):

julia> model = UNIFAC(["methanol","ethanol"])
UNIFAC{PR{BasicIdeal, PRAlpha, NoTranslation, vdW1fRule}} with 2 components:
 "methanol": "CH3OH" => 1
 "ethanol": "CH2" => 1, "CH3" => 1, "OH (P)" => 1
Group Type: UNIFACDortmund
Contains parameters: A, B, C, R, Q

julia> model = PSRKUNIFAC(["methanol","ethanol"])
UNIFAC{BasicIdeal} with 2 components:
 "methanol": "CH3OH" => 1
 "ethanol": "CH2" => 1, "CH3" => 1, "OH" => 1
Group Type: PSRK
Contains parameters: A, B, C, R, Q

The models are the same (UNIFAC), but the group parametrizations are different. this is specified with the grouptype keyword. for example, if we see UNIFAC_groups.csv, it starts with:

Clapeyron Database File,
modified UNIFAC (Dortmund) Groups [csvtype = groups,grouptype = UNIFACDortmund]
species,groups
ethane,"[""CH3"" => 2]"
propane,"[""CH3"" => 2, ""CH2"" => 1]"
butane,"[""CH3"" => 2, ""CH2"" => 2]"
...

For compatibility reasons, if you pass a CSV without grouptype, it will be accepted, but two CSV with different specified group types cannot be merged:

x1 = """Clapeyron Database File
       paramterization 1 [csvtype = like,grouptype = param1]
       species,a,b
       sp1,1000,0.05
       sp2,700,0.41
       """
x2 = """Clapeyron Database File
       fitted to data [csvtype = like,grouptype = fitted]
       species,a,b
       sp1,912,0.067
       sp2,616,0.432
       """

If we pass the same parameters, with different group types, the parser will fail

julia> Clapeyron.getparams(["sp1","sp2"],userlocations = [x1,x2])
ERROR: cannot join two databases with different group types:
current group type: param1
incoming group type: fitted

Note, that the parser will not fail if you pass different parameters with different group types (For example if a has param1 group type and b has fit group type)

Creating Files from Parameters

Clapeyron.ParamTableFunction
ParamTable(type::Symbol,table;
location = nothing,
name = nothing,
grouptype = :unknown,
options = ParamOptions())

Creates a clapeyron CSV file and returns the location of that file. the type determines the table type:

  • :single creates a table with single parameters
  • :pair creates a table with pair parameters
  • :assoc creates a table with association parameters
  • :group creates a table with association parameters

By default, the name is generated randomly, and the table is stored as a temporary scratch space (provided by Scratch.jl). You can clean said scratch space by using Clapeyron.cleartemp!().

Examples:

julia> data = (species = ["water"],Mw = [18.03]) #it could be a Dict, a named tuple, or any Tables.jl compatible table
(species = ["water"], Mw = [18.9])
julia> file = ParamTable(:single,data ,name="water_new_mw")
"C:\Users\user\.julia\scratchspaces\7c7805af-46cc-48c9-995b-ed0ed2dc909a\ParamTables\singledata_water_new_mw.csv"
julia> model = PCSAFT(["water","methanol"],userlocations = [file])
PCSAFT{BasicIdeal} with 2 components:
 "water"
 "methanol"
Contains parameters: Mw, segment, sigma,
epsilon, epsilon_assoc, bondvol
julia> model.params.Mw
SingleParam{Float64}("Mw") with 2 components:
 "water" => 18.9
 "methanol" => 32.042
Clapeyron.cleartemp!Function
cleartemp!()

Deletes all files in the temporary Clapeyron scratch space, used to store the csvs created by ParamTable.

Parameter types

Clapeyron.SingleParamType
SingleParam{T}

Struct designed to contain single parameters. Basically a vector with some extra info.

Creation:

julia> mw = SingleParam("molecular weight",["water","ammonia"],[18.01,17.03])
SingleParam{Float64}("molecular weight") with 2 components:
 "water" => 18.01
 "ammonia" => 17.03
julia> mw.values
2-element Vector{Float64}:
 18.01
 17.03
julia> mw.components
2-element Vector{String}:
 "water"
 "ammonia"
julia> mw2 = SingleParam(mw,"new name")
SingleParam{Float64}("new name") with 2 components:
 "water" => 18.01
 "ammonia" => 17.03
julia> has_oxigen = [true,false]; has_o = SingleParam(mw2,has_oxigen)
SingleParam{Bool}("new name") with 2 components:
 "water" => true
 "ammonia" => false

Example usage in models:

function molecular_weight(model,molar_frac)
    mw = model.params.mw.values
    res = zero(eltype(molarfrac))
    for i in @comps #iterating through all components
        res += molar_frac[i]*mw[i]
    end
    return res
end
Clapeyron.PairParamType
PairParam{T}

Struct designed to contain pair data. used a matrix as underlying data storage.

Creation:

julia> kij = PairParam("interaction params",["water","ammonia"],[0.1 0.0;0.1 0.0])
PairParam{Float64}["water", "ammonia"]) with values:
2×2 Matrix{Float64}:
 0.1  0.0
 0.1  0.0
julia> kij.values
2×2 Matrix{Float64}:
 0.1  0.0
 0.1  0.0
julia> diagvalues(kij)
2-element view(::Vector{Float64}, 
1:3:4) with eltype Float64:
 0.1
 0.0

Example usage in models:

#lets compute ∑xᵢxⱼkᵢⱼ
function alpha(model,x)
    kij = model.params.kij.values
    ki = diagvalues(model.params.kij)
    res = zero(eltype(molarfrac))
    for i in @comps 
        @show ki[i] #diagonal values
        for j in @comps 
            res += x[i]*x[j]*kij[i,j]
        end
    end
    return res
end
Clapeyron.GroupParamType
GroupParam

Struct holding group parameters.contains:

  • components: a list of all components
  • groups: a list of groups names for each component
  • grouptype: used to differenciate between different group models.
  • i_groups: a list containing the number of groups for each component
  • n_groups: a list of the group multiplicity of each group corresponding to each group in i_groups
  • flattenedgroups: a list of all unique groups–the parameters correspond to this list
  • n_flattenedgroups: the group multiplicities corresponding to each group in flattenedgroups

You can create a group param by passing a `Vector{Tuple{String, Vector{Pair{String, Int64}}}}. For example:

julia> grouplist = [
           ("ethanol", ["CH3"=>1, "CH2"=>1, "OH"=>1]),
           ("nonadecanol", ["CH3"=>1, "CH2"=>18, "OH"=>1]),
           ("ibuprofen", ["CH3"=>3, "COOH"=>1, "aCCH"=>1, "aCCH2"=>1, "aCH"=>4])];
julia> groups = GroupParam(grouplist)
GroupParam(:unknown) with 3 components:
 "ethanol": "CH3" => 1, "CH2" => 1, "OH" => 1
 "nonadecanol": "CH3" => 1, "CH2" => 18, "OH" => 1
 "ibuprofen": "CH3" => 3, "COOH" => 1, "aCCH" => 1, "aCCH2" => 1, "aCH" => 4
julia> groups.flattenedgroups
7-element Vector{String}:
 "CH3"
 "CH2"
 "OH"
 "COOH"
 "aCCH"
 "aCCH2"
 "aCH"
julia> groups.i_groups
3-element Vector{Vector{Int64}}:
 [1, 2, 3]
 [1, 2, 3]
 [1, 4, 5, 6, 7]
julia> groups.n_groups
3-element Vector{Vector{Int64}}:
 [1, 1, 1]
 [1, 18, 1]
 [3, 1, 1, 1, 4]
julia> groups.n_flattenedgroups
 3-element Vector{Vector{Int64}}:
 [1, 1, 1, 0, 0, 0, 0]
 [1, 18, 1, 0, 0, 0, 0]
 [3, 0, 0, 1, 1, 1, 4]

if you have CSV with group data, you can also pass those, to automatically query the missing groups in your input vector:

julia> grouplist = [
           "ethanol",
           ("nonadecanol", ["CH3"=>1, "CH2"=>18, "OH"=>1]),
           ("ibuprofen", ["CH3"=>3, "COOH"=>1, "aCCH"=>1, "aCCH2"=>1, "aCH"=>4])];
           julia> groups = GroupParam(grouplist, ["SAFT/SAFTgammaMie/SAFTgammaMie_groups.csv"])
           GroupParam with 3 components:
            "ethanol": "CH2OH" => 1, "CH3" => 1
            "nonadecanol": "CH3" => 1, "CH2" => 18, "OH" => 1
            "ibuprofen": "CH3" => 3, "COOH" => 1, "aCCH" => 1, "aCCH2" => 1, "aCH" => 4

In this case, SAFTGammaMie files support the second order group CH2OH.

Clapeyron.SiteParamType
SiteParam

Struct holding site parameters. Is built by parsing all association parameters in the input CSV files. It has the following fields:

  • components: a list of all components (or groups in Group Contribution models)
  • sites: a list containing a list of all sites corresponding to each component (or group) in the components field
  • n_sites: a list of the site multiplicities corresponding to each site in flattenedsites
  • flattenedsites: a list of all unique sites
  • i_sites: an iterator that goes through the indices corresponding to each site in flattenedsites
  • n_flattenedsites: the site multiplicities corresponding to each site in flattenedsites
  • i_flattenedsites: an iterator that goes through the indices for each flattened site

Let's explore the sites in a 3-component SAFTGammaMie model:

julia> model3 = SAFTgammaMie([    
                "ethanol",
                ("nonadecanol", ["CH3"=>1, "CH2"=>18, "OH"=>1]),     
                ("ibuprofen", ["CH3"=>3, "COOH"=>1, "aCCH"=>1, "aCCH2"=>1, "aCH"=>4])
                               ])
SAFTgammaMie{BasicIdeal} with 3 components:
 "ethanol"
 "nonadecanol"
 "ibuprofen"
Contains parameters: segment, shapefactor, lambda_a, lambda_r, sigma, epsilon, epsilon_assoc, bondvol 
julia> model3.sites
SiteParam with 8 sites:
 "CH2OH": "H" => 1, "e1" => 2     
 "CH3": (no sites)
 "CH2": (no sites)
 "OH": "H" => 1, "e1" => 2        
 "COOH": "e2" => 2, "H" => 1, "e1" => 2
 "aCCH": (no sites)
 "aCCH2": (no sites)
 "aCH": (no sites)
julia> model3.sites.flattenedsites
3-element Vector{String}:
 "H"
 "e1"
 "e2"
julia> model3.sites.i_sites       
8-element Vector{Vector{Int64}}:
 [1, 2]
 []
 []
 [1, 2]
 [1, 2, 3]
 []
 []
 []
julia> model3.sites.n_sites       
8-element Vector{Vector{Int64}}:
 [1, 2]
 []
 []
 [1, 2]
 [2, 1, 2]
 []
 []
 []
Clapeyron.AssocOptionsType
AssocOptions(;rtol = 1e-12,atol = 1e-12,max_iters = 1000,dampingfactor = 0.5,combining =:nocombining,dense = true)

Struct containing iteration parameters for the solver of association sites. the combining option controls the type of combining rule applied to the association strength:

  • nocombining (default). Does not perform any combination rules.
  • :cr1: "combining rule - 1": ε[i,j][a,b] = (ε[i,i][a,b] + ε[j,j][a,b])/2 β[i,j][a,b] = √(β[i,i][a,b] * β[j,j][a,b])
  • :esd,:elliott: Elliott–Suresh–Donohue combining rule: ε[i,j][a,b] = (ε[i,i][a,b] + ε[j,j][a,b])/2 β[i,j][a,b] = √(β[i,i][a,b] * β[j,j][a,b]) * (σ[i]*σ[j]/σ[i,j])^3
  • :esd_runtime,:elliott_runtime: combining rule, performed at runtime. Dense matrix only. slowest, but it is the only way to handle temperature-dependent σ: Δ[i,j][a,b] = (ε[i,i][a,b] + ε[j,j][a,b])/2 Δ[i,j][a,b] = √(Δ[i,i][a,b] * Δ[j,j][a,b])
Association Scheme matters

all combining rules implicitly requires that both Δ(i,i,a,b) and Δ(j,j,a,b) are non-zero, that means that components that don't self associate will not be combined.

Combining Rules

Clapeyron.kij_mixFunction
kij_mix(f,p::ClapeyronParam,k::PairParam)::PairParam
kij_mix(f,p::ClapeyronParam)::PairParam

General combining rule for pair parameter with a kᵢⱼ interaction parameter. returns a pair parameter with non diagonal entries equal to:

pᵢⱼ = f(pᵢ,pⱼ,kᵢⱼ)

Where f is a 'combining' function that follows the rules:

f(pᵢ,pⱼ,0) = f(pⱼ,pᵢ,0)
f(pᵢ,pᵢ,0) = pᵢ

and k must follow:

kᵢᵢ = 0 

Ignores non-diagonal entries already set.

If a Single Parameter is passed as input, it will be converted to a Pair Parameter with pᵢᵢ = pᵢ.

Clapeyron.pair_mixFunction
pair_mix(g,P::ClapeyronParam,Q::ClapeyronParam)::PairParam
pair_mix(g,P::ClapeyronParam,Q::ClapeyronParam)::PairParam

General combining rule for a pair and a single parameter. returns a pair parameter P with non diagonal entries equal to:

Pᵢⱼ = g(Pᵢ,Pⱼ,Qᵢ,Qⱼ,Qᵢⱼ)

Where f is a 'combining' function that follows the rules:

Pᵢⱼ = Pⱼᵢ = g(Pᵢ,Pⱼ,Qᵢ,Qⱼ,Qᵢⱼ) = g(Pⱼ,Pᵢ,Qⱼ,Qᵢ,Qᵢⱼ)
g(Pᵢ,Pᵢ,Qᵢ,Qᵢ,Qᵢ) = Pᵢ

it is a more general form of kij_mix, where kij_mix(f,P,Q) == pair_mix(g,P,Q) is correct if:

f(Pᵢ,Pⱼ,Qᵢⱼ) = g(Pᵢ,Pⱼ,_,_,Qᵢⱼ)

If you pass a `SingleParam` or a vector as input for `Q`, then `Qᵢⱼ` will be considered 0.
Clapeyron.mirror_pairFunction
mirror_pair(param::PairParam,f = identity)

performs an operation f over the indices of p such as p[j,i] = f(p[i,j]). by default, f = identity (a symmetric matrix). One key difference is that it sets the ismissingvalues field for each modified index to false

Clapeyron.sigma_LorentzBerthelotFunction
sigma_LorentzBerthelot(σ::SingleOrPair,ζ::PairParam)::PairParam
sigma_LorentzBerthelot(σ::SingleOrPair)::PairParam

Combining rule for a single or pair parameter. returns a pair parameter with non diagonal entries equal to:

σᵢⱼ = (1 - ζᵢⱼ)*(σᵢ + σⱼ)/2

If ζᵢⱼ is not defined, the definition is reduced to a simple arithmetic mean:

σᵢⱼ = (σᵢ + σⱼ)/2

Ignores non-diagonal entries already set.

If a Single Parameter is passed as input, it will be converted to a Pair Parameter with σᵢᵢ = σᵢ.

Clapeyron.epsilon_LorentzBerthelotFunction
epsilon_LorentzBerthelot(ϵ::SingleOrPair,k::PairParam)::PairParam
epsilon_LorentzBerthelot(ϵ::SingleOrPair)::PairParam

Combining rule for a single or pair parameter. returns a pair parameter with non diagonal entries equal to:

ϵᵢⱼ = (1 - kᵢⱼ)*√(ϵᵢϵⱼ)

If kᵢⱼ is not defined, the definition is reduced to a simple geometric mean:

ϵᵢⱼ = √(ϵᵢϵⱼ)

Ignores non-diagonal entries already set.

If a Single Parameter is passed as input, it will be converted to a Pair Parameter with ϵᵢᵢ = ϵᵢ.

Clapeyron.epsilon_HudsenMcCoubreyFunction
epsilon_HudsenMcCoubrey(ϵ::SingleOrPair,σ::PairParam)::PairParam
epsilon_HudsenMcCoubrey(ϵ::SingleOrPair)::PairParam

Combining rule for a single or pair parameter. returns a pair parameter with non diagonal entries equal to:

ϵᵢⱼ = √(ϵᵢϵⱼ)*(σᵢᵢ^3 * σⱼⱼ^3)/σᵢⱼ^6 

If σᵢⱼ is not defined, the definition is reduced to a simple geometric mean:

ϵᵢⱼ = √(ϵᵢϵⱼ)

Ignores non-diagonal entries already set.

If a Single Parameter is passed as input, it will be converted to a Pair Parameter with ϵᵢᵢ = ϵᵢ.

Clapeyron.lambda_LorentzBerthelotFunction
lambda_LorentzBerthelot(λ::ClapeyronParam,k = 3)::PairParam

Combining rule for a single or pair parameter. returns a pair parameter with non diagonal entries equal to:

λᵢⱼ = k + √((λᵢᵢ - k)(λⱼⱼ - k))

with k = 0 the definition is reduced to a simple geometric mean:

ϵᵢⱼ = √(ϵᵢϵⱼ)

Ignores non-diagonal entries already set.

If a Single Parameter is passed as input, it will be converted to a Pair Parameter with λᵢᵢ = λᵢ.

Clapeyron.lambda_squarewellFunction
lambda_squarewell(λ::SingleOrPair,σ::PairParam)::PairParam

Combining rule for a single or pair parameter. returns a pair parameter with non diagonal entries equal to:

λᵢⱼ = (σᵢᵢλᵢᵢ + σⱼⱼλⱼⱼ)/(σᵢᵢ + σⱼⱼ)

Ignores non-diagonal entries already set.

If a Single Parameter is passed as input, it will be converted to a Pair Parameter with λᵢᵢ = λᵢ.

Group Combining Rules

Clapeyron.group_sumFunction
group_sum(groups::GroupParam,P::SingleParameter)

Given a GroupParam and a Single Parameter P for group data, it will return a single parameter p of component data, where:

pᵢ = ∑Pₖνᵢₖ

where νᵢₖ is the number of groups k at component i.

group_sum(groups::GroupParam,P::AbstractVector)

Given a GroupParam and a Vector P for group data, it will return a Vector p of component data, where:

pᵢ = ∑Pₖνᵢₖ

where νᵢₖ is the number of groups k at component i.

group_sum(groups::GroupParam,::Nothing)

Given a GroupParam, it will return a vector p of component data, where:

pᵢ = ∑νᵢₖ

where νᵢₖ is the number of groups k at component i.

Clapeyron.group_pairmeanFunction
group_pairmean(groups::GroupParam,param::PairParam)
group_pairmean(f,groups::GroupParam,param::SingleParam)

Given a GroupParam and a parameter P it will return a single parameter p of component data, where:

pᵢ = ∑νᵢₖ(∑(νᵢₗ*P(i,j))) / ∑νᵢₖ(∑νᵢₗ)

where νᵢₖ is the number of groups k at component i and P(i,j) depends on the type of P:

  • if P is a single paremeter, then P(i,j) = f(P[i],P[j])
  • if P is a pair paremeter, then P(i,j) = p[i,j]
Clapeyron.group_matrixFunction
group_matrix(groups::GroupParam)

returns a matrix of size (k,i) with values νₖᵢ. when multiplied with a molar amount, it returns the amount of moles of each group.

Model Splitting

Clapeyron.split_modelFunction
split_model(model::EoSModel)

Takes in a model for a multi-component system and returns a vector of model for each pure system.

Example:

julia> gerg2 = GERG2008(["propane","pentane"])
GERG008 model with 2 components:
"propane"
"pentane"
julia> split_model(gerg2)
2-element Vector{GERG2008}:
 GERG2008("propane")
 GERG2008("pentane")
Clapeyron.is_splittableFunction
is_splittable(model)::Bool

Trait to determine if a EoSModel should be splitted by itself or can be simply filled into a vector. This is useful in the case of models without any parameters, as those models are impossible by definition to split, because they don't have any underlying data. The Default is is_splittable(model) = true.

Clapeyron.index_reductionFunction
index_reduction(model::EoSModel,z,zmin = sum(z)*4*eps(eltype(z)))
index_reduction(model::EoSModel,bools <: AbstractVector{Bool})

Removes any component with composition z[i] < zmin. returns a reduced model model_r and a vector of indices idx_r, such as:

model_r,idx_r = index_reduction(model,z)
eos(model,V,T,z) ≈ eos(model_r,V,T,z[idx_r])

if the model does not have empty compositions, it will just return the input model.

The function will error if the reduction results in an empty model.

you can pass an arbitrary boolean vector (bools) to perform the reduction. ```

Model Exporting

Clapeyron.export_modelFunction
export_model(model::EoSModel,name="";location=".")

Exports model parameters to CSVs. Unless the name kwarg is specified, the name of the files will follow the convention singledata_EoS, pairdata_EoS and assocdata_EoS. Files will be saved within the current directory unless the location argument is specified.

Note that it will export all submodel parameters (e.g. Alpha function parameters for cubic EoS).

Model Citing

Clapeyron.doiFunction
doi(model)

Returns a Vector of strings containing the top-level bibliographic references of the model, in DOI format. if there isn't a references field, it defaults to default_references(model)

julia> umr = UMRPR(["water"],idealmodel = WalkerIdeal);Clapeyron.doi(umr)
1-element Vector{String}:
 "10.1021/I160057A011"
Clapeyron.citeFunction
cite(model,out = :doi)

Returns a Vector of strings containing all bibliographic references of the model, in the format indicated by the out argument. this includes any nested models.

julia> umr = UMRPR(["water"],idealmodel = WalkerIdeal);Clapeyron.cite(umr) #should cite UMRPR, UNIFAC, WalkerIdeal
4-element Vector{String}:
 "10.1021/I160057A011"
 "10.1021/ie049580p"
 "10.1021/i260064a004"
 "10.1021/acs.jced.0c00723"

the out argument supports two values:

  • :doi: returns the stored values on each EoS. by default those are DOI identifiers.
  • :bib: returns BibTeX entries. to use this, an internet connection is required.
julia> model = SAFTVRQMie(["helium"])
SAFTVRQMie{BasicIdeal} with 1 component:
 "helium"
Contains parameters: Mw, segment, sigma, lambda_a, lambda_r, epsilon

julia> Clapeyron.cite(model,:bib)
2-element Vector{String}:
 "@article{Aasen_2019,
	doi = {10" ⋯ 463 bytes ⋯ "Journal of Chemical Physics}
}"
 "@article{Aasen_2020,
	doi = {10" ⋯ 452 bytes ⋯ "Journal of Chemical Physics}
}"

This list will displayed by each EoSModel on future versions. you can enable/disable this by setting ENV["CLAPEYRON_SHOW_REFERENCES"] = "TRUE"/"FALSE"

Clapeyron.doi2bibFunction
doi2bib(doi::String)

given a DOI identifier, returns a BibTeX entry. requires an internet connection. if the value is not found, returns an empty string. results are cached in Clapeyron.DOI2BIB_CACHE

Example

julia> Clapeyron.doi2bib("10.1063/1.5136079")
"@article{Aasen_2020,
	doi = {10.1063/1.5136079},
	url = {https://doi.org/10.1063%2F1.5136079},
	year = 2020,
	month = {feb},
	publisher = {{AIP} Publishing},
	volume = {152},
	number = {7},
	pages = {074507},
	author = {Ailo Aasen and Morten Hammer and Erich A. Müller and {\O}ivind Wilhelmsen},
	title = {Equation of state and force fields for Feynman{\textendash}Hibbs-corrected Mie fluids. {II}. Application to mixtures of helium, neon, hydrogen, and deuterium},
	journal = {The Journal of Chemical Physics}
}"