Stereopermutators

Currently, this library models two types of stereocenters: The first models the atom-centric relative placement of substituents across geometries ranging from linear to square-antiprismatic. The second type models configurational permutations owing to bond-centric rotational barriers.

In order to determine whether a particular atom in a molecule is an atom-centric stereocenter, its substituents are ranked according to a nearly-compliant implementation of the IUPAC sequence rules as laid out in the 2013 Blue Book. The rankings are transformed into an abstract substituent case (e.g. octahedral (A-A)BBCD) and a symbolic computation is carried out to determine the number of permutations that are not superimposable via spatial rotations within the idealized shape. The set of resulting permutations of the substituent symbols is called the set of stereopermutations. If this set contains more than one stereopermutation, then the atom is an atom-centric stereocenter under that shape.

If substituents are haptic or multidentate, an additional algorithm removes stereopermutations it deems clearly impossible. All bridge lengths between pairs of chelating atoms of a multidentate ligand are checked against the atom pair’s bite angle within the idealized shape. Additionally, haptic ligands’ cones are checked to ensure they do not overlap. Indices within the set of not clearly impossible stereopermutations are called assignments.

Bond-centric stereocenters model rotational configurations of arbitrary combinations of two shapes and their respective fused shape positions. The fused shape positions of each side affect the overall permutations if the shape has multiple position groups. For instance, this is the case in square pyramidal shapes, where there are axial and equatorial shape positions.

class scine_molassembler.StereopermutatorList

Manages all stereopermutators that are part of a molecule.

>>> # A sample molecule with one stereogenic atom and bond stereopermutator each
>>> mol = io.experimental.from_smiles("CC=C[C@](F)(Cl)[H]")
>>> permutators = mol.stereopermutators
>>> is_stereogenic = lambda p: p.num_assignments > 1
>>> atom_permutators = permutators.atom_stereopermutators()
>>> bond_permutators = permutators.bond_stereopermutators()
>>> stereogenic_atom_permutators = [p for p in atom_permutators if is_stereogenic(p)]
>>> stereogenic_bond_permutators = [p for p in bond_permutators if is_stereogenic(p)]
>>> # Atom stereopermutators are instantiated on every non-terminal atom
>>> permutators.A() > len(stereogenic_atom_permutators)
True
>>> len(stereogenic_atom_permutators) # But only one of them is stereogenic here
1
>>> # Bond stereopermutators are instantiated only where they are stereogenic
>>> # or conserve information on relative spatial orientation
>>> permutators.B() > len(stereogenic_bond_permutators)
False
>>> len(stereogenic_bond_permutators)
1
>>> permutators.has_unassigned_permutators() # The stereo of the double bond is unspecified
True
>>> permutators.has_zero_assignment_permutators()
False
>>> double_bond_index = BondIndex(1, 2)
>>> assert mol.graph.bond_type(double_bond_index) == BondType.Double
>>> # When looking up a stereopermutator, remember that you can get None back
>>> bond_stereopermutator = permutators.option(double_bond_index)
>>> bond_stereopermutator is None
False
>>> is_stereogenic(bond_stereopermutator)
True
A(self: scine_molassembler.StereopermutatorList) int

Returns the number of AtomStereopermutator instances stored

B(self: scine_molassembler.StereopermutatorList) int

Returns the number of BondStereopermutator instances stored

__init__(*args, **kwargs)
atom_stereopermutators(self: scine_molassembler.StereopermutatorList) Iterator

Returns a range of all AtomStereopermutator

bond_stereopermutators(self: scine_molassembler.StereopermutatorList) Iterator

Returns a range of all BondStereopermutator

empty(self: scine_molassembler.StereopermutatorList) bool

Whether the list is empty or not. Remember that since atom stereopermutators are instantiated on all non-terminal atoms, this is rare.

>>> h2 = Molecule() # Default constructor makes the diatomic hydrogen molecule
>>> h2.stereopermutators.empty() # Both atoms are terminal here, so no permutators
True
has_unassigned_permutators(self: scine_molassembler.StereopermutatorList) bool

Returns whether the list contains any stereopermutators that are unassigned

has_zero_assignment_permutators(self: scine_molassembler.StereopermutatorList) bool

Returns whether the list contains any stereopermutators that have no possible stereopermutations

option(*args, **kwargs)

Overloaded function.

  1. option(self: scine_molassembler.StereopermutatorList, atom: int) -> Optional[scine_molassembler.AtomStereopermutator]

    Fetches a read-only option to an AtomStereopermutator, if present on this atom index

  2. option(self: scine_molassembler.StereopermutatorList, bond_index: scine_molassembler.BondIndex) -> Optional[scine_molassembler.BondStereopermutator]

    Fetches a read-only option to a BondStereopermutator, if present on this atom index

Atom stereopermutators

class scine_molassembler.AtomStereopermutator

This class handles the permutation of ranked ligands around a central atom. It models its haptic ligands’ binding sites and bridges in multidentate ligands in order to decide which stereopermutations are feasible. A stereopermutation may be infeasible, i.e. not realizable in three-dimensional space, if either haptic ligands would intersect due to too close ligand angles for their spatial heft, or if a multidentate ligand’s bridge length between binding sites were too short to match the angle. The list of stereopermutations reduced by infeasible stereopermutations is then re-indexed and those indices into the list are called assignments.

A Stereopermutator can be unassigned, i.e. the distinct stereopermutation that the substituents are can be indeterminate. If you choose to generate conformers for a molecule that includes unassigned stereopermutators, every conformer will choose an assignment from the pool of feasible assignments randomly, but consistent with relative statistical occurrence weights.

Stereopermutator instances themselves are nonmodifiable. To change them, you have to make changes at the molecule level.

>>> # Enantiomeric pair of asymmetric tetrahedral carbon atoms
>>> import scine_utilities as utils
>>> asym_carbon = io.experimental.from_smiles("N[C@](Br)(O)F")
>>> carbon_index = asym_carbon.graph.atoms_of_element(utils.ElementType.C)[0]
>>> carbon_stereopermutator = asym_carbon.stereopermutators.option(carbon_index)
>>> assert carbon_stereopermutator is not None
>>> carbon_stereopermutator.shape == shapes.Shape.Tetrahedron
True
>>> carbon_stereopermutator.assigned is not None
True
>>> enantiomer = io.experimental.from_smiles("N[C@@](Br)(O)F")
>>> assert enantiomer.graph.element_type(carbon_index) == utils.ElementType.C
>>> enantiomer_stereopermutator = enantiomer.stereopermutators.option(carbon_index)
>>> enantiomer_stereopermutator.assigned is not None
True
>>> carbon_stereopermutator.assigned != enantiomer_stereopermutator.assigned
True
__init__(self: scine_molassembler.AtomStereopermutator, graph: scine_molassembler.Graph, shape: scine_molassembler.shapes.Shape, placement: int, ranking: scine_molassembler.RankingInformation) None

Instantiate an atom stereopermutator at a particular position in a graph

angle(self: scine_molassembler.AtomStereopermutator, site_index_i: int, site_index_j: int) float

Fetches the angle between substituent site indices in radians

>>> # The tetrahedron angle
>>> import math
>>> tetrahedron_angle = 2 * math.atan(math.sqrt(2))
>>> methane = io.experimental.from_smiles("C")
>>> a = methane.stereopermutators.option(0)
>>> math.isclose(a.angle(0, 1), tetrahedron_angle)
True
property assigned

The assignment integer if assigned, None otherwise.

>>> # A stereo-unspecified tetrahedral asymmetric carbon atom
>>> import scine_utilities as utils
>>> asymmetric_carbon = io.experimental.from_smiles("NC(Br)(O)F")
>>> carbon_index = asymmetric_carbon.graph.atoms_of_element(utils.ElementType.C)[0]
>>> stereopermutator = asymmetric_carbon.stereopermutators.option(carbon_index)
>>> assert stereopermutator is not None
>>> stereopermutator.num_assignments == 2
True
>>> stereopermutator.assigned is None
True
property index_of_permutation

The index of permutation if assigned, otherwise None. Indices of permutation are the abstract index of permutation within the set of permutations that do not consider feasibility. This is not necessarily equal to the assignment index.

>>> # The shipscrew lambda and delta isomers where trans-ligation is impossible
>>> shipscrew_smiles = "[Fe@OH1+3]123(OC(=O)C(=O)O1)(OC(=O)C(=O)O2)OC(=O)C(=O)O3"
>>> shipscrew = io.experimental.from_smiles(shipscrew_smiles)
>>> permutator = shipscrew.stereopermutators.option(0)
>>> assert permutator is not None
>>> permutator.num_stereopermutations # Number of abstract permutations
4
>>> permutator.num_assignments # Number of spatially feasible permutations
2
>>> permutator.index_of_permutation
2
>>> permutator.assigned
1
>>> shipscrew.assign_stereopermutator(0, None) # Dis-assign the stereopermutator
>>> permutator = shipscrew.stereopermutators.option(0)
>>> assert permutator is not None
>>> permutator.index_of_permutation is None
True
>>> permutator.assigned is None
True
property num_assignments

The number of feasible assignments. See index_of_permutation.

property num_stereopermutations

The number of stereopermutations. See index_of_permutation.

property placement

The central atom this permutator is placed on

property ranking

Get the underlying ranking state of substituents

Return type

RankingInformation

property shape

Returns the underlying shape

Return type

shapes.Shape

site_groups(self: scine_molassembler.AtomStereopermutator) List[List[int]]

Site indices grouped by whether their shape vertices are interconvertible by rotation. Can be used to distinguish e.g. apical / equatorial ligands in a square pyramid shape.

Raises

RuntimeError if the stereopermutator is not assigned

property thermalized

Whether the stereopermutations are thermalized

property vertex_map

Bond stereopermutators

class scine_molassembler.BondStereopermutator

Handles specific relative arrangements of two atom stereopermutators joined by a bond. This includes, importantly, E/Z stereocenters at double bonds.

>>> # The bond stereopermutator in but-2z-ene
>>> z_butene = io.experimental.from_smiles("C/C=C\C")
>>> bond_index = BondIndex(1, 2)
>>> assert z_butene.graph.bond_type(bond_index) == BondType.Double
>>> permutator = z_butene.stereopermutators.option(bond_index)
>>> assert permutator is not None
>>> permutator.assigned is not None
True
>>> permutator.num_assignments
2
class Alignment

How dihedrals are aligned in the generation of stereopermutations

Members:

Eclipsed : At least two shape vertices eclipse one another along the axis

Staggered : At least one pair of substituents are staggered along the axis

EclipsedAndStaggered : Both eclipsed and staggered alignments are generated

BetweenEclipsedAndStaggered : Offset exactly halfway between eclipsed and staggered alignments

BetweenEclipsedAndStaggered = <Alignment.BetweenEclipsedAndStaggered: 3>
Eclipsed = <Alignment.Eclipsed: 0>
EclipsedAndStaggered = <Alignment.EclipsedAndStaggered: 2>
Staggered = <Alignment.Staggered: 1>
__init__(self: scine_molassembler.BondStereopermutator.Alignment, value: int) None
property name
property value
class FittingMode

Differentiates how viable assignments are chosen during fitting

Members:

Thresholded : Positions must be close to the idealized assignment geometry

Nearest : The assignment closest to the idealized geometry is chosen

Nearest = <FittingMode.Nearest: 1>
Thresholded = <FittingMode.Thresholded: 0>
__init__(self: scine_molassembler.BondStereopermutator.FittingMode, value: int) None
property name
property value
__init__(*args, **kwargs)
property assigned

An integer indicating the assignment of the stereopermutator or None if the stereopermutator is unassigned.

>>> # An unassigned bond stereopermutator
>>> butene = io.experimental.from_smiles("CC=CC")
>>> bond_index = BondIndex(1, 2)
>>> assert butene.graph.bond_type(bond_index) == BondType.Double
>>> permutator = butene.stereopermutators.option(bond_index)
>>> assert permutator is not None
>>> permutator.assigned is None
True
>>> permutator.num_assignments
2
property composite

The underlying stereopermutation generating shape composite

dihedral(self: scine_molassembler.BondStereopermutator, stereopermutator_a: scine_molassembler.AtomStereopermutator, site_index_a: int, stereopermutator_b: scine_molassembler.AtomStereopermutator, site_index_b: int) float

Returns the dihedral angle between two sites of the constituting atom stereopermutators in radians

You can glean site indices from the individual constituting atom stereopermutators’ rankings.

Parameters
  • stereopermutator_a – One constituting AtomStereopermutator

  • site_index_a – The site index of stereopermutator_a starting the dihedral sequence

  • stereopermutator_b – The other constituting AtomStereopermutator

  • site_index_b – The site index of stereopermutator_b ending the dihedral sequence

property index_of_permutation

Returns an integer indicating the index of permutation if the stereopermutator is assigned or None if the stereopermutator is unassigned.

>>> # A case in which the number of abstract and feasible permutations
>>> # differ: bond stereopermutators in small cycles (<= 6)
>>> benzene = io.experimental.from_smiles("C1=CC=CC=C1")
>>> permutators = benzene.stereopermutators.bond_stereopermutators()
>>> has_two_stereopermutations = lambda p: p.num_stereopermutations == 2
>>> has_one_assignment = lambda p: p.num_assignments == 1
>>> all(map(has_two_stereopermutations, permutators))
True
>>> all(map(has_one_assignment, permutators))
True
property num_assignments

The number of assignments. Valid assignment indices range from 0 to this number minus one.

property num_stereopermutations

Returns the number of stereopermutations.

property placement

The edge this stereopermutator is placed on.