Molassembler  1.0.0
Molecule graph and conformer library
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Scine::Molassembler::AtomStereopermutator Class Reference

Handles the steric permutation of substituents of a non-terminal central atom. More...

#include <AtomStereopermutator.h>

Data Structures

class  Impl
 

Public Types

using ShapeMap = Temple::StrongIndexFlatMap< SiteIndex, Shapes::Vertex >
 
using PropagatedState = std::tuple< RankingInformation, Stereopermutators::Abstract, Stereopermutators::Feasible, ShapeMap >
 Old state dumped upon propagation.
 
using MinimalChiralConstraint = std::array< boost::optional< SiteIndex >, 4 >
 Site index sequence defining a chiral constraint. If a site index is None, then it denotes the position of the central index.
 

Public Member Functions

Special member functions
 AtomStereopermutator (AtomStereopermutator &&other) noexcept
 Construct an AtomStereopermutator. More...
 
AtomStereopermutatoroperator= (AtomStereopermutator &&other) noexcept
 Construct an AtomStereopermutator. More...
 
 AtomStereopermutator (const AtomStereopermutator &other)
 Construct an AtomStereopermutator. More...
 
AtomStereopermutatoroperator= (const AtomStereopermutator &other)
 Construct an AtomStereopermutator. More...
 
 ~AtomStereopermutator ()
 Construct an AtomStereopermutator. More...
 
 AtomStereopermutator (const Graph &graph, Shapes::Shape shape, AtomIndex centerAtom, RankingInformation ranking)
 Construct an AtomStereopermutator. More...
 
Modifiers
void assign (boost::optional< unsigned > assignment)
 Changes the assignment of the stereopermutator. More...
 
void assignRandom (Random::Engine &engine)
 Assign the Stereopermutator randomly using relative statistical weights. More...
 
void applyPermutation (const std::vector< AtomIndex > &permutation)
 Applies an atom index permutation. More...
 
boost::optional< ShapeMapfit (const Graph &graph, const AngstromPositions &angstromWrapper)
 Determine the shape and assignment realized in positions. More...
 
boost::optional< PropagatedStatepropagate (const Graph &graph, RankingInformation newRanking, boost::optional< Shapes::Shape > shapeOption)
 Propagate the stereocenter state through a possible ranking change. More...
 
void propagateVertexRemoval (AtomIndex removedIndex)
 Adapts atom indices in the internal state to the removal of an atom. More...
 
void setShape (Shapes::Shape shape, const Graph &graph)
 Change the underlying shape of the permutator. More...
 
Observers
double angle (SiteIndex i, SiteIndex j) const
 Fetches angle between binding sites in the idealized shape. More...
 
boost::optional< unsigned > assigned () const
 Returns the permutation index within the set of feasible permutations, if set. More...
 
AtomIndex placement () const
 Returns the atom this permutator is placed on. More...
 
boost::optional< unsigned > indexOfPermutation () const
 Returns IOP within the set of symbolic ligand permutations. More...
 
std::vector
< MinimalChiralConstraint
minimalChiralConstraints (bool enforce=false) const
 Returns a minimal representation of chiral constraints. More...
 
std::string info () const
 Returns an information string for diagnostic purposes. More...
 
std::string rankInfo () const
 Returns an information string for ranking equality checking purposes. More...
 
std::vector< std::vector
< SiteIndex > > 
siteGroups () const
 Returns site indices grouped by rotational interconversion. More...
 
const Stereopermutators::AbstractgetAbstract () const
 Returns the underlying feasible stereopermutations object. More...
 
const Stereopermutators::FeasiblegetFeasible () const
 Returns the underlying feasible stereopermutations object. More...
 
const RankingInformationgetRanking () const
 Returns the underlying ranking. More...
 
Shapes::Shape getShape () const
 Returns the underlying shape. More...
 
const ShapeMapgetShapePositionMap () const
 Yields the mapping from site indices to shape positions. More...
 
unsigned numAssignments () const
 Returns the number of possible assignments. More...
 
unsigned numStereopermutations () const
 Returns the number of possible stereopermutations. More...
 
Operators
bool operator== (const AtomStereopermutator &other) const
 Checks whether the underlying shape, central atom index, number of stereopermutations and current assignment match. More...
 
bool operator!= (const AtomStereopermutator &other) const
 Inverts operator ==.
 
bool operator< (const AtomStereopermutator &other) const
 Lexicographically compares the central atom index, the shape, the number of stereopermutations, and the current assignment. More...
 

Static Public Member Functions

Shape picking
static Shapes::Shape up (Shapes::Shape shape)
 Picks a shape retaining as much chiral state as possible on a shape size increase. More...
 
static Shapes::Shape down (Shapes::Shape shape, Shapes::Vertex removedVertex)
 Picks a shape retaining as much chiral state as possible on a shape size decrease. More...
 

Private Attributes

std::unique_ptr< ImplpImpl_
 

Detailed Description

Handles the steric permutation of substituents of a non-terminal central atom.

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.

E.g. a square planar AABC ligand set will have an A-A cis stereopermutation that occurs twice as often as the A-A trans stereopermutation.

Note
An instance of this class on a given central atom does not indicate that that atom is a stereocenter. That is only the case if there are multiple stereopermutations of the ranked substituents / ligands.

Constructor & Destructor Documentation

Scine::Molassembler::AtomStereopermutator::AtomStereopermutator ( AtomStereopermutator &&  other)
noexcept

Construct an AtomStereopermutator.

Parameters
graphThe molecule's graph. This information is needed to model haptic ligands.
shapeThe local idealized shape to model. Typically the result of Molecule's inferShape.
centerAtomThe atom index within the molecule that is the center of the local idealized shape
rankingThe ranking of the central atom's substituents and ligand sites. Typically the result of Molecule's rankPriority.

Complexity \(L\cdot S!\) where \(L\) is the number of links and \(S\) is the size of shape

Postcondition
The stereopermutator is unassigned
Scine::Molassembler::AtomStereopermutator::AtomStereopermutator ( const AtomStereopermutator other)

Construct an AtomStereopermutator.

Parameters
graphThe molecule's graph. This information is needed to model haptic ligands.
shapeThe local idealized shape to model. Typically the result of Molecule's inferShape.
centerAtomThe atom index within the molecule that is the center of the local idealized shape
rankingThe ranking of the central atom's substituents and ligand sites. Typically the result of Molecule's rankPriority.

Complexity \(L\cdot S!\) where \(L\) is the number of links and \(S\) is the size of shape

Postcondition
The stereopermutator is unassigned
Scine::Molassembler::AtomStereopermutator::~AtomStereopermutator ( )

Construct an AtomStereopermutator.

Parameters
graphThe molecule's graph. This information is needed to model haptic ligands.
shapeThe local idealized shape to model. Typically the result of Molecule's inferShape.
centerAtomThe atom index within the molecule that is the center of the local idealized shape
rankingThe ranking of the central atom's substituents and ligand sites. Typically the result of Molecule's rankPriority.

Complexity \(L\cdot S!\) where \(L\) is the number of links and \(S\) is the size of shape

Postcondition
The stereopermutator is unassigned
Scine::Molassembler::AtomStereopermutator::AtomStereopermutator ( const Graph graph,
Shapes::Shape  shape,
AtomIndex  centerAtom,
RankingInformation  ranking 
)

Construct an AtomStereopermutator.

Parameters
graphThe molecule's graph. This information is needed to model haptic ligands.
shapeThe local idealized shape to model. Typically the result of Molecule's inferShape.
centerAtomThe atom index within the molecule that is the center of the local idealized shape
rankingThe ranking of the central atom's substituents and ligand sites. Typically the result of Molecule's rankPriority.

Complexity \(L\cdot S!\) where \(L\) is the number of links and \(S\) is the size of shape

Postcondition
The stereopermutator is unassigned

Member Function Documentation

double Scine::Molassembler::AtomStereopermutator::angle ( SiteIndex  i,
SiteIndex  j 
) const

Fetches angle between binding sites in the idealized shape.

Precondition
The stereopermutator must be assigned
Parameters
iSite index one
jSite index two

Complexity \(\Theta(1)\)

Exceptions
std::runtime_errorIf the stereopermutator is unassigned
std::out_of_rangeIf any site index is out of range
See Also
getRanking()
void Scine::Molassembler::AtomStereopermutator::applyPermutation ( const std::vector< AtomIndex > &  permutation)

Applies an atom index permutation.

Complexity \(\Theta(1)\)

Parameters
permutationThe permutation to apply
void Scine::Molassembler::AtomStereopermutator::assign ( boost::optional< unsigned >  assignment)

Changes the assignment of the stereopermutator.

Parameters
assignmentThe new assignment of the stereopermutator. May be boost::none, which sets the chiral state as indeterminate. Must be less than the number of assignments if not None.

Complexity \(\Theta(1)\) if assignment is boost::none. \(\Theta(S)\) otherwise

boost::optional<unsigned> Scine::Molassembler::AtomStereopermutator::assigned ( ) const

Returns the permutation index within the set of feasible permutations, if set.

Returns the information of whether the stereopermutator is assigned or not, and if so, which assignment it is.

Complexity \(\Theta(1)\)

void Scine::Molassembler::AtomStereopermutator::assignRandom ( Random::Engine engine)

Assign the Stereopermutator randomly using relative statistical weights.

Stereopermutations are generated with relative statistical occurrence weights. The assignment is then chosen from the possible stereopermutations with a discrete distribution whose weights are the corresponding relative statistical occurrences.

Complexity \(\Theta(1)\) if assignment is boost::none. \(\Theta(S)\) otherwise

Note
If the stereocenter is already assigned, it is reassigned.
The state of the passed PRNG is advanced.
static Shapes::Shape Scine::Molassembler::AtomStereopermutator::down ( Shapes::Shape  shape,
Shapes::Vertex  removedVertex 
)
static

Picks a shape retaining as much chiral state as possible on a shape size decrease.

Complexity \(O(S!)\) if uncached, \(\Theta(1)\) otherwise

Exceptions
std::logic_errorIf there are no smaller shapes
boost::optional<ShapeMap> Scine::Molassembler::AtomStereopermutator::fit ( const Graph graph,
const AngstromPositions angstromWrapper 
)

Determine the shape and assignment realized in positions.

The shape and assignment are determined from Cartesian coordinates by the probability that the continuous shape measure with respect to an idealized shape is least likely to be drawn from a sample of randomly generated vertex positions. For shapes with many vertices, generating enough samples for a statistical argument is too expensive, and the shape whose continuous shape measure is lowest is chosen instead.

Parameters
graphThe molecule's graph which this permutator helps model
angstromWrapperThe wrapped positions
Returns
A mapping of site indices to shape vertices if a stereopermutation could be found, None otherwise

Complexity \(\Theta(S!)\)

const Stereopermutators::Abstract& Scine::Molassembler::AtomStereopermutator::getAbstract ( ) const

Returns the underlying feasible stereopermutations object.

Complexity \(\Theta(1)\)

Note
This is library-internal and not part of the public API
const Stereopermutators::Feasible& Scine::Molassembler::AtomStereopermutator::getFeasible ( ) const

Returns the underlying feasible stereopermutations object.

Complexity \(\Theta(1)\)

Note
This is library-internal and not part of the public API
const RankingInformation& Scine::Molassembler::AtomStereopermutator::getRanking ( ) const

Returns the underlying ranking.

Complexity \(\Theta(1)\)

Shapes::Shape Scine::Molassembler::AtomStereopermutator::getShape ( ) const

Returns the underlying shape.

Complexity \(\Theta(1)\)

const ShapeMap& Scine::Molassembler::AtomStereopermutator::getShapePositionMap ( ) const

Yields the mapping from site indices to shape positions.

Complexity \(\Theta(1)\)

Exceptions
std::logic_errorif the stereopermutator is unassigned.
boost::optional<unsigned> Scine::Molassembler::AtomStereopermutator::indexOfPermutation ( ) const

Returns IOP within the set of symbolic ligand permutations.

This is different to the assignment. The assignment denotes the index within the set of possible (more specifically, not obviously infeasible) stereopermutations.

Complexity \(\Theta(1)\)

std::string Scine::Molassembler::AtomStereopermutator::info ( ) const

Returns an information string for diagnostic purposes.

Complexity \(\Theta(1)\)

std::vector<MinimalChiralConstraint> Scine::Molassembler::AtomStereopermutator::minimalChiralConstraints ( bool  enforce = false) const

Returns a minimal representation of chiral constraints.

Every minimal representation consists only of site indices. If no site index is present, this position is the location of the central atom.

The minimal representation assumes that all shape tetrahedra are defined as Positive targets, which is checked in the shapes tests.

Complexity \(\Theta(T)\) where \(T\) is the number of tetrahedra defined for the modeled shape

Parameters
enforceEmit minimal representations of chiral constraints even if the stereopermutator does not have any chiral state, i.e. numStereopermutators() <= 1, as long as it is assigned.
unsigned Scine::Molassembler::AtomStereopermutator::numAssignments ( ) const

Returns the number of possible assignments.

The number of possible assignments is the number of non-superposable arrangements of the abstract ligand case reduced by trans-arranged multidentate pairs where the bridge length is too short or overlapping haptic cones.

For instance, if octahedral M[(A-A)3], there are four abstract arrangements

  • trans-trans-trans
  • trans-cis-cis
  • 2x cis-cis-cis (Δ and Λ isomers, ship propeller-like)

However, the number of stereopermutations for a concrete case in which the bridges are too short to allow trans bonding is reduced by all arrangements containing a trans-bonded bidentate ligand, i.e. only Δ and Λ remain. The number of assignments is then only two.

This is the upper exclusive bound on Some-type arguments to assign().

Complexity \(\Theta(1)\)

unsigned Scine::Molassembler::AtomStereopermutator::numStereopermutations ( ) const

Returns the number of possible stereopermutations.

The number of possible stereopermutations is the number of non-superposable arrangements of the abstract ligand case without removing trans-arranged multidentate pairs or overlapping haptic cones.

For instance, if octahedral M[(A-A)3], there are four abstract arrangements

  • trans-trans-trans
  • trans-cis-cis
  • 2x cis-cis-cis (Δ and Λ isomers, ship propeller-like chirality)

However, the number of assignments for a concrete case in which the bridges are too short to allow trans binding is reduced by all arrangements containing a trans-bonded bidentate ligand, i.e. only Δ and Λ remain.

Fetches the number of permutations determined by symbolic ligand calculation, not considering linking or haptic ligand cones.

Complexity \(\Theta(1)\)

bool Scine::Molassembler::AtomStereopermutator::operator< ( const AtomStereopermutator other) const

Lexicographically compares the central atom index, the shape, the number of stereopermutations, and the current assignment.

Parameters
otherThe other atom stereopermutator to compare against
Returns
Which atom stereopermutator is lexicographically smaller
AtomStereopermutator& Scine::Molassembler::AtomStereopermutator::operator= ( AtomStereopermutator &&  other)
noexcept

Construct an AtomStereopermutator.

Parameters
graphThe molecule's graph. This information is needed to model haptic ligands.
shapeThe local idealized shape to model. Typically the result of Molecule's inferShape.
centerAtomThe atom index within the molecule that is the center of the local idealized shape
rankingThe ranking of the central atom's substituents and ligand sites. Typically the result of Molecule's rankPriority.

Complexity \(L\cdot S!\) where \(L\) is the number of links and \(S\) is the size of shape

Postcondition
The stereopermutator is unassigned
AtomStereopermutator& Scine::Molassembler::AtomStereopermutator::operator= ( const AtomStereopermutator other)

Construct an AtomStereopermutator.

Parameters
graphThe molecule's graph. This information is needed to model haptic ligands.
shapeThe local idealized shape to model. Typically the result of Molecule's inferShape.
centerAtomThe atom index within the molecule that is the center of the local idealized shape
rankingThe ranking of the central atom's substituents and ligand sites. Typically the result of Molecule's rankPriority.

Complexity \(L\cdot S!\) where \(L\) is the number of links and \(S\) is the size of shape

Postcondition
The stereopermutator is unassigned
bool Scine::Molassembler::AtomStereopermutator::operator== ( const AtomStereopermutator other) const

Checks whether the underlying shape, central atom index, number of stereopermutations and current assignment match.

Parameters
otherThe other atom stereopermutator to compare against
Returns
Whether the underlying shape, central atom index, number of stereopermutations and current assignment match
AtomIndex Scine::Molassembler::AtomStereopermutator::placement ( ) const

Returns the atom this permutator is placed on.

Complexity \(\Theta(1)\)

boost::optional<PropagatedState> Scine::Molassembler::AtomStereopermutator::propagate ( const Graph graph,
RankingInformation  newRanking,
boost::optional< Shapes::Shape shapeOption 
)

Propagate the stereocenter state through a possible ranking change.

In case a graph modification changes the ranking of this stereopermutator's substituents, it must be redetermined whether the new configuration is a stereopermutator and if so, which assignment corresponds to the previous one.

Complexity \(L\cdot S!\) where \(L\) is the number of links and \(S\) is the size of shape

void Scine::Molassembler::AtomStereopermutator::propagateVertexRemoval ( AtomIndex  removedIndex)

Adapts atom indices in the internal state to the removal of an atom.

Atom indices are adapted to a graph-level removal of an atom. The removed index is changed to a placeholder value.

Complexity \(\Theta(1)\)

std::string Scine::Molassembler::AtomStereopermutator::rankInfo ( ) const

Returns an information string for ranking equality checking purposes.

Complexity \(\Theta(1)\)

void Scine::Molassembler::AtomStereopermutator::setShape ( Shapes::Shape  shape,
const Graph graph 
)

Change the underlying shape of the permutator.

Complexity \(L\cdot S!\) where \(L\) is the number of links and \(S\) is the size of shape

Postcondition
The permutator is unassigned (chiral state is discarded)
std::vector<std::vector<SiteIndex> > Scine::Molassembler::AtomStereopermutator::siteGroups ( ) const

Returns site indices grouped by rotational interconversion.

Precondition
The stereopermutator must be assigned

Complexity \(\Theta(S^2)\)

Exceptions
std::runtime_errorIf the stereopermutator is unassigned
static Shapes::Shape Scine::Molassembler::AtomStereopermutator::up ( Shapes::Shape  shape)
static

Picks a shape retaining as much chiral state as possible on a shape size increase.

Complexity \(O(S!)\) if uncached, \(\Theta(1)\) otherwise

Exceptions
std::logic_errorIf there are no larger shapes

The documentation for this class was generated from the following file: