swiftpol
A Python package for building in silico polymer systems
swiftpol.build
- swiftpol.build.build_polymer(sequence, monomer_list, reaction, terminal='hydrogen', chain_num=1, chainID='A')[source]
Constructs a polymer from a given sequence of monomers.
- Parameters:
sequence (str) – A string representing the sequence of monomers (e.g., ‘ABAB’).
monomer_list (list) – A list of SMILES strings representing the monomers.
reaction (rdkit.Chem.rdChemReactions.ChemicalReaction) – An RDKit reaction object used to link monomers.
terminal (str, optional) – The terminal group to be added to the polymer. Options are ‘hydrogen’, ‘carboxyl’, ‘ester’, or a canonical smiles string to insert as the terminal. Default is ‘hydrogen’.
chain_number (int, optional) – The number of polymer chains to construct. Default is 1. Input used for ensemble build.
- Returns:
The constructed polymer as an RDKit molecule object.
- Return type:
rdkit.Chem.rdchem.Mol
- swiftpol.build.build_linear_copolymer(sequence, monomer_a_smiles, monomer_b_smiles, reaction=<rdkit.Chem.rdChemReactions.ChemicalReaction object>)[source]
Constructs a linear co-polymer from the provided sequence of monomers.
This function takes a sequence of monomers represented as ‘A’ and ‘B’, and the SMILES strings of two monomers. It constructs a co-polymer based on the sequence, using the provided reaction SMARTS for joining the monomers. The function returns the sanitized polymer and the percentage composition of each monomer in the polymer.
- Parameters:
sequence (str) – A string representing the sequence of monomers. ‘A’ represents monomer_a and ‘B’ represents monomer_b.
monomer_a_smiles (str) – The SMILES string of monomer A.
monomer_b_smiles (str) – The SMILES string of monomer B.
reaction (rdkit.Chem.rdChemReactions.ChemicalReaction, optional) – The reaction SMARTS used for joining the monomers. Defaults to ‘[C:1][HO:2].[HO:3][C:4]>>[C:1][O:2][C:4].[O:3]’, representing a condensation polymerisation.
- Returns:
A tuple containing the following elements: - sanitized_polymer (rdkit.Chem.rdchem.Mol): The constructed and sanitized polymer. - percentage_monomer_a (float): The percentage composition of monomer A in the polymer. - percentage_monomer_b (float): The percentage composition of monomer B in the polymer.
- Return type:
tuple
- swiftpol.build.PDI(chains)[source]
Calculates the Polydispersity Index (PDI), number-average molecular weight (Mn), and weight-average molecular weight (Mw) of a list of chains.
- Parameters:
chains (list) – A list of molecular chains. Each chain is represented as an RDKit molecule object.
- Returns:
A tuple containing the following elements: - PDI (float): The Polydispersity Index, which is the ratio of Mw to Mn. - Mn (float): The number-average molecular weight. - Mw (float): The weight-average molecular weight.
- Return type:
tuple
- swiftpol.build.blockiness_gen(sequence, wrt='A')[source]
Calculate the blockiness and average block length of a co-polymer sequence.
This function calculates the blockiness of a co-polymer sequence by counting the occurrences of ‘BB’ and ‘BA’ or ‘AB’ in the sequence. It also calculates the average block length of ‘A’ and ‘B’ monomers in the sequence.
- Parameters:
sequence (str) – A string representing the co-polymer sequence. ‘A’ represents one type of monomer and ‘B’ represents another type.
- Returns:
A tuple containing the following elements: - blockiness (float): The blockiness of the co-polymer sequence. Calculated as the ratio of ‘BB’ to ‘BA’ or ‘AB’. - block_length_A (float): The average block length of ‘A’ in the sequence. - block_length_B (float): The average block length of ‘B’ in the sequence.
- Return type:
tuple
Notes
If the sequence does not contain both ‘A’ and ‘B’, the function returns a string indicating that the molecule is not a co-polymer.
- swiftpol.build.calculate_box_components(chains, monomers, sequence, salt_concentration=<Quantity(0.0, 'mole / liter')>, residual_monomer=0.0)[source]
Calculates the components required to construct a simulation box for a given set of molecular chains.
This function determines the quantity of each molecule type required, considering the salt concentration and residual monomer concentration. It is adapted from the OpenFF Toolkit Packmol wrapper’s solvate_topology function.
- Parameters:
chains (list) – A list of molecular chains to be included in the simulation box.
sequence (str) – A string representing the sequence of the molecular chains. ‘G’ and ‘L’ represent different types of monomers.
salt_concentration (float, optional) – The desired salt concentration in the simulation box. Defaults to 0 M.
residual_monomer (float, optional) – The desired residual monomer concentration in the simulation box. Defaults to 0.00%.
- Returns:
tuple – A tuple containing the following elements: - molecules (list): A list of molecules to be included in the simulation box. - number_of_copies (list): A list indicating the quantity of each molecule to be included in the simulation box. - topology (openff.toolkit.topology.Topology): The topology of the simulation box. - box_vectors (numpy.ndarray): The vectors defining the dimensions of the simulation box.
Notes
——
This function is adapted from the OpenFF Toolkit Packmol wrapper’s solvate_topology function.
- swiftpol.build.introduce_stereoisomers(stereo_monomer, instance, seq)[source]
Introduce stereoisomers into a polymer sequence by replacing a specified percentage of ‘A’ monomers with ‘S’.
This function replaces a specified percentage of ‘A’ monomers with ‘S’ in the given sequence. The replacements are made in pairs of ‘A’ monomers to ensure stereoisomerism.
- Parameters:
stereo_monomer (str) – The monomer to be replaced with its stereoisomer (e.g., ‘A’).
instance (float) – The fraction of ‘stereo_monomer’ monomers to be replaced with stereoisomer (e.g., 0.5 for 50%).
seq (str) – The original polymer sequence.
- Returns:
The modified sequence with stereoisomers introduced.
- Return type:
str
swiftpol.build.polymer_system
- class swiftpol.build.polymer_system(monomer_list, reaction, length_target, num_chains, stereoisomerism_input=None, terminals='hydrogen', perc_A_target=100, blockiness_target=[1.0, 'A'], copolymer=False, diblock=False, acceptance=10)[source]
Bases:
object- string = <module 'string' from '/usr/share/miniconda/envs/swiftpol/lib/python3.12/string.py'>
- __init__(monomer_list, reaction, length_target, num_chains, stereoisomerism_input=None, terminals='hydrogen', perc_A_target=100, blockiness_target=[1.0, 'A'], copolymer=False, diblock=False, acceptance=10)[source]
Initialize the polymer system and build the polymer chains.
- Parameters:
monomer_list (list) (List of monomers to be used in the polymerization.)
reaction (str) (The type of reaction to be used for polymerization.)
length_target (float) (The target length of the polymer chains.)
num_chains (int) (The number of polymer chains to be generated.)
stereoisomerism_input (tuple, optional) (A tuple containing the monomer, instance fraction (e.g. 0.5 for 50% stereoisomer), and SMILES string of the stereoisomer to be introduced. Default is None.)
terminals (str, optional) (The type of terminal groups to be used. Default is ‘hydrogen’, adds a hydrogen atom.)
perc_A_target (float, optional) (The target percentage of monomer A in the copolymer. Default is 100.)
blockiness_target (float, optional) (The target blockiness of the copolymer, and indication of method calculation. Default is 1.0, with reference to ‘A’ monomer linkages.)
copolymer (bool, optional) (Flag to indicate if the system is a copolymer. Default is False.)
diblock (bool, optional) (Flag to indicate if the system is a diblock copolymer. Default is False.)
acceptance (float, optional) (The percentage deviation of blockiness and A percentage from target values. Default is 10%.)
- Variables:
(float) (max_length)
(str) (terminals)
(float)
(float)
(list) (blockiness_list)
(list)
(list)
(list)
(float)
(float)
(list)
(float)
(float)
(float)
(float)
(float)
(int) (num_chains)
(float)
(float)
(float)
- generate_conformers(rough=False)[source]
Generate conformers for each polymer chain in the system.
This method uses the OpenFF toolkit OpenEye Wrapper to generate conformers for each polymer chain in the system. It first checks if the OpenEye toolkit is licensed and available. If it is, it uses the OpenEyeToolkitWrapper to generate conformers. Otherwise, it falls back to using the RDKitToolkitWrapper. Each chain is processed to generate a single conformer, and unique atom names are assigned to each chain.
- Parameters:
rough (bool, optional) – If True, generates rough initial conformers using RDKit’s EmbedMolecule function with the following settings: - useExpTorsionAnglePrefs = False - useBasicKnowledge = False - useRandomCoords = True - maxIterations = 500 - numThreads = 1 Default is False.
- Raises:
ImportError – If neither RDKit nor OpenEye toolkits are available.
- charge_system(charge_scheme)[source]
Assign partial charges to each polymer chain in the system.
This method uses one of AM1-BCC, Espaloma, or OpenFF NAGL to assign partial charges to each polymer chain in the system. It iterates over each chain in the self.chains list and assigns partial charges to the chain.
- Parameters:
charge_scheme (str) – The charge assignment scheme to use. Options are ‘AM1_BCC’, ‘espaloma’, or ‘NAGL’.
- Raises:
ImportError – If the selected toolkit is not available.
- export_to_csv(filename)[source]
Export all the instances in polymer_system.__init__() into a pandas DataFrame and save it as a CSV file.
- Parameters:
filename (str) – The name of the CSV file to save the data.
- pack_solvated_system(salt_concentration=<Quantity(0.0, 'mole / liter')>, residual_monomer=0.0)[source]
Pack a solvated system using Packmol functions, and the OpenFF Packmol wrapper.
This method uses Packmol to build a solvated system by packing molecules into a simulation box. It considers the salt concentration and residual monomer concentration to determine the quantity of each molecule type required.
- Parameters:
salt_concentration (openff.units.Quantity, optional) – The desired salt concentration in the simulation box. Default is 0.0 M.
residual_monomer (float, optional) – The desired residual monomer concentration in the simulation box. Default is 0.00.
- Returns:
An Interchange object representing the packed solvated system.
- Return type:
openff.toolkit.topology.topology.Topology
Notes
This function uses the OpenFF Interchange Packmol wrapper to pack the molecules into the simulation box. It removes any molecules with zero copies before packing, to avoid packmol errors. Assigns % residual monomer value to ensemble under self.residual_monomer_actual
- generate_polyply_files(residual_monomer=0.0, residual_oligomer=0.0)[source]
Generate input files for Polyply from the system.
This method generates the input files required for Polyply (https://github.com/marrink-lab/polyply_1.0/) from the system.
- Parameters:
residual_monomer (float, optional) – The desired residual monomer concentration in the simulation box. Default is 0.00.
residual_oligomer (float, optional) – The desired residual oligomer concentration in the simulation box. Default is 0.00.
- Returns:
A tuple containing the paths to the .gro, .top, and .pdb files generated for Polyply and GROMACS insert-molecules input.
- Return type:
tuple
Notes
This function uses OpenFF Interchange to generate the input files for Polyply.
- Raises:
UserWarning – If partial charges are not assigned to the system, processing large systems may raise errors from OpenFF-Interchange.
- calculate_residuals(residual_monomer=0, residual_oligomer=0, return_rdkit=False)[source]
Generate residual monomer and oligomer molecules, and molecule counts.
- Parameters:
residual_monomer (float) – The desired residual monomer concentration in the simulation box. Default is 0.00.
residual_oligomer (float) – The desired residual oligomer concentration in the simulation box. Default is 0.00.
return_rdkit (bool, optional) – If True, return RDKit molecule objects instead of OpenFF Molecule objects. Default is False.
- Returns:
tuple
A tuple containing the following elements –
A list of OpenFF Molecule objects representing the residual monomer and oligomer molecules.
A list of molecule counts, corresponding with the molecule list
Actual value for % residual monomer
Actual value for % residual oligomer
Notes
This function works independently of, and is unrelated to, the build.calculate_box_components function, which will be deprecated in future releases.
- Raises:
UserWarning - If the residual monomer concentration is close to or lower than the weight of a single monomer. –
swiftpol.build.polymer_system_from_PDI
- class swiftpol.build.polymer_system_from_PDI(monomer_list, reaction, length_target, num_chains, PDI_target, stereoisomerism_input=None, terminals='hydrogen', perc_A_target=100, blockiness_target=[1.0, 'A'], copolymer=False, diblock=False, acceptance=10, verbose=None)[source]
Bases:
object- string = <module 'string' from '/usr/share/miniconda/envs/swiftpol/lib/python3.12/string.py'>
- __init__(monomer_list, reaction, length_target, num_chains, PDI_target, stereoisomerism_input=None, terminals='hydrogen', perc_A_target=100, blockiness_target=[1.0, 'A'], copolymer=False, diblock=False, acceptance=10, verbose=None)[source]
Initialize the polymer system and build the polymer chains.
This class initializes a polymer system and generates polymer chains based on the specified parameters. It supports homopolymers, copolymers, and diblock copolymers, with options for stereoisomerism, terminal groups, and blockiness control. The class also allows for fine-tuning of the polydispersity index (PDI) and other polymer properties.
- Parameters:
monomer_list (list) – List of monomers to be used in the polymerization.
reaction (str) – The type of reaction to be used for polymerization.
length_target (float) – The target length of the polymer chains.
num_chains (int) – The number of polymer chains to be generated.
PDI_target (float) – The target polydispersity index (PDI) of the polymer chains. A PDI value of up to 3.0 is supported. Smaller systems may not achieve high PDI values.
stereoisomerism_input (tuple, optional) –
- A tuple containing:
The monomer (str).
The instance fraction (float, e.g., 0.5 for 50% stereoisomer).
The SMILES string (str) of the stereoisomer to be introduced.
Default is None.
terminals (str, optional) – The type of terminal groups to be used. Default is ‘hydrogen’, which adds a hydrogen atom.
perc_A_target (float, optional) – The target percentage of monomer A in the copolymer. Default is 100.
blockiness_target (float, optional) – The target blockiness of the copolymer, indicating the degree of clustering of monomers A and B. Default is 1.0, with reference to ‘A’ monomer linkages.
copolymer (bool, optional) – Flag to indicate if the system is a copolymer. Default is False.
diblock (bool, optional) – Flag to indicate if the system is a diblock copolymer. Default is False.
acceptance (float, optional) – The percentage deviation allowed for blockiness and A percentage from target values. Default is 10%.
verbose (bool, optional) – If True, enables logging messages for debugging and progress tracking. Default is False.
- Variables:
length_target (float) – The target length of the polymer chains.
terminals (str) – The type of terminal groups used.
blockiness_target (float) – The target blockiness of the copolymer.
A_target (float) – The target percentage of monomer A in the copolymer.
chains (list) – List of polymer chains as OpenFF Molecule objects.
chain_rdkit (list) – List of polymer chains as RDKit molecule objects.
lengths (list) – List of lengths of the polymer chains.
perc_A_actual (list) – List of actual percentages of monomer A in the polymer chains.
B_block_length (float) – The average block length of monomer B in the copolymer.
A_block_length (float) – The average block length of monomer A in the copolymer.
blockiness_list (list) – List of blockiness values for the polymer chains.
mean_blockiness (float) – The mean blockiness of the polymer chains.
mol_weight_average (float) – The average molecular weight of the polymer chains.
PDI (float) – The polydispersity index of the polymer chains.
Mn (float) – The number-average molecular weight of the polymer chains.
Mw (float) – The weight-average molecular weight of the polymer chains.
num_chains (int) – The number of polymer chains generated.
length_average (float) – The average length of the polymer chains.
min_length (float) – The minimum length of the polymer chains.
max_length (float) – The maximum length of the polymer chains.
Notes
This class is under development and may change significantly in future versions of SwiftPol.
For debugging and progress tracking, set verbose=True.
- generate_conformers(rough=False)[source]
Generate conformers for each polymer chain in the system.
This method uses the OpenFF toolkit OpenEye Wrapper to generate conformers for each polymer chain in the system. It first checks if the OpenEye toolkit is licensed and available. If it is, it uses the OpenEyeToolkitWrapper to generate conformers. Otherwise, it falls back to using the RDKitToolkitWrapper. Each chain is processed to generate a single conformer, and unique atom names are assigned to each chain.
- Parameters:
rough (bool, optional) – If True, generates rough initial conformers using RDKit’s EmbedMolecule function with the following settings: - useExpTorsionAnglePrefs = False - useBasicKnowledge = False - useRandomCoords = True - maxIterations = 500 - numThreads = 1 Default is False.
- Raises:
ImportError – If neither RDKit nor OpenEye toolkits are available.
- charge_system(charge_scheme)[source]
Assign partial charges to each polymer chain in the system.
This method uses one of AM1-BCC, Espaloma, or OpenFF NAGL to assign partial charges to each polymer chain in the system. It iterates over each chain in the self.chains list and assigns partial charges to the chain.
- Parameters:
charge_scheme (str) – The charge assignment scheme to use. Options are ‘AM1_BCC’, ‘espaloma’, or ‘NAGL’.
- Raises:
ImportError – If the selected toolkit is not available.
- export_to_csv(filename)[source]
Export all the instances in polymer_system.__init__() into a pandas DataFrame and save it as a CSV file.
- Parameters:
filename (str) – The name of the CSV file to save the data.
- pack_solvated_system(salt_concentration=0.0, residual_monomer=0.0)[source]
Pack a solvated system using Packmol functions, and the OpenFF Packmol wrapper.
This method uses Packmol to build a solvated system by packing molecules into a simulation box with water. It considers the salt concentration and residual monomer concentration to determine the quantity of each molecule type required.
- Parameters:
salt_concentration (openff.units.Quantity, optional) – The desired salt concentration in the simulation box. Default is 0.0 M.
residual_monomer (float, optional) – The desired residual monomer concentration in the simulation box. Default is 0.00.
- Returns:
An Interchange object representing the packed solvated system.
- Return type:
openff.toolkit.topology.topology.Topology
Notes
This function uses the OpenFF Interchange Packmol wrapper to pack the molecules into the simulation box. It removes any molecules with zero copies before packing, to avoid packmol errors. Assigns % residual monomer value to ensemble under self.residual_monomer_actual
- generate_polyply_files(residual_monomer=0.0, residual_oligomer=0.0)[source]
Generate input files for Polyply from the system.
This method generates the input files required for Polyply (https://github.com/marrink-lab/polyply_1.0/) from the system.
- Parameters:
residual_monomer (float, optional) – The desired residual monomer concentration in the simulation box. Default is 0.00.
residual_oligomer (float, optional) – The desired residual oligomer concentration in the simulation box. Default is 0.00.
- Returns:
A tuple containing the paths to the .gro, .top, and .pdb files generated for Polyply and GROMACS insert-molecules input.
- Return type:
tuple
Notes
This function uses OpenFF Interchange to generate the input files for Polyply.
- Raises:
UserWarning – If partial charges are not assigned to the system, processing large systems may raise errors from OpenFF-Interchange.
- calculate_residuals(residual_monomer=0, residual_oligomer=0, return_rdkit=False)[source]
Generate residual monomer and oligomer molecules, and molecule counts.
- Parameters:
residual_monomer (float) – The desired residual monomer concentration in the simulation box. Default is 0.00.
residual_oligomer (float) – The desired residual oligomer concentration in the simulation box. Default is 0.00.
return_rdkit (bool, optional) – If True, return RDKit molecule objects instead of OpenFF Molecule objects. Default is False.
- Returns:
tuple
A tuple containing the following elements –
A list of OpenFF Molecule objects representing the residual monomer and oligomer molecules.
A list of molecule counts, corresponding with the molecule list
Actual value for % residual monomer
Actual value for % residual oligomer
Notes
This function works independently of, and is unrelated to, the build.calculate_box_components function, which will be deprecated in future releases.
- Raises:
UserWarning - If the residual monomer concentration is close to or lower than the weight of a single monomer. –
swiftpol.parameterize
- swiftpol.parameterize.charge_polymer(polymer, charge_scheme)[source]
Calculate and return the partial charges of a polymer chain based on the specified charge scheme.
- Parameters:
polymer (rdkit.Chem.rdchem.Mol) – A polymer chain for which the charges are to be calculated.
charge_scheme (str) – A string that specifies the charge scheme to be used. It can be either ‘AM1_BCC’, ‘espaloma’, or ‘NAGL’.
- Returns:
The partial charges of the polymer chain according to the specified charge scheme.
- Return type:
numpy.ndarray
- Raises:
AttributeError – If the charge_scheme input is not ‘AM1_BCC’, ‘NAGL’, or ‘espaloma’.
- swiftpol.parameterize.charge_openff_polymer(openff_chain, charge_scheme, overwrite=True)[source]
Assign partial charges to a polymer chain based on the specified charge scheme.
This function assigns partial charges to a polymer chain using one of the following charge schemes: ‘AM1_BCC’, ‘NAGL’, or ‘espaloma’. It can optionally overwrite existing charges.
- Parameters:
openff_chain (openff.toolkit.topology.Molecule) – The original polymer chain to which partial charges will be assigned.
charge_scheme (str) – The charge scheme to use for assigning partial charges. Options are ‘AM1_BCC’, ‘NAGL’, or ‘espaloma’.
overwrite (bool, optional) – Whether to overwrite existing partial charges. Default is True.
- Returns:
The partial charges assigned to the polymer chain.
- Return type:
openff.toolkit.topology.molecule.PartialChargeArray
- Raises:
ModuleNotFoundError – If the installed version of openff-toolkit is below what is required to use NAGL.
ImportError – If the required package for the specified charge scheme is not installed.
AttributeError – If the charge_scheme input is not ‘AM1_BCC’, ‘NAGL’, or ‘espaloma’.
swiftpol.crosslink
- swiftpol.crosslink.replace_halogens_with_hydrogens(mol)[source]
Replaces all halogen atoms (F, Cl, Br, I) in an RDKit molecule with hydrogen atoms.
This function iterates through the atoms in the given molecule, identifies halogen atoms (fluorine, chlorine, bromine, iodine etc.), removes them, and replaces them with explicit hydrogens.
- Parameters:
mol (rdkit.Chem.Mol) – The input RDKit molecule from which halogens will be replaced with hydrogens.
- Returns:
A new RDKit molecule with all halogens replaced by explicit hydrogens.
- Return type:
rdkit.Chem.Mol
- swiftpol.crosslink.validate_linear_reaction(starting_polymer, linear_activate, linear_activate_reactants, linear_react)[source]
Validates the linear activation and polymerization reactions for compatibility with the starting polymer and reactants.
Note
This function is not recommended for standalone use.
- Parameters:
starting_polymer (rdkit.Chem.Mol) – The initial polymer molecule to validate against the reactions.
linear_activate (rdkit.Chem.rdChemReactions.ChemicalReaction) – The reaction template for linear activation.
linear_activate_reactants (list of rdkit.Chem.Mol) – The reactants required for the linear activation reaction.
linear_react (rdkit.Chem.rdChemReactions.ChemicalReaction) – The reaction template for linear polymerization.
- Raises:
ValueError – If the reactions are invalid or incompatible with the starting polymer.
- swiftpol.crosslink.validate_branched_reaction(starting_polymer, branched_activate, branched_activate_reactants, branched_react, linear_activate, linear_activate_reactants, linear_react)[source]
Validates the branched and linear activation and polymerization reactions for compatibility with the starting polymer and reactants.
Note
This function is not recommended for standalone use.
- Parameters:
starting_polymer (rdkit.Chem.Mol) – The initial polymer molecule to validate against the reactions.
branched_activate (rdkit.Chem.rdChemReactions.ChemicalReaction) – The reaction template for branched activation.
branched_activate_reactants (list of rdkit.Chem.Mol) – The reactants required for the branched activation reaction.
branched_react (rdkit.Chem.rdChemReactions.ChemicalReaction) – The reaction template for branched polymerization.
linear_activate (rdkit.Chem.rdChemReactions.ChemicalReaction) – The reaction template for linear activation.
linear_activate_reactants (list of rdkit.Chem.Mol) – The reactants required for the linear activation reaction.
linear_react (rdkit.Chem.rdChemReactions.ChemicalReaction) – The reaction template for linear polymerization.
- Raises:
ValueError – If the reactions are invalid or incompatible with the starting polymer.
- swiftpol.crosslink.iterative_chainID_update(polymer, last_polymer_added)[source]
Updates the chain IDs of atoms in a polymer molecule iteratively based on the last polymer added.
This function modifies the polymer molecule by updating the chain IDs of its atoms. The chain IDs are incremented sequentially, starting from the chain ID of the last_polymer_added. The function ensures that the chain IDs wrap around using the English alphabet (A-Z).
- Parameters:
polymer (rdkit.Chem.Mol) – The polymer molecule whose chain IDs will be updated.
last_polymer_added (rdkit.Chem.Mol) – The last polymer molecule added to the network. The chain IDs of this molecule are used as the starting point for updating the chain IDs in the polymer.
- Returns:
polymer – The updated polymer molecule with modified chain IDs.
- Return type:
rdkit.Chem.Mol
Notes
Not recommended for standalone use.
The function uses the AtomPDBResidueInfo object to access and modify chain IDs during the network building process.
Chain IDs are updated using the English alphabet (A-Z) and wrap around if they exceed ‘Z’.
Atoms without MonomerInfo are skipped, and a warning is printed for each skipped atom. For polymers build with SwiftPol, this should not occur.
- swiftpol.crosslink.build_branched_polymer(starting_polymer, reaction_templates, num_iterations=None, target_mw=None, probability_of_branched_addition=0.5, probability_of_linear_addition=0.5)[source]
Builds a branched polymer.
- Parameters:
starting_polymer (rdkit.Chem.Mol) – The initial polymer molecule to which reactions will be applied.
reaction_templates (dict) – A dictionary containing reaction templates for both linear and branched additions. Keys should include:
“branched_activate”: Activation reaction for branching.
“branched_react”: Reaction for branching.
“linear_activate”: Activation reaction for linear addition.
“linear_react”: Reaction for linear addition.
num_iterations (int, optional) – The number of iterations to perform. If specified, the process will stop after this many successful reactions, regardless of the target molecular weight.
target_mw (float, optional) – The target molecular weight for the polymer. If specified, the process will stop once the molecular weight of the polymer reaches or exceeds this value.
probability_of_branched_addition (float, optional, default=0.5) – The probability of performing a branched addition at each iteration.
probability_of_linear_addition (float, optional, default=0.5) – The probability of performing a linear addition at each iteration.
- Returns:
polymer_network – The resulting polymer molecule after the specified number of iterations or upon reaching the target molecular weight. Random stereochemsitry assigned to each stereocenter.
- Return type:
rdkit.Chem.Mol
Notes
The function alternates between linear and branched additions based on the specified probabilities.
- If both num_iterations and target_mw are specified, the process will stop as soon as either
condition is met.
- The reaction templates must be compatible with the starting polymer and follow RDKit’s reaction
SMARTS format.
- Raises:
ValueError – If the reaction templates are invalid or incompatible with the starting polymer.
AssertionError – If the sum of probabilities does not equal 1.
ValueError – If neither num_iterations nor target_mw is specified.
- swiftpol.crosslink.crosslink_polymer(mol, percentage_to_crosslink=50)[source]
Crosslinks a specified percentage of terminal chlorine atoms in a polymer and replaces the remaining halogens with hydrogens.
- Parameters:
mol (rdkit.Chem.Mol) – The input molecule representing the polymer. Must contain terminal chlorine atoms at the point of crosslinking. The Chlorine dependency will be removed in future releases.
percentage_to_crosslink (float, optional, default=50) – The percentage (0-100) of terminal chlorine atoms to crosslink. The remaining halogens will be replaced with hydrogens.
- Returns:
mol (rdkit.Chem.Mol) – The modified molecule with the specified percentage of crosslinks and the remaining halogens replaced by hydrogens.
Notes
——
- The function identifies terminal chlorine atoms bonded to carbon atoms.
- Crosslinking is performed by creating single bonds between carbon atoms – of selected chlorine sites.
- Remaining halogens are replaced with hydrogens using the – replace_halogens_with_hydrogens function.
The percentage_to_crosslink parameter determines the proportion of – chlorine atoms to crosslink.
- Raises:
ValueError – If percentage_to_crosslink is not between 0 and 100.