Available Models¶
Fragmentation models implemented.
All models can be imported directly with:
from ugropy import joback, psrk, unifac
You can check the group list and their SMARTS representation with:
joback.subgroups psrk.subgroups unifac.subgroups
In the case of a PropertiesEstimator like joback, you can check the contribution of each group to the properties with:
joback.properties_contributions
Core¶
Core module.
FragmentationModel subgroups detection functions.
- check_can_fit_atoms(mol_object: Mol, mol_subgroups: dict, model: FragmentationModel) bool [source]¶
Check if a solution can be fitted in the mol_object atoms.
- Parameters:
mol_object (Chem.rdchem.Mol) – RDKit Mol object.
mol_subgroups (dict) – Subgroups of mol_object.
model (FragmentationModel) – FragmentationModel object.
- Returns:
True if the solution can be fitted.
- Return type:
bool
- check_has_composed(mol_subgroups: dict, model: FragmentationModel) tuple[bool, ndarray] [source]¶
Check if the molecule has composed structures.
A composed structure is a subgroup of FragmentationModel that can be decomposed into two or more FragmentationModel subgroups. For example, ACCH2 can be decomposed into the AC and CH2 groups.
- Parameters:
mol_subgroups (dict) – Dictionary with the detected subgroups.
model (FragmentationModel) – FragmentationModel object.
- Returns:
True if the molecule has composed structures.
- Return type:
bool
- draw(mol_object: Mol, subgroups: dict, model: FragmentationModel, title: str = '', width: float = 400, height: float = 200, title_font_size: float = 12, legend_font_size: float = 12, font: str = 'Helvetica') str [source]¶
Create a svg representation of the fragmentation result.
- Parameters:
mol_object (Chem.rdchem.Mol) – RDKit Mol object.
mol_subgroups (Union[dict, List[dict]]) – Subgroups of mol_object.
model (FragmentationModel) – FragmentationModel object.
title (str, optional) – Graph title, by default “”
width (int, optional) – Graph width, by default 400
height (int, optional) – Graph height, by default 200
title_font_size (int, optional) – Font size of graph’s title, by default 12
legend_font_size (int, optional) – Legend font size, by default 12
font (str, optional) – Text font, by default “Helvetica”
- Returns:
SVG string.
- Return type:
str
- check_has_hiden(mol_object: Mol, mol_subgroups: dict, model: FragmentationModel) bool [source]¶
Check for hidden subgroups in composed structures.
The principal subgroups that can be hidden in composed structures for the models UNIFAC, PSRK and Dortmund are CH2 and CH. The algorithm checks that the number of CH2 and CH groups in mol_subgroups dictionary is equal to the number of free CH2 and CH. If these numbers are not equal reveals that the CH2 and CH are hidden in composed structures, eg: ACCH2, ACCH. This phenomenon occurs when two subgroups fight for the same CH2 or CH. For example the molecule:
CCCC1=CC=C(COC(C)(C)C)C=C1
Here an ACCH2 and a CH2O are fighting to have the same CH2. But since there is a free CH2 in the molecule, the algorithm prefers to keep both ACCH2 and CH2O groups without any free CH2 subgroup. This check counts all the CH2 that are participating in a CH2 hideout (ACCH2 and CH2O are examples of hideouts). The algorithm notices that there is one free CH2 and there are zero free CH2 groups in the mol_subgroups dictionary and returns ‘True’ (mol_object has a hidden CH2).
- Parameters:
mol_object (Chem.rdchem.Mol) – RDKit Mol object.
mol_subgroups (dict) – Subgroups of mol_object.
model (FragmentationModel) – FragmentationModel object.
- Returns:
True if has hidden subgroups.
- Return type:
bool
- check_has_molecular_weight_right(mol_object: Mol, mol_subgroups: dict, model: FragmentationModel) bool [source]¶
Check the molecular weight of the molecule using its functional groups.
Compares the RDKit molecular weight of the molecule to the computed molecular weight from the functional groups. Returns True if both molecular weights are equal with 0.5 u (half hydrogen atom) as atol of numpy.allclose(). Also, the method will check if the molecule has negative occurrences on its functional groups, also returning False in that case.
- Parameters:
mol_object (Chem.rdchem.Mol) – RDKit Chem object
mol_subgroups (dict) – FragmentationModel subgroups of the mol_object
model (FragmentationModel) – FragmentationModel object.
- Returns:
True if RDKit and ugropy molecular weight are equal with a tolerance.
- Return type:
bool
- check_has_composed_overlapping(mol_object: Mol, mol_subgroups: dict, model: FragmentationModel) bool [source]¶
Check if in the solution are composed structures overlapping.
- Parameters:
mol_object (Chem.rdchem.Mol) – RDKit Mol object.
mol_subgroups (dict) – Subgroups of mol_object.
model (FragmentationModel) – FragmentationModel object.
- Returns:
Treu if the solution has overlapping composed structures.
- Return type:
bool
- correct_composed(mol_object: Mol, mol_subgroups: dict, model: FragmentationModel) dict [source]¶
Correct composed structures.
A priori is not easy to recognize what composed structures in mol_subgroups need to be decomposed to correct the solution. By that, all the combinations are tried. For example, a molecule that can’t be solved has one ACCH2 and two ACCH composed structures. The decomposition combinatory will be:
[[ACCH2], [ACCH], [ACCH2, ACCH], [ACCH, ACCH], [ACCH2, ACCH, ACCH]]
- Parameters:
mol_object (Chem.rdchem.Mol) – RDKit Mol object.
mol_subgroups (dict) – Molecule’s FragmentationModel subgroups.
model (FragmentationModel) – FragmentationModel object.
- Returns:
Corrected subgroups due to decomposing composed structures.
- Return type:
dict or list[dict]
- fit_atoms(mol_object: Mol, mol_subgroups: dict, model: FragmentationModel) dict [source]¶
Assign the atoms indexes for each mol_subgroup.
- Parameters:
mol_object (Chem.rdchem.Mol) – RDKit Mol object.
mol_subgroups (dict) – Subgroups of mol_object.
model (FragmentationModel) – FragmentationModel object.
- Returns:
Atom indexes in mol_object of each subgroup.
- Return type:
dict
- get_groups(model: FragmentationModel, identifier: str | Mol, identifier_type: str = 'name') Fragmentation [source]¶
Obtain the FragmentationModel’s subgroups of an RDkit Mol object.
- Parameters:
model (FragmentationModel) – FragmentationModel object.
identifier (str or rdkit.Chem.rdchem.Mol) – Identifier of a molecule (name, SMILES or Chem.rdchem.Mol). Example: hexane or CCCCCC.
identifier_type (str, optional) – Use ‘name’ to search a molecule by name, ‘smiles’ to provide the molecule SMILES representation or ‘mol’ to provide a rdkit.Chem.rdchem.Mol object, by default “name”.
- Returns:
FragmentationModel’s subgroups
- Return type:
Fragmentation
- instantiate_mol_object(identifier: str | Mol, identifier_type: str = 'name') Mol [source]¶
Instantiate a RDKit Mol object from molecule’s name or SMILES.
- Parameters:
identifier (str or rdkit.Chem.rdchem.Mol) – Identifier of a molecule (name, SMILES or rdkit.Chem.rdchem.Mol). Example: hexane or CCCCCC for name or SMILES respectively.
identifier_type (str, optional) – Use ‘name’ to search a molecule by name, ‘smiles’ to provide the molecule SMILES representation or ‘mol’ to provide a rdkit.Chem.rdchem.Mol object, by default “name”.
- Returns:
RDKit Mol object.
- Return type:
rdkit.Chem.rdchem.Mol
- correct_problematics(mol_object: Mol, mol_subgroups: dict, model: FragmentationModel) dict [source]¶
Correct problematic structures in mol_object.
- Parameters:
mol_object (Chem.rdchem.Mol) – RDKit Chem object
mol_subgroups (dict) – Dictionary with the subgroups not problematics corrected in mol_object.
model (FragmentationModel) – FragmentationModel object.
- Returns:
Molecule’s subrgoups corrected with the problematic structures list.
- Return type:
dict
Constants¶
constants module.
- R¶
Gas constant [J/mol/K]
- Type:
float
Groups¶
Groups module.
- class Groups(identifier: str, identifier_type: str = 'name', normal_boiling_temperature: float = None)[source]¶
Bases:
object
Group class.
Stores the solved FragmentationModels subgroups of a molecule.
- Parameters:
identifier (str or rdkit.Chem.rdchem.Mol) – Identifier of a molecule (name, SMILES or Chem.rdchem.Mol). Example: hexane or CCCCCC.
identifier_type (str, optional) – Use ‘name’ to search a molecule by name, ‘smiles’ to provide the molecule SMILES representation or ‘mol’ to provide a rdkit.Chem.rdchem.Mol object, by default “name”.
normal_boiling_temperature (float, optional) – If provided, will be used to estimate critical temperature, acentric factor, and vapor pressure instead of the estimated normal boiling point in the Joback group contribution model, by default None.
- identifier¶
Identifier of a molecule. Example: hexane or CCCCCC.
- Type:
str
- identifier_type¶
Use ‘name’ to search a molecule by name or ‘smiles’ to provide the molecule SMILES representation, by default “name”.
- Type:
str, optional
- mol_object¶
RDKit Mol object.
- Type:
rdkit.Chem.rdchem.Mol
- molecular_weight¶
Molecule’s molecular weight from rdkit.Chem.Descriptors.MolWt [g/mol].
- Type:
float
- unifac¶
Classic LV-UNIFAC subgroups.
- Type:
Fragmentation
- psrk¶
Predictive Soave-Redlich-Kwong subgroups.
- Type:
Fragmentation
- joback¶
JobackProperties object that contains the Joback subgroups and the estimated properties of the molecule.
- Type:
Properties¶
Properties module.
- class JobackProperties(identifier: str, identifier_type: str = 'name', normal_boiling_point: float = None)[source]¶
Bases:
object
Joback group contribution properties estimator.
The class recieves either the Joback and Reid model’s [10, 11] groups, name or smiles of a molecule and estimates the its properties.
- Parameters:
identifier (str or rdkit.Chem.rdchem.Mol) – Identifier of a molecule (name, SMILES, groups, or Chem.rdchem.Mol). Example: you can use hexane, CCCCCC, {“-CH3”: 2, “-CH2-”: 4} for name, SMILES and groups respectively.
identifier_type (str, optional) – Use ‘name’ to search a molecule by name, ‘smiles’ to provide the molecule SMILES representation, ‘groups’ to provide Joback groups or ‘mol’ to provide a rdkit.Chem.rdchem.Mol object, by default “name”.
normal_boiling_point (float, optional) – If provided, will be used to estimate critical temperature, acentric factor, and vapor pressure instead of the estimated normal boiling point, by default None.
- subgroups¶
Joback functional groups of the molecule.
- Type:
dict
- experimental_boiling_temperature¶
User provided experimental normal boiling point [K].
- Type:
float
- critical_temperature¶
Joback estimated critical temperature [K].
- Type:
float
- critical_pressure¶
Joback estimated critical pressure [bar].
- Type:
float
- critical_volume¶
Joback estimated critical volume [cm³/mol].
- Type:
float
- normal_boiling_point¶
Joback estimated normal boiling point [K].
- Type:
float
- fusion_temperature¶
Joback estimated fusion temperature [K].
- Type:
float
- h_formation¶
Joback estimated enthalpy of formation ideal gas at 298 K [kJ/mol].
- Type:
float
- g_formation¶
Joback estimated Gibbs energy of formation ideal gas at 298 K [K].
- Type:
float
- heat_capacity_ideal_gas_params¶
Joback estimated Reid’s ideal gas heat capacity equation parameters [J/mol/K].
- Type:
dict
- h_fusion¶
Joback estimated fusion enthalpy [kJ/mol].
- Type:
float
- h_vaporization¶
Joback estimated vaporization enthalpy at the normal boiling point [kJ/mol].
- Type:
float
- sum_na¶
Joback n_A contribution to liquid viscosity [N/s/m²].
- Type:
float
- sum_nb¶
Joback n_B contribution to liquid viscosity [N/s/m²].
- Type:
float
- molecular_weight¶
Molecular weight from Joback’s subgroups [g/mol].
- Type:
float
- vapor_pressure_params¶
Vapor pressure G and k parameters for the Riedel-Plank-Miller [10] equation [bar].
- Type:
dict
- heat_capacity_ideal_gas(temperature: float | ndarray[Any, dtype[_ScalarType_co]]) float | ndarray[Any, dtype[_ScalarType_co]] [source]¶
Calculate the ideal gas heat capacity [J/mol/K].
Uses the Joback estimated Reid’s ideal gas heat capacity equation parameters [J/mol/K].
- Parameters:
temperature (Union[float, NDArray]) – Temperature [K]
- Returns:
Ideal gas heat capacity [J/mol/K].
- Return type:
Union[float, NDArray]
- heat_capacity_liquid(temperature: float | ndarray[Any, dtype[_ScalarType_co]]) float | ndarray[Any, dtype[_ScalarType_co]] [source]¶
Calculate the liquid heat capacity [J/mol/K].
Uses the Rowlinson-Bondi [10] equation with the Joback estimated properties.
- Parameters:
temperature (Union[float, NDArray]) – Temperature [K]
- Returns:
Ideal gas heat capacity [J/mol/K].
- Return type:
Union[float, NDArray]
- viscosity_liquid(temperature: float | ndarray[Any, dtype[_ScalarType_co]]) float | ndarray[Any, dtype[_ScalarType_co]] [source]¶
Calculate the Joback estimated liquid viscosity [N/s/m²].
- Parameters:
temperature (Union[float, NDArray]) – Temperature [K]
- Returns:
Liquid viscosity [N/s/m²].
- Return type:
Union[float, NDArray]
- vapor_pressure(temperature: float | ndarray[Any, dtype[_ScalarType_co]]) float | ndarray[Any, dtype[_ScalarType_co]] [source]¶
Calculate the vapor pressure [bar].
Uses the Riedel-Plank-Miller [10] equation with the Joback estimated properties.
- Parameters:
temperature (Union[float, NDArray]) – Temperature [K]
- Returns:
Vapor pressure [bar]
- Return type:
Union[float, NDArray]
Fragmentation Models¶
FragmentationModel module.
All ugropy models (joback, unifac, psrk) are instances of the FragmentationModule class.
- class FragmentationModel(subgroups: DataFrame, split_detection_smarts: List[str] = [], problematic_structures: DataFrame | None = None, hideouts: DataFrame | None = None)[source]¶
Bases:
object
FragmentationModel class.
All ugropy supported models are an instance of this class. This class can be used by the user to create their own FragmentationModels.
- Parameters:
subgroups (pd.DataFrame) – Model’s subgroups. Index: ‘group’ (groups names). Mandatory columns: ‘detection_smarts’ (SMARTS representations of the group to detect its precense in the molecule), ‘smarts’ (true SMARTS of the group witouht additional atom detections), ‘contribute’ (dictionary as a string with the group contribution), ‘composed’ (y or n if it is or is not a composed structure), ‘molecular_weight’ (molecular weight of the subgroups).
split_detection_smarts (List[str], optional) – List of subgroups that have different SMARTS representations. by default []
problematic_structures (Union[pd.DataFrame, None], optional) – Model’s problematic/ambiguous structures. Index: ‘smarts’ (SMARTS of the problematic structure). Mandatory columns: ‘contribute’ (dictionary as a string with the structure contribution), by default None
hideouts (Union[pd.DataFrame, None], optional) – Hideouts for each group. Index: ‘group’ (Group of the model that can be hiden). Mandatory columns: ‘hideout’ (other subgroups to find the hiden subgroup), by defautl None
- subgroups¶
Model’s subgroups. Index: ‘group’ (groups names). Mandatory columns: ‘detection_smarts’ (SMARTS representations of the group to detect its precense in the molecule), ‘smarts’ (true SMARTS of the group witouht additional atom detections. If a value is missing uses the corresponding detection_smarts), ‘contribute’ (dictionary as a string with the group contribution), ‘composed’ (y or n if it is or is not a composed structure), ‘molecular_weight’ (molecular weight of the subgroups).
- Type:
pd.DataFrame
- split_detection_smarts¶
List of subgroups that have different SMARTS representations.
- Type:
List[str]
- problematic_structures¶
Model’s problematic/ambiguous structures. Index: ‘smarts’ (SMARTS of the problematic structure). Mandatory columns: ‘contribute’ (dictionary as a string with the structure contribution)
- Type:
pd.Dataframe
- hideouts¶
Hideouts for each group. Index: ‘group’ (Group of the model that can be hiden). Mandatory columns: ‘hideout’ (other subgroups to find the hiden subgroup)
- Type:
pd.DataFrame
- detection_mols¶
Dictionary cotaining all the rdkit Mol object from the detection_smarts subgroups column.
- Type:
dict
- fit_mols¶
Dictionary cotaining all the rdkit Mol object from the smarts subgroups column.
- Type:
dict
- contribution_matrix¶
Contribution matrix of the model built from the subgroups contribute.
- Type:
pd.Dataframe
Gibbs Models¶
GibbsModel module.
- class GibbsModel(subgroups: DataFrame, split_detection_smarts: List[str] = [], problematic_structures: DataFrame | None = None, hideouts: DataFrame | None = None, subgroups_info: DataFrame | None = None, main_groups: DataFrame | None = None)[source]¶
Bases:
FragmentationModel
GibbsModel it’s a fragmentation model dedicated to Gibbs excess models.
unifac, psrk, dortmund are instances of this class.
- Parameters:
subgroups (pd.DataFrame) – Model’s subgroups. Index: ‘group’ (groups names). Mandatory columns: ‘detection_smarts’ (SMARTS representations of the group to detect its precense in the molecule), ‘smarts’ (true SMARTS of the group witouht additional atom detections), ‘contribute’ (dictionary as a string with the group contribution), ‘composed’ (y or n if it is or is not a composed structure), ‘molecular_weight’ (molecular weight of the subgroups).
split_detection_smarts (List[str], optional) – List of subgroups that have different SMARTS representations. by default []
problematic_structures (Union[pd.DataFrame, None], optional) – Model’s problematic/ambiguous structures. Index: ‘smarts’ (SMARTS of the problematic structure). Mandatory columns: ‘contribute’ (dictionary as a string with the structure contribution), by default None
hideouts (Union[pd.DataFrame, None], optional) – Hideouts for each group. Index: ‘group’ (Group of the model that can be hiden). Mandatory columns: ‘hideout’ (other subgroups to find the hiden subgroup), by defautl None
subgroups_info (Union[pd.DataFrame, None], optional) – Information of the model’s subgroups (R, Q, subgroup_number, main_group), by default None
main_groups (Union[pd.DataFrame, None], optional) – Main groups information (no., maingroup_name, subgroups), by default None
- subgroups¶
Model’s subgroups. Index: ‘group’ (groups names). Mandatory columns: ‘detection_smarts’ (SMARTS representations of the group to detect its precense in the molecule), ‘smarts’ (true SMARTS of the group witouht additional atom detections. If a value is missing uses the corresponding detection_smarts), ‘contribute’ (dictionary as a string with the group contribution), ‘composed’ (y or n if it is or is not a composed structure), ‘molecular_weight’ (molecular weight of the subgroups).
- Type:
pd.DataFrame
- split_detection_smarts¶
List of subgroups that have different SMARTS representations.
- Type:
List[str]
- problematic_structures¶
Model’s problematic/ambiguous structures. Index: ‘smarts’ (SMARTS of the problematic structure). Mandatory columns: ‘contribute’ (dictionary as a string with the structure contribution)
- Type:
pd.Dataframe
- hideouts¶
Hideouts for each group. Index: ‘group’ (Group of the model that can be hiden). Mandatory columns: ‘hideout’ (other subgroups to find the hiden subgroup)
- Type:
pd.DataFrame
- detection_mols¶
Dictionary cotaining all the rdkit Mol object from the detection_smarts subgroups column.
- Type:
dict
- fit_mols¶
Dictionary cotaining all the rdkit Mol object from the smarts subgroups column.
- Type:
dict
- contribution_matrix¶
Contribution matrix of the model built from the subgroups contribute.
- Type:
pd.Dataframe
- subgroups_info¶
Information of the model’s subgroups (R, Q, subgroup_number, main_group).
- Type:
pd.DataFrame
- main_groups¶
Main groups information (no., maingroup_name, subgroups).
- Type:
pd.DataFrame
Properties Estimators¶
PropertiesEstimator module.
- class PropertiesEstimator(subgroups: DataFrame, split_detection_smarts: List[str] = [], problematic_structures: DataFrame | None = None, hideouts: DataFrame | None = None, properties_contributions: DataFrame | None = None)[source]¶
Bases:
FragmentationModel
Fragmentation model dedicated to properties estimation models.
joback is a instance of this class.
- Parameters:
subgroups (pd.DataFrame) – Model’s subgroups. Index: ‘group’ (groups names). Mandatory columns: ‘detection_smarts’ (SMARTS representations of the group to detect its precense in the molecule), ‘smarts’ (true SMARTS of the group witouht additional atom detections), ‘contribute’ (dictionary as a string with the group contribution), ‘composed’ (y or n if it is or is not a composed structure), ‘molecular_weight’ (molecular weight of the subgroups).
split_detection_smarts (List[str], optional) – List of subgroups that have different SMARTS representations. by default []
problematic_structures (Union[pd.DataFrame, None], optional) – Model’s problematic/ambiguous structures. Index: ‘smarts’ (SMARTS of the problematic structure). Mandatory columns: ‘contribute’ (dictionary as a string with the structure contribution), by default None
hideouts (Union[pd.DataFrame, None], optional) – Hideouts for each group. Index: ‘group’ (Group of the model that can be hiden). Mandatory columns: ‘hideout’ (other subgroups to find the hiden subgroup), by defautl None
properties_contributions (pd.DataFrame, optional) – Group’s properties contributions, by default None
- subgroups¶
Model’s subgroups. Index: ‘group’ (groups names). Mandatory columns: ‘detection_smarts’ (SMARTS representations of the group to detect its precense in the molecule), ‘smarts’ (true SMARTS of the group witouht additional atom detections. If a value is missing uses the corresponding detection_smarts), ‘contribute’ (dictionary as a string with the group contribution), ‘composed’ (y or n if it is or is not a composed structure), ‘molecular_weight’ (molecular weight of the subgroups).
- Type:
pd.DataFrame
- split_detection_smarts¶
List of subgroups that have different SMARTS representations.
- Type:
List[str]
- problematic_structures¶
Model’s problematic/ambiguous structures. Index: ‘smarts’ (SMARTS of the problematic structure). Mandatory columns: ‘contribute’ (dictionary as a string with the structure contribution)
- Type:
pd.Dataframe
- hideouts¶
Hideouts for each group. Index: ‘group’ (Group of the model that can be hiden). Mandatory columns: ‘hideout’ (other subgroups to find the hiden subgroup)
- Type:
pd.DataFrame
- detection_mols¶
Dictionary cotaining all the rdkit Mol object from the detection_smarts subgroups column.
- Type:
dict
- fit_mols¶
Dictionary cotaining all the rdkit Mol object from the smarts subgroups column.
- Type:
dict
- contribution_matrix¶
Contribution matrix of the model built from the subgroups contribute.
- Type:
pd.Dataframe
- properties_contributions¶
Group’s properties contributions.
- Type:
pd.DataFrame
Writers¶
Writers module.
The witers module contains functions to construct the neccesary inputs for other thermodynamics libraries.
Supported:
Clapeyron.jl: https://github.com/ClapeyronThermo/Clapeyron.jl
- to_clapeyron(molecules_names: List[str], unifac_groups: List[dict] = [], psrk_groups: List[dict] = [], joback_objects: List[JobackProperties] = [], path: str = 'database', batch_name: str = '') None [source]¶
Write the .csv input files for Clapeyron.jl.
The provided lists must have the same length. If one of the model lists is left empty, that model will not be writen in the .csv.
- Parameters:
molecules_names (List[str]) – List of names for each chemical to write in the .csv files.
unifac_groups (List[dict], optional) – List of classic liquid-vapor UNIFAC groups, by default [].
psrk_groups (List[dict], optional) – List of Predictive Soave-Redlich-Kwong groups, by default [].
joback_objects (List[JobackProperties], optional) – List of ugropy.properties.JobackProperties objects, by default [].
path (str, optional) – Path to the directory to store de .csv files, by default “./database”.
batch_name (str, optional) – Name of the writing batch. For example, if you name the batch with “batch1”, the output of the UNIFAC groups will be: “batch1_ogUNIFAC_groups.csv”. With the default value will be “ogUNIFAC_groups.csv”, by default “”.
- to_thermo(mol_subgroups: dict, model: GibbsModel) dict [source]¶
Obtain the subgroups dictionary to the Caleb Bell’s Thermo library.
Thermo: https://github.com/CalebBell/thermo
- Parameters:
mol_subgroups (dict) – ugropy subgroups.
model (GibbsModel) – Gibbs excess FragmentationModel (unifac or psrk).
- Returns:
Thermo fragmentation subgroups.
- Return type:
dict
Clapeyron writers¶
Clapeyron writers functions module.
- write_critical(path: Path, batch_name: str, molecules_names: List[str], joback_objects: List[JobackProperties] = []) None [source]¶
Create the DataFrame with the critical properties for Clapeyron.jl.
Uses the Joback to estimate the critical properties of the molecules.
- Parameters:
path (pathlib.Path, optional) – Path to the directory to store de .csv files, by default “./database”.
batch_name (str, optional) – Name of the writing batch. For example, if you name the batch with “batch1”, the output of the UNIFAC groups will be: “batch1_ogUNIFAC_groups.csv”. With the default value will be “ogUNIFAC_groups.csv”, by default “”.
molecules_names (List[str]) – List of names for each chemical to write in the .csv files.
joback_objects (List[Joback], optional) – List of ugropy.properties.JobackProperties objects, by default [].
- Returns:
DataFrame with the molecular weights for Clapeyron.jl
- Return type:
pd.DataFrame
- write_molar_mass(path: Path, batch_name: str, molecules_names: List[str], unifac_groups: List[dict] = [], psrk_groups: List[dict] = [], joback_objects: List[JobackProperties] = []) None [source]¶
Create the DataFrame with the molecular weights for Clapeyron.jl.
- Parameters:
path (pathlib.Path) – Path to the directory to store de .csv files, by default “./database”.
batch_name (str, optional) – Name of the writing batch. For example, if you name the batch with “batch1”, the output of the UNIFAC groups will be: “batch1_ogUNIFAC_groups.csv”. With the default value will be “ogUNIFAC_groups.csv”, by default “”.
molecules_names (List[str]) – List of names for each chemical to write in the .csv files.
unifac_groups (List[dict], optional) – List of classic liquid-vapor UNIFAC groups, by default [].
psrk_groups (List[dict], optional) – List of Predictive Soave-Redlich-Kwong groups, by default [].
joback_objects (List[Joback], optional) – List of ugropy.properties.JobackProperties objects, by default [].
- Returns:
DataFrame with the molecular weights for Clapeyron.jl
- Return type:
pd.DataFrame
- write_psrk(path: Path, batch_name: str, molecules_names: List[str], psrk_groups: List[dict]) None [source]¶
Create the DataFrame with the PSRK groups for Clapeyron.jl.
- Parameters:
path (pathlib.Path) – Path to the directory to store de .csv files, by default “./database”.
batch_name (str, optional) – Name of the writing batch. For example, if you name the batch with “batch1”, the output of the UNIFAC groups will be: “batch1_ogUNIFAC_groups.csv”. With the default value will be “ogUNIFAC_groups.csv”, by default “”.
molecules_names (List[str]) – List of names for each chemical to write in the .csv files.
psrk_groups (List[dict], optional) – List of Predictive Soave-Redlich-Kwong groups.
- Returns:
DataFrame with the LV-UNIFAC groups for Clapeyron.jl
- Return type:
pd.DataFrame
- write_unifac(path: Path, batch_name: str, molecules_names: List[str], unifac_groups: List[dict]) None [source]¶
Create the DataFrame with the classic LV-UNIFAC groups for Clapeyron.jl.
- Parameters:
path (pathlib.Path) – Path to the directory to store de .csv files, by default “./database”.
batch_name (str, optional) – Name of the writing batch. For example, if you name the batch with “batch1”, the output of the UNIFAC groups will be: “batch1_ogUNIFAC_groups.csv”. With the default value will be “ogUNIFAC_groups.csv”, by default “”.
molecules_names (List[str]) – List of names for each chemical to write in the .csv files.
unifac_groups (List[dict], optional) – List of classic liquid-vapor UNIFAC groups.