Contents
Index
- Clapeyron.AssocOptions
- Clapeyron.AssocParam
- Clapeyron.ECS
- Clapeyron.FluidCorrelation
- Clapeyron.GammaPhi
- Clapeyron.GroupParam
- Clapeyron.PairParam
- Clapeyron.ParamOptions
- Clapeyron.SingleParam
- Clapeyron.SiteParam
- Clapeyron.EoSFunctions.bmcs_hs
- Clapeyron.EoSFunctions.xlogx
- Clapeyron.ParamTable
- Clapeyron.cite
- Clapeyron.cleartemp!
- Clapeyron.doi
- Clapeyron.doi2bib
- Clapeyron.epsilon_HudsenMcCoubrey
- Clapeyron.epsilon_HudsenMcCoubrey!
- Clapeyron.epsilon_HudsenMcCoubreysqrt
- Clapeyron.epsilon_HudsenMcCoubreysqrt!
- Clapeyron.epsilon_LorentzBerthelot
- Clapeyron.epsilon_LorentzBerthelot!
- Clapeyron.export_model
- Clapeyron.getparams
- Clapeyron.group_matrix
- Clapeyron.group_pairmean
- Clapeyron.group_sum
- Clapeyron.index_reduction
- Clapeyron.is_splittable
- Clapeyron.kij_mix
- Clapeyron.lambda_LorentzBerthelot
- Clapeyron.lambda_LorentzBerthelot!
- Clapeyron.lambda_squarewell
- Clapeyron.lambda_squarewell!
- Clapeyron.mirror_pair
- Clapeyron.mirror_pair!
- Clapeyron.pair_mix
- Clapeyron.sigma_LorentzBerthelot
- Clapeyron.sigma_LorentzBerthelot!
- Clapeyron.split_model
- Clapeyron.split_model_binaries
Parsing Parameters from Files
Clapeyron.ParamOptions — TypeParamOptions(;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- getparamsdisplayed 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#1and- #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#1and- #site#2, where- #site#is the value of this option.
- 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- ewill 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.getparams — Functionparams = 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 true, getparams will add a "sites" value in the params result, containing a SiteParam built with the input parameters.
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,850my_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  700This 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- @REPLACEwill 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, QThe 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: fittedNote, 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.ParamTable — FunctionParamTable(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:
- :singlecreates a table with single parameters
- :paircreates a table with pair parameters
- :assoccreates a table with association parameters
- :groupcreates 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.042Clapeyron.cleartemp! — Functioncleartemp!()Deletes all files in the temporary Clapeyron scratch space, used to store the csvs created by ParamTable.
Parameter types
Clapeyron.SingleParam — TypeSingleParam{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" => falseExample 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
endClapeyron.PairParam — TypePairParam{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.0Example 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
endClapeyron.AssocParam — TypeAssocParam{T}Struct holding association parameters.
Clapeyron.GroupParam — TypeGroupParamStruct 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
- n_intragroups: a list containining the connectivity graph (as a matrix) between each group for each component.
- 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" => 4In this case, SAFTGammaMie files support the second order group CH2OH.
Clapeyron.SiteParam — TypeSiteParamStruct 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.AssocOptions — TypeAssocOptions(;rtol = 1e-12,atol = 1e-12,max_iters = 1000,dampingfactor = 0.5,combining =:nocombining)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:- Δ[i,j][a,b] = √(Δ[i,i][a,b] * Δ[j,j][a,b])
- :dufal,- mie15: combining rule used for- SAFTVRMie15- ε[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])^3 / 8
Combining Rules, out-of-place methods
Clapeyron.kij_mix — Functionkij_mix(f,p::ClapeyronParam,k::PairParam)::PairParam
kij_mix(f,p::ClapeyronParam)::PairParamGeneral 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_mix — Functionpair_mix(g,P::ClapeyronParam,Q::ClapeyronParam)::PairParam
pair_mix(g,P::ClapeyronParam,Q::ClapeyronParam)::PairParamGeneral 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_pair — Functionmirror_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_LorentzBerthelot — Functionsigma_LorentzBerthelot(σ::SingleOrPair,ζ)::PairParam
sigma_LorentzBerthelot(σ::SingleOrPair)::PairParam
sigma_LorentzBerthelot(σ::Union{AbstractVector,AbstractMatrix},ζ)::AbstractMatrix
sigma_LorentzBerthelot(σ::Union{AbstractVector,AbstractMatrix})::AbstractMatrixCombining rule for a single or pair parameter. Returns a pair parameter with non diagonal entries equal to:
σᵢⱼ = (1 - ζᵢⱼ)*(σᵢ + σⱼ)/2If ζᵢⱼ is not defined, the definition is reduced to a simple arithmetic mean:
σᵢⱼ = (σᵢ + σⱼ)/2Ignores non-diagonal entries already set.
If a Single Parameter (or vector) is passed as input, it will be converted to a Pair Parameter with σᵢᵢ = σᵢ.
Clapeyron.epsilon_LorentzBerthelot — Functionepsilon_LorentzBerthelot(ϵ::SingleOrPair,k)::PairParam
epsilon_LorentzBerthelot(ϵ::SingleOrPair)::PairParam
epsilon_LorentzBerthelot(ϵ::Union{AbstractVector,AbstractMatrix},k)::AbstractMatrix
epsilon_LorentzBerthelot(ϵ::Union{AbstractVector,AbstractMatrix})::AbstractMatrixCombining 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_HudsenMcCoubrey — Functionepsilon_HudsenMcCoubrey(ϵ::SingleOrPair,σ)::PairParam
epsilon_HudsenMcCoubrey(ϵ::SingleOrPair)::PairParam
epsilon_HudsenMcCoubrey(ϵ::Union{AbstractVector,AbstractMatrix},σ)::AbstractMatrix
epsilon_HudsenMcCoubrey(ϵ::Union{AbstractVector,AbstractMatrix})::AbstractMatrixCombining 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.epsilon_HudsenMcCoubreysqrt — Functionepsilon_HudsenMcCoubreysqrt(ϵ::SingleOrPair,σ)::PairParam
epsilon_HudsenMcCoubreysqrt(ϵ::SingleOrPair)::PairParam
epsilon_HudsenMcCoubreysqrt(ϵ::Union{AbstractVector,AbstractMatrix},σ)::AbstractMatrix
epsilon_HudsenMcCoubreysqrt(ϵ::Union{AbstractVector,AbstractMatrix})::AbstractMatrixCombining rule for a single or pair parameter. Returns a pair parameter with non diagonal entries equal to:
ϵᵢⱼ = √(ϵᵢϵⱼ * σᵢᵢ^3 * σⱼⱼ^3)/σᵢⱼ^3If σᵢⱼ 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_LorentzBerthelot — Functionlambda_LorentzBerthelot(λ::SingleOrPair,k = 3)::PairParam
lambda_LorentzBerthelot(λ::Union{AbstractVector,AbstractMatrix},k = 3)::AbstractMatrixCombining 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_squarewell — Functionlambda_squarewell(λ::SingleOrPair,σ)::PairParam
lambda_squarewell(λ::Union{AbstractVector,AbstractMatrix},σ)::AbstractMatrixCombining 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 λᵢᵢ = λᵢ.
Combining Rules, in-place methods
Missing docstring for Clapeyron.kij_mix!. Check Documenter's build log for details.
Missing docstring for Clapeyron.pair_mix!. Check Documenter's build log for details.
Clapeyron.mirror_pair! — Functionmirror_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_LorentzBerthelot! — Functionsigma_LorentzBerthelot!(σ::PairParameter,ζ)::PairParam
sigma_LorentzBerthelot!(σ::AbstractMatrix,ζ)::AbstractMatrixCombining rule for a single or pair parameter. Returns a pair parameter with non diagonal entries equal to:
σᵢⱼ = (1 - ζᵢⱼ)*(σᵢ + σⱼ)/2If ζᵢⱼ is not defined, the definition is reduced to a simple arithmetic mean:
σᵢⱼ = (σᵢ + σⱼ)/2The method overwrites the entries in σ, with the exception of diagonal entries.
Clapeyron.epsilon_LorentzBerthelot! — Functionepsilon_LorentzBerthelot!(ϵ::PairParameter,k)::PairParam
epsilon_LorentzBerthelot!(ϵ::PairParameter)::PairParam
epsilon_LorentzBerthelot!(ϵ::AbstractMatrix,k)::AbstractMatrix
epsilon_LorentzBerthelot!(ϵ::AbstractMatrix)::AbstractMatrixCombining 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:
ϵᵢⱼ = √(ϵᵢϵⱼ)The method overwrites the entries in ϵ, with the exception of diagonal entries.
Clapeyron.epsilon_HudsenMcCoubrey! — Functionepsilon_HudsenMcCoubrey!(ϵ::PairParameter,σ)::PairParam
epsilon_HudsenMcCoubrey!(ϵ::PairParameter)::PairParam
epsilon_HudsenMcCoubrey!(ϵ::AbstractMatrix,σ)::AbstractMatrix
epsilon_HudsenMcCoubrey!(ϵ::AbstractMatrix)::AbstractMatrixCombining 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:
ϵᵢⱼ = √(ϵᵢϵⱼ)The method overwrites the entries in ϵ, with the exception of diagonal entries.
Clapeyron.epsilon_HudsenMcCoubreysqrt! — Functionepsilon_HudsenMcCoubreysqrt!(ϵ::PairParameter,σ)::PairParam
epsilon_HudsenMcCoubreysqrt!(ϵ::PairParameter)::PairParam
epsilon_HudsenMcCoubreysqrt!(ϵ::AbstractMatrix,σ)::AbstractMatrix
epsilon_HudsenMcCoubreysqrt!(ϵ::AbstractMatrix)::AbstractMatrixCombining rule for a single or pair parameter. Returns a pair parameter with non diagonal entries equal to:
ϵᵢⱼ = √(ϵᵢϵⱼ * σᵢᵢ^3 * σⱼⱼ^3)/σᵢⱼ^3If σᵢⱼ is not defined, the definition is reduced to a simple geometric mean:
ϵᵢⱼ = √(ϵᵢϵⱼ)The method overwrites the entries in ϵ, with the exception of diagonal entries.
Clapeyron.lambda_LorentzBerthelot! — Functionlambda_LorentzBerthelot!(λ::PairParameter,k = 3)::PairParam
lambda_LorentzBerthelot!(λ::AbstractMatrix,k = 3)::AbstractMatrixCombining rule for a pair parameters. Returns a pair parameter or matrix with non diagonal entries equal to:
λᵢⱼ = k + √((λᵢᵢ - k)(λⱼⱼ - k))with k = 0 the definition is reduced to a simple geometric mean:
λᵢⱼ = √(λᵢλⱼ)The method overwrites the entries in λ, with the exception of diagonal entries.
Clapeyron.lambda_squarewell! — Functionlambda_squarewell!(λ::PairParameter,σ)::PairParam
lambda_squarewell!(λ::AbstractMatrix,σ)::AbstractMatrixCombining rule for pair parameters. Returns a pair parameter with non diagonal entries equal to:
λᵢⱼ = (σᵢᵢλᵢᵢ + σⱼⱼλⱼⱼ)/(σᵢᵢ + σⱼⱼ)Ignores non-diagonal entries already set.
The method overwrites the entries in λ, with the exception of diagonal entries.
Group Combining Rules
Clapeyron.group_sum — Functiongroup_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_pairmean — Functiongroup_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 Pis a single paremeter, thenP(i,j) = f(P[i],P[j])
- if Pis a pair paremeter, thenP(i,j) = p[i,j]
Clapeyron.group_matrix — Functiongroup_matrix(groups::MixedGCSegmentParam)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_model — Functionsplit_model(model::EoSModel)
split_model(model::EoSModel, splitter)Takes in a model for a multi-component system and returns a vector of models. The result depends on the splitter used.
A model can be splitted in a list of submodels, where each submodel is of the same type as the original model, but it has a different combination of components.
Group-Contribution models are also splitted in a component basis, split_model takes care of converting between Group and Component basis automatically.
A splitter is just a list of indices for each submodel, valid splitters are:
- A list of integers: split_model(model,[1,5,2])will return three pure models.
- An integer: split_model(model,1)will return a list with one model corresponding to the first component
- A list of lists: split_model(model,[[1,2],[3])will return a list with two models, the first one will contain two components, the second one will be a pure model.
The default splitter is 1:length(model), that will return a list with all pure models.
Examples
julia> model = MonomerIdeal(["methane","propane","butane"])
MonomerIdeal with 3 components:
 "methane"
 "propane"
 "butane"
Contains parameters: Mw, reference_state
julia> split_model(model)
3-element Vector{MonomerIdeal}:
 MonomerIdeal("methane")
 MonomerIdeal("propane")
 MonomerIdeal("butane")
 
julia> split_model(model,[[1,2],[3,1]])
2-element Vector{MonomerIdeal}:
 MonomerIdeal("methane", "propane")
 MonomerIdeal("butane", "methane")
julia> split_model(model,2)
1-element Vector{MonomerIdeal}:
 MonomerIdeal("propane")Clapeyron.is_splittable — Functionis_splittable(model)::BoolTrait 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_reduction — Functionindex_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.
Clapeyron.split_model_binaries — Functionsplit_model_binaries(model::EoSModel)::Vector{EoSModel}Given a multicomponent EoSModel, returns a list with the combination of all binary models.
Example
julia> model = PCPSAFT(["water", "methanol", "propyleneglycol","methyloxirane"])
PCPSAFT{BasicIdeal} with 4 components:
 "water"
 "methanol"
 "propyleneglycol"
 "methyloxirane"
Contains parameters: Mw, segment, sigma, epsilon, dipole, dipole2, epsilon_assoc, bondvol
julia> split_model_binaries(model)
6-element Vector{PCPSAFT{BasicIdeal}}:
 PCPSAFT{BasicIdeal}("water", "methanol")
 PCPSAFT{BasicIdeal}("water", "propyleneglycol")
 PCPSAFT{BasicIdeal}("water", "methyloxirane")
 PCPSAFT{BasicIdeal}("methanol", "propyleneglycol")
 PCPSAFT{BasicIdeal}("methanol", "methyloxirane")
 PCPSAFT{BasicIdeal}("propyleneglycol", "methyloxirane")Math/Utility Functions
Clapeyron.EoSFunctions.bmcs_hs — Functionbmcs_hs(ζ0,ζ1,ζ2,ζ3)Boublík-Mansoori-Leeland-Carnahan-Starling hard sphere term:
1/ζ₀ * (3ζ₁*ζ₂/(1-ζ₃) + ζ₂^3/(ζ₃*(1-ζ₃)^2) + (ζ₂^3/ζ₃^2-ζ₀)*log(1-ζ₃))Clapeyron.EoSFunctions.xlogx — Functionxlogx(x::Real,k = one(x))Returns x * log(k*x) for x ≥ 0, handling $x = 0$ by taking the downward limit.
copied from LogExpFunctions.jl
Model Exporting
Clapeyron.export_model — Functionexport_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).
Clapeyron.FluidCorrelation — TypeFluidCorrelation{V,L,Sat,Cp} <: RestrictedEquilibriaModelWrapper struct to signal that a CompositeModel uses correlations for calculation of saturation points, vapour and liquid phase volumes.
Clapeyron.GammaPhi — TypeGammaPhi{γ,Φ} <: RestrictedEquilibriaModelWrapper struct to signal that a CompositeModel uses an activity model in conjunction with a fluid.
Missing docstring for Clapeyron.ExtendedCorrespondingStates. Check Documenter's build log for details.
Clapeyron.ECS — TypeExtendedCorrespondingStates <: EoSModel
function ECS(components,
    refmodel=PropaneRef(),
    shapemodel=SRK(components),
    shaperef = SRK(refmodel.components))Input Models
- shape_model: shape model
- shape_ref: shape reference. Is the same type of EoS that- shape_model.
- model_ref: Reference model
Description
A Extended Corresponding states method.
The idea is to use a "shape model" that provides a corresponding states parameters and a "reference model" that implements a Helmholtz energy function, so that:
eos(shape_model,v,T,x)/RT = eos(model_ref,v₀,T₀)/RT₀where:
T₀ = T/f
v₀ = v/h
f,h = shape_factors(model::ECS,shape_ref::EoSModel,V,T,z)shape_factors can be used to create custom Extended Corresponding state models.
References
.1 Mollerup, J. (1998). Unification of the two-parameter equation of state and the principle of corresponding states. Fluid Phase Equilibria, 148(1–2), 1–19. doi:10.1016/s0378-3812(98)00230-1
Model Citing
Clapeyron.doi — Functiondoi(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.cite — Functioncite(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.doi2bib — Functiondoi2bib(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}
}"