Molecule¶
-
class
scine_molassembler.
Molecule
¶ Models a molecule as a
Graph
and aStereopermutatorList
.-
__init__
(*args, **kwargs)¶ Overloaded function.
__init__(self: scine_molassembler.Molecule) -> None
Initialize a hydrogen molecule
>>> h2 = Molecule() >>> h2.graph.V 2 >>> h2.graph.E 1
__init__(self: scine_molassembler.Molecule, arg0: Scine::Utils::ElementType) -> None
Initialize a single-atom molecule.
This is a bit of a paradox, yes, and it might have been preferable for the concept of a molecule to contain at least two bonded atoms, but unfortunately single atoms occur everywhere and enforcing the concept would complicate many interfaces.
>>> import scine_utilities as utils >>> f = Molecule(utils.ElementType.F) >>> f.graph.V 1 >>> f.graph.E 0
__init__(self: scine_molassembler.Molecule, first_element: Scine::Utils::ElementType, second_element: Scine::Utils::ElementType, bond_type: scine_molassembler.BondType = <BondType.Single: 0>) -> None
Initialize a molecule from two element types and a mutual
BondType
>>> # Make H-F >>> import scine_utilities as utils >>> hf = Molecule(utils.ElementType.H, utils.ElementType.F) >>> hf.graph.V == 2 True
__init__(self: scine_molassembler.Molecule, graph: scine_molassembler.Graph) -> None
Initialize a molecule from connectivity alone, inferring shapes and stereopermutators from the graph.
>>> # Rebuild a molecule with an assigned stereopermutator from just the graph >>> a = io.experimental.from_smiles("[C@](F)(Cl)(C)[H]") >>> a.stereopermutators.has_unassigned_permutators() False >>> b = Molecule(a.graph) >>> b.stereopermutators.has_unassigned_permutators() True
-
add_atom
(self: scine_molassembler.Molecule, element: Scine::Utils::ElementType, adjacent_to: int, bond_type: scine_molassembler.BondType = <BondType.Single: 0>) → int¶ Add an atom to the molecule, attaching it to an existing atom by a specified bond type.
- Parameters
element – Element type of the new atom
adjacent_to – Atom to which the new atom is added
bond_type –
BondType
with which the new atom is attached
>>> # Let's make linear H3 >>> import scine_utilities as utils >>> mol = Molecule() # Default constructor makes H2 >>> _ = mol.add_atom(utils.ElementType.H, 0) # Make linear H3
-
add_bond
(self: scine_molassembler.Molecule, first_atom: int, second_atom: int, bond_type: scine_molassembler.BondType = <BondType.Single: 0>) → scine_molassembler.BondIndex¶ Adds a bond between two existing atoms.
- Parameters
first_atom – First atom to bond
second_atom – Second atom to bond
bond_type –
BondType
with which to bond the atoms
>>> # Let's make triangular H3 >>> import scine_utilities as utils >>> mol = Molecule() # Default constructor makes H2 >>> _ = mol.add_atom(utils.ElementType.H, 0) # Make linear H3 >>> _ = mol.add_bond(1, 2) # Make triangular H3
-
add_permutator
(self: scine_molassembler.Molecule, bond: scine_molassembler.BondIndex, alignment: scine_molassembler.BondStereopermutator.Alignment = scine_molassembler.BondStereopermutator.Alignment.Eclipsed) → scine_molassembler.BondStereopermutator¶ Add a BondStereopermutator to the molecule
Note
You can’t add AtomStereopermutators to the molecule manually. These are automatically present on non-terminal atoms.
- Parameters
bond – Bond to place the permutator at
alignment – Alignment with which to construct the permutator
- Returns
A reference to the added stereopermutator
- Raises
RuntimeError – If there is already a permutator present at this bond
-
static
apply_canonicalization_map
(canonicalization_index_map: List[int], atom_collection: Scine::Utils::AtomCollection) → Scine::Utils::AtomCollection¶ Reorders an atom collection according to an index mapping from canonicalization.
- Parameters
canonicalization_index_map – Index mapping saved from previous canonicalization
atom_collection – Atom collection to reorder
- Returns
Reordered atom collection
-
assign_stereopermutator
(*args, **kwargs)¶ Overloaded function.
assign_stereopermutator(self: scine_molassembler.Molecule, atom: int, assignment_option: Optional[int]) -> None
Sets the atom stereopermutator assignment at a particular atom
- param atom
Atom index of the
AtomStereopermutator
to set- param assignment_option
An assignment integer if the stereopermutator is to be assigned or
None
if the stereopermutator is to be dis-assigned.
>>> # Assign an unspecified asymmetric carbon atom and then dis-assign it >>> mol = io.experimental.from_smiles("F[CH1](Br)C") >>> asymmetric_carbon_index = 1 >>> mol.assign_stereopermutator(asymmetric_carbon_index, 0) >>> mol.stereopermutators.option(asymmetric_carbon_index).assigned 0 >>> mol.assign_stereopermutator(asymmetric_carbon_index, None) >>> mol.stereopermutators.option(asymmetric_carbon_index).assigned is None True
assign_stereopermutator(self: scine_molassembler.Molecule, bond_index: scine_molassembler.BondIndex, assignment_option: Optional[int]) -> None
Sets the bond stereopermutator assignment at a particular bond
- param bond_index
BondIndex
of theBondStereopermutator
to set- param assignment_option
An assignment integer if the stereopermutator is to be assigned or
None
if the stereopermutator is to be dis-assigned.
>>> # Dis-assign an assigned bond stereopermutator >>> ethene = io.experimental.from_smiles("C/C=C\C") >>> double_bond_index = BondIndex(1, 2) >>> assert ethene.graph.bond_type(double_bond_index) == BondType.Double >>> ethene.stereopermutators.option(double_bond_index).assigned is not None True >>> ethene.assign_stereopermutator(double_bond_index, None) >>> ethene.stereopermutators.option(double_bond_index).assigned is not None False
-
assign_stereopermutator_randomly
(*args, **kwargs)¶ Overloaded function.
assign_stereopermutator_randomly(self: scine_molassembler.Molecule, atom: int) -> None
Assigns an
AtomStereopermutator
at random (assignments are weighted by relative statistical occurence).- param atom
Atom index of the stereopermutator to assign randomly.
Note
This function advances
molassembler
’s global PRNG state.>>> # Assign an unspecified chiral center >>> mol = io.experimental.from_smiles("S[As](F)(Cl)(Br)(N)[H]") >>> as_index = 1 >>> mol.stereopermutators.option(as_index).assigned is None True >>> mol.assign_stereopermutator_randomly(1) >>> mol.stereopermutators.option(as_index).assigned is None False
assign_stereopermutator_randomly(self: scine_molassembler.Molecule, bond_index: scine_molassembler.BondIndex) -> None
Assigns a
BondStereopermutator
at random.- param bond_index
BondIndex
of the stereopermutator to assign randomly.
Note
This function advances
molassembler
’s global PRNG state.>>> # Assign an unspecified double bond randomly >>> mol = io.experimental.from_smiles("CC=CC") >>> double_bond_index = BondIndex(1, 2) >>> assert mol.graph.bond_type(double_bond_index) == BondType.Double >>> mol.stereopermutators.option(double_bond_index).assigned is None True >>> mol.assign_stereopermutator_randomly(double_bond_index) >>> mol.stereopermutators.option(double_bond_index).assigned is None False
-
canonical_compare
(self: scine_molassembler.Molecule, other: scine_molassembler.Molecule, components_bitmask: scine_molassembler.AtomEnvironmentComponents = <AtomEnvironmentComponents.All: 15>) → bool¶ Modular comparison of this Molecule with another, assuming that both are in some (possibly partial) canonical form.
For comparisons of fully canonical molecule pairs, regular equality comparison will just call this function with all environment components considered instead of performing a full isomorphism.
This function is similar to modular_isomorphism, but faster, since if both molecules are in a canonical form, comparison does not require an isomorphism, but merely a same-graph test over the components used.
- Parameters
other – The other (canonical) molecule to compare against
components_bitmask – The components of an atom’s environment to include in the comparison. You should use the same bitmask as when canonicalizing the molecules you are comparing here. It may be possible to use a bitmask with fewer components, but certainly not one with more.
>>> # Bring two molecules into a partial canonical form and compare them >>> a = io.experimental.from_smiles("OCC") >>> b = io.experimental.from_smiles("SCC") >>> a == b False >>> # A and B are identical when considered purely by their graph >>> part = AtomEnvironmentComponents.Connectivity >>> _ = a.canonicalize(part) >>> _ = b.canonicalize(part) >>> a.canonical_compare(b, part) True >>> a == b # Partial canonicalization does not change the meaning of strict equality False >>> # Another pair that is identical save for a stereopermutation >>> c = io.experimental.from_smiles("N[C@](Br)(O)C") >>> d = io.experimental.from_smiles("N[C@@](Br)(O)C") >>> c == d # Strict equality includes stereopermutation False >>> part = AtomEnvironmentComponents.ElementsBondsAndShapes >>> _ = c.canonicalize(part) >>> _ = d.canonicalize(part) >>> c.canonical_compare(d, part) # Limited comparison yields equality True
-
property
canonical_components
¶ Yields the components of the molecule that were used in a previous canonicalization. Can be
None
if the molecule was never canonicalized.- Return type
AtomEnvironmentComponents
orNone
>>> # Canonicalize something and retrieve its canonical components >>> mol = io.experimental.from_smiles("C12(CCC1)COCC2") >>> mol.canonical_components is None True >>> _ = mol.canonicalize() >>> mol.canonical_components == AtomEnvironmentComponents.All True
-
canonicalize
(self: scine_molassembler.Molecule, components_bitmask: scine_molassembler.AtomEnvironmentComponents = <AtomEnvironmentComponents.All: 15>) → List[int]¶ Transform the molecule to a canonical form.
- Warning
Invalidates all external atom and bond indices.
Molecule instances can be canonicalized. Graph canonicalization is an algorithm that reduces all isomorphic forms of an input graph into a canonical form. After canonicalization, isomorphism tests are reduced to mere identity tests.
The canonicalization itself, however, is computationally at least as expensive as an isomorphism itself. Therefore, no expense is saved if an isomorphism test is to be computed only once for two molecules by canonizing both. Only if a molecule instance is to be a repeated candidate for isomorphism tests is there value in canonizing it.
This library takes the approach of adding a tag to molecules that identifies which components of the graph and stereocenters have been used in the generation of the canonical form. This tag is voided with the use of any non-const member function. Pay close attention to the documentation of comparison member functions and operators to ensure that you are making good use of the provided shortcuts.
Note that canonicalization information is only retained across IO boundaries using the JSON serialization variations.
- Parameters
components_bitmask – The components of the molecular graph to include in the canonicalization procedure.
- Returns
Flat index mapping/permutation from old indices to new
>>> # Create two different representations of the same molecule >>> a = io.experimental.from_smiles("N[C@](Br)(O)C") >>> b = io.experimental.from_smiles("Br[C@](O)(N)C") >>> # a and be represent the same molecule, but have different vertex order >>> a == b # Equality operators perform an isomorphism for non-canonical pairs True >>> amap = a.canonicalize() >>> bmap = b.canonicalize() >>> amap == bmap # This shows the vertex order was different False >>> a == b # Equality operators perform a same-graph test for canonical pairs (faster) True
-
dump_graphviz
(self: scine_molassembler.Molecule) → str¶ Returns a graphviz string representation of the molecule
-
hash
(self: scine_molassembler.Molecule) → int¶ Calculates a convoluted hash of a molecule. The molecule must be at least partially canonical. Hashes between molecules of different canonicity are not comparable.
>>> # Show that hash values differ at various levels of canonicity >>> from copy import copy >>> spiro = io.experimental.from_smiles("C12(CCC1)CCC2") >>> # We make two variants of the molecule that have different canonicalization states >>> # to demonstrate that their hashes are unequal. We discard the mappings >>> # we get from canonicalize() >>> partially_canonical = copy(spiro) >>> _ = partially_canonical.canonicalize(AtomEnvironmentComponents.ElementsAndBonds) >>> fully_canonical = copy(spiro) >>> _ = fully_canonical.canonicalize() >>> partially_canonical == fully_canonical True >>> partially_canonical.hash() == fully_canonical.hash() False
-
modular_isomorphism
(self: scine_molassembler.Molecule, other: scine_molassembler.Molecule, components_bitmask: scine_molassembler.AtomEnvironmentComponents) → Optional[List[int]]¶ Modular comparison of this Molecule with another.
This permits detailed specification of which elements of the molecular information you want to use in the comparison.
Equality comparison is performed in several stages: First, at each atom position, a hash is computed that encompasses all local information that is specified to be used by the
components_bitmask
parameter. This hash is then used during graph isomorphism calculation to avoid finding an isomorphism that does not consider the specified factors.If an isomorphism is found, it is then validated. Bond orders and stereopermutators across both molecules are compared using the found isomorphism as an index map.
Shortcuts to
canonical_compare
ifcomponents_bitmask
matches the canonical components of both molecules (seecanonical_components
).- Parameters
other – The molecule to compare against
components_bitmask – The components of the molecule to use in the comparison
- Returns
None if the molecules are not isomorphic, a List[int] index mapping from self to other if the molecules are isomorphic.
>>> a = io.experimental.from_smiles("OCC") >>> b = io.experimental.from_smiles("SCC") >>> a == b False >>> # A and B are identical when considered purely by their graph >>> a.modular_isomorphism(b, AtomEnvironmentComponents.Connectivity) is not None True >>> # Another pair that is identical save for a stereopermutation >>> c = io.experimental.from_smiles("N[C@](Br)(O)C") >>> d = io.experimental.from_smiles("N[C@@](Br)(O)C") >>> c == d # Strict equality includes stereopermutation False >>> c.modular_isomorphism(d, AtomEnvironmentComponents.ElementsBondsAndShapes) is not None True
-
remove_atom
(self: scine_molassembler.Molecule, atom: int) → None¶ Remove an atom from the graph, including bonds to it, after checking that removing it is safe, i.e. the removal does not disconnect the graph.
- Warning
Invalidates all external atom and bond indices.
- Parameters
atom – Atom to remove
>>> m = Molecule() # Make H2 >>> [a for a in m.graph.atoms()] [0, 1] >>> m.graph.can_remove(0) # We can remove a hydrogen from H2 True >>> m.remove_atom(0) >>> m.graph.V # We are left with just a hydrogen atom 1 >>> m.graph.E 0
-
remove_bond
(*args, **kwargs)¶ Overloaded function.
remove_bond(self: scine_molassembler.Molecule, first_atom: int, second_atom: int) -> None
Remove a bond from the graph, after checking that removing it is safe, i.e. the removal does not disconnect the graph.
- warning
Invalidates all external atom and bond indices.
- param first_atom
First atom of the bond to be removed
- param second_atom
Second atom of the bond to be removed
>>> cyclopropane = io.experimental.from_smiles("C1CC1") >>> # In cyclopropane, we can remove a C-C bond without disconnecting the graph >>> cyclopropane.graph.can_remove(BondIndex(0, 1)) True >>> V_before = cyclopropane.graph.V >>> E_before = cyclopropane.graph.E >>> cyclopropane.remove_bond(BondIndex(0, 1)) >>> V_before - cyclopropane.graph.V # The number of atoms is unchanged 0 >>> E_before - cyclopropane.graph.E # We really only removed a bond 1 >>> # Note that now the valence of the carbon atoms where we removed >>> # the bond is... funky >>> cyclopropane.graph.degree(0) 3 >>> expected_bonds = [BondType.Single, BondType.Single, BondType.Single] >>> g = cyclopropane.graph >>> [g.bond_type(b) for b in g.bonds(0)] == expected_bonds True >>> cyclopropane.stereopermutators.option(0).shape == shapes.Shape.VacantTetrahedron True
remove_bond(self: scine_molassembler.Molecule, bond_index: scine_molassembler.BondIndex) -> None
Remove a bond from the graph, after checking that removing it is safe, i.e. the removal does not disconnect the graph.
- param bond_index
BondIndex
of the bond to be removed
-
remove_permutator
(self: scine_molassembler.Molecule, bond_index: scine_molassembler.BondIndex) → bool¶ Remove a bond stereopermutator from the molecule, if present
- Parameters
bond_index – Bond from which to remove the stereopermutator
- Returns
Whether a stereopermutator was removed
-
set_bond_type
(self: scine_molassembler.Molecule, first_atom: int, second_atom: int, bond_type: scine_molassembler.BondType) → bool¶ Change the bond type between two atoms. Inserts the bond if it doesn’t yet exist.
- Parameters
first_atom – First atom of the bond to be changed
second_atom – Second atom of the bond to be changed
bond_type – The new
BondType
- Returns
Whether the bond already existed
>>> # You really do have full freedom when it comes to your graphs: >>> h2 = Molecule() >>> _ = h2.set_bond_type(0, 1, BondType.Double) # Double bonded hydrogen atoms!
-
set_element_type
(self: scine_molassembler.Molecule, atom: int, element: Scine::Utils::ElementType) → None¶ Change the element type of an atom.
- Parameters
atom – Atom index of the atom to alter
element – New element type to set
>>> # Transform H2 into HF >>> import scine_utilities as utils >>> from copy import copy >>> H2 = Molecule() >>> HF = copy(H2) >>> HF.set_element_type(0, utils.ElementType.F) >>> HF == H2 False
-
set_shape_at_atom
(self: scine_molassembler.Molecule, atom: int, shape: scine_molassembler.shapes.Shape) → None¶ Change the local shape at an atom.
This sets the local shape at a specific atom index. There are a number of cases that this function treats differently, besides faulty arguments: If there is already a AtomStereopermutator instantiated at this atom index, its underlying shape is altered. If there is no AtomStereopermutator at this index, one is instantiated. In all cases, new or modified stereopermutators are default-assigned if there is only one possible assignment.
>>> # Make methane square planar >>> from copy import copy >>> methane = io.experimental.from_smiles("C") >>> square_planar_methane = copy(methane) >>> square_planar_methane.set_shape_at_atom(0, shapes.Shape.Square) >>> methane == square_planar_methane False
-
property
stereopermutators
¶ Read only access to the list of stereopermutators
- Return type
-
thermalize_stereopermutator
(self: scine_molassembler.Molecule, atom_index: int, thermalization: bool = True) → None¶ Change the thermalization at an atom stereopermutator
Alters the thermalization of stereopermutations at an atom stereopermutator.
- Parameters
atom_index – Atom whose atom stereopermutator’s thermalization to change
thermalization – New status of thermalization to set
-
-
class
scine_molassembler.
AtomEnvironmentComponents
¶ Denotes information parts of molecules. Relevant for molecular comparison and hashing.
Members:
Connectivity : Consider only the graph
ElementTypes : Element types
BondOrders : Bond orders
Shapes : Shapes
Stereopermutations : Stereopermutations
ElementsAndBonds : Consider element types and bond orders
ElementsBondsAndShapes : Consider element types, bond orders and shapes
All : Consider element types, bond orders, shapes and stereopermutations
-
All
= <AtomEnvironmentComponents.All: 15>¶
-
BondOrders
= <AtomEnvironmentComponents.BondOrders: 2>¶
-
Connectivity
= <AtomEnvironmentComponents.Connectivity: 0>¶
-
ElementTypes
= <AtomEnvironmentComponents.ElementTypes: 1>¶
-
ElementsAndBonds
= <AtomEnvironmentComponents.ElementsAndBonds: 3>¶
-
ElementsBondsAndShapes
= <AtomEnvironmentComponents.ElementsBondsAndShapes: 7>¶
-
Shapes
= <AtomEnvironmentComponents.Shapes: 4>¶
-
Stereopermutations
= <AtomEnvironmentComponents.Stereopermutations: 8>¶
-
__init__
(self: scine_molassembler.AtomEnvironmentComponents, value: int) → None¶
-
property
name
¶
-
property
value
¶
-