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

Class performing spatial modeling of molecules. More...

#include <SpatialModel.h>

Data Structures

struct  BoundsMatrixHelper
 

Public Types

Member types
template<unsigned size>
using BoundsMapType = std::unordered_map< std::array< AtomIndex, size >, ValueBounds, boost::hash< std::array< AtomIndex, size > > >
 ValueBounds container for various internal coordinates.
 
using FixedPositionsMapType = std::unordered_map< AtomIndex, Utils::Position >
 Type used to store fixed positions in angstrom.
 
using BoundsMatrix = Eigen::MatrixXd
 Type used to represent atom-pairwise distance bounds constructed from bounds on internal coordinates. More...
 

Public Member Functions

Special member functions
 SpatialModel (const Molecule &molecule, const Configuration &configuration)
 Model a molecule into internal coordinate bounds stored internally. More...
 
Modifiers
void setBondBoundsIfEmpty (std::array< AtomIndex, 2 > bondIndices, ValueBounds bounds)
 Sets bond bounds to exact value bounds if unset. More...
 
void setAngleBoundsIfEmpty (std::array< AtomIndex, 3 > angleIndices, ValueBounds bounds)
 Adds angle bounds to the model if unset. More...
 
void setDihedralBoundsIfEmpty (std::array< AtomIndex, 4 > dihedralIndices, ValueBounds bounds)
 Adds the dihedral bounds to the model if unset. More...
 
void addAtomStereopermutatorInformation (const AtomStereopermutator &permutator, const PrivateGraph &graph, double looseningMultiplier, const std::unordered_map< AtomIndex, Utils::Position > &fixedAngstromPositions, bool forceChiralConstraintEmission)
 Adds angle information to the internal coordinate bounds and collects chiral constraints. More...
 
void addBondStereopermutatorInformation (const BondStereopermutator &permutator, const AtomStereopermutator &stereopermutatorA, const AtomStereopermutator &stereopermutatorB, double looseningMultiplier, const std::unordered_map< AtomIndex, Utils::Position > &fixedAngstromPositions)
 Adds dihedral information to the internal coordinate bounds and collects dihedral constraints. More...
 
bool modelPartiallyFixedBond (const BondStereopermutator &permutator, const std::pair< const AtomStereopermutator &, const AtomStereopermutator & > &atomStereopermutators, const std::unordered_map< AtomIndex, Utils::Position > &fixedAngstromPositions)
 Sets bond bounds to exact value bounds if unset. More...
 
Information
std::vector< ChiralConstraintgetChiralConstraints () const
 Yields all collected chiral constraints. More...
 
std::vector< DihedralConstraintgetDihedralConstraints () const
 Yields all collected dihedral constraints. More...
 
BoundsMatrix makePairwiseBounds () const
 Generate unsmoothed atom-pairwise distance bounds matrix. More...
 
std::string dumpGraphviz () const
 Generates a string graphviz representation of the modeled molecule. More...
 
void writeGraphviz (const std::string &filename) const
 Writes a graphviz representation of a modeled molecule to a file. More...
 

Private Member Functions

Private methods
void addDefaultAngles_ ()
 Adds [0, π] default angle bounds for all bonded atom triples.
 
void addDefaultDihedrals_ ()
 Adds [0, π] default dihedrals to the model (of connected sequences) More...
 
void instantiateMissingAtomStereopermutators_ (Random::Engine &engine)
 
void modelBondDistances_ (const FixedPositionsMapType &fixedAngstromPositions, double looseningFactor)
 Add bond distances to the underlying model.
 
void modelFlatCycles_ (const FixedPositionsMapType &fixedAngstromPositions, double looseningFactor)
 Add precise information to the model for cycles that are flat.
 
void modelSpirocenters_ (const FixedPositionsMapType &fixedAngstromPositions)
 Adds spirocenter modelling to the data set. More...
 

Private Attributes

const Moleculemolecule_
 
BoundsMapType< 2 > constraints_
 Constraints by fixed positions.
 
BoundsMapType< 2 > bondBounds_
 Bond distance value bounds on index pairs.
 
BoundsMapType< 3 > angleBounds_
 Angle value bounds on index triples.
 
BoundsMapType< 4 > dihedralBounds_
 Dihedral value bounds on index quadruplets.
 
std::vector< ChiralConstraintchiralConstraints_
 Chiral constraints.
 
std::vector< DihedralConstraintdihedralConstraints_
 

Static members

static constexpr double bondRelativeVariance = 0.01
 Relative bond distance variance 0.0x mean x% variance. Must fulfill 0 < x << 1.
 
static constexpr double angleRelativeVariance = 0.02
 Relative angle variance. More...
 
static constexpr double dihedralAbsoluteVariance = M_PI / 90
 Absolute dihedral angle variance in radians.
 
static constexpr ValueBounds angleClampBounds {0.0, M_PI}
 Defines clamping bounds on angles.
 
static const ValueBounds defaultDihedralBounds
 Defines the range (-π,π] using std::nextafter (not constexpr)
 
static double modelDistance (AtomIndex i, AtomIndex j, const PrivateGraph &graph)
 Models the equilibrium distance between two bonded atoms. More...
 
static double modelDistance (const BondIndex &bond, const PrivateGraph &graph)
 
static std::vector< BondIndexcycleConsistingOfExactly (const std::vector< AtomIndex > &atoms, const PrivateGraph &graph)
 Tries to find a cycle that consists of a particular set of atoms. More...
 
static boost::optional
< ValueBounds
coneAngle (const std::vector< AtomIndex > &baseConstituents, const ValueBounds &coneHeightBounds, const Graph &graph)
 Tries to calculate the cone angle for a possibly haptic ligand site. More...
 
static double spiroCrossAngle (double alpha, double beta)
 Calculates the cross angle between opposite cycle atoms in a spirocenter. More...
 
static ValueBounds siteDistanceFromCenter (const std::vector< AtomIndex > &siteAtomList, AtomIndex placement, const Graph &graph)
 Calculates bounds on the distance of a (possibly haptic) ligand site to a central index. More...
 
static ValueBounds makeBoundsFromCentralValue (double centralValue, double absoluteVariance)
 Yields central value plus/minus some absolute variance bounds. More...
 
static double smallestCycleDistortionMultiplier (const AtomIndex i, const Cycles &cycles)
 Idealization strictness loosening multiplier due to small cycles. More...
 
static BoundsMatrix makePairwiseBounds (unsigned N, const BoundsMapType< 2 > &fixedPositionBounds, const BoundsMapType< 2 > &bondBounds, const BoundsMapType< 3 > &angleBounds, const BoundsMapType< 4 > &dihedralBounds)
 Relative bond distance variance 0.0x mean x% variance. Must fulfill 0 < x << 1.
 
static double siteCentralAngle (AtomIndex placement, const Shapes::Shape &shape, const RankingInformation &ranking, const AtomStereopermutator::ShapeMap &shapeVertexMap, const std::pair< SiteIndex, SiteIndex > &sites, const PrivateGraph &inner)
 Determines the central value of the angle between AtomStereopermutator sites. More...
 
static double siteCentralAngle (const AtomStereopermutator &permutator, const std::pair< SiteIndex, SiteIndex > &sites, const PrivateGraph &inner)
 
static ValueBounds modelSiteAngleBounds (const AtomStereopermutator &permutator, const std::pair< SiteIndex, SiteIndex > &sites, const double looseningMultiplier, const PrivateGraph &inner)
 Models bounds on the angle between AtomStereopermutator sites. More...
 
static ChiralConstraint makeChiralConstraint (const AtomStereopermutator::MinimalChiralConstraint &minimalConstraint, const AtomStereopermutator &permutator, const double looseningMultiplier)
 Creates a volume-specific chiral constraint from a list of. More...
 
static ValueBounds clamp (ValueBounds bounds, const ValueBounds &clampBounds)
 Analogous to C++17's clamp, reduces ValueBounds to the maximal extent specified by some clamping bounds. More...
 
static void checkFixedPositionsPreconditions (const Molecule &molecule, const Configuration &configuration)
 Ensures that the preconditions for fixed positions in Distance Geometry detailed in the Configuration object are met. More...
 

Detailed Description

Class performing spatial modeling of molecules.

Keeps a record of the internal coordinate bounds that a molecular graph is interpreted as and generates a list of atom-pairwise bounds from which classes performing smoothing can be constructed.

Member Typedef Documentation

Type used to represent atom-pairwise distance bounds constructed from bounds on internal coordinates.

Diagonal entries are undefined. Strict lower triangle contains lower bounds, strict upper triangle contains upper bounds.

Constructor & Destructor Documentation

Scine::Molassembler::DistanceGeometry::SpatialModel::SpatialModel ( const Molecule molecule,
const Configuration configuration 
)

Model a molecule into internal coordinate bounds stored internally.

Main space where modeling is performed. This is where a molecule's graph and list of stereopermutators are transformed into explicit bounds on 1-2 (bonds), 1-3 (angles) and 1-4 (dihedral) internal coordinates. This determines which conformations are accessible.

Complexity At least \(O(P_2 + P_3 + P_4)\) where \(P_i\) is the number of distinct paths of length \(i\) in the graph. That should scale at least linearly in the number of vertices.

Parameters
moleculeThe molecule that is to be modeled. This may not contain stereopermutators with zero assignments or unassigned stereopermutators.
configurationThe Distance Geometry configuration object. Relevant for this stage of the process are the loosening multiplier and fixed positions, if set.

Member Function Documentation

void Scine::Molassembler::DistanceGeometry::SpatialModel::addAtomStereopermutatorInformation ( const AtomStereopermutator permutator,
const PrivateGraph graph,
double  looseningMultiplier,
const std::unordered_map< AtomIndex, Utils::Position > &  fixedAngstromPositions,
bool  forceChiralConstraintEmission 
)

Adds angle information to the internal coordinate bounds and collects chiral constraints.

Complexity \(O(S^2)\) where \(S\) is the size of the modeled shape

Parameters
permutatorThe AtomStereopermutator to collect information from
graphThe graph context
looseningMultiplierA loosening factor for the overall model
fixedAngstromPositionsA mapping between atom indices and fixed spatial positions
forceChiralConstraintEmissionIf set, the stereopermutator emits chiral constraints even if the stereopermutator does not have chiral state, i.e. numAssignments() <= 1, as long as it is assigned.
void Scine::Molassembler::DistanceGeometry::SpatialModel::addBondStereopermutatorInformation ( const BondStereopermutator permutator,
const AtomStereopermutator stereopermutatorA,
const AtomStereopermutator stereopermutatorB,
double  looseningMultiplier,
const std::unordered_map< AtomIndex, Utils::Position > &  fixedAngstromPositions 
)

Adds dihedral information to the internal coordinate bounds and collects dihedral constraints.

Complexity \(O(S^2)\) where \(S\) is the size of the larger modeled shape

Parameters
permutatorThe BondStereopermutator to collect information from
stereopermutatorAOne AtomStereopermutator constituting the BondStereopermutator
stereopermutatorBThe other AtomStereopermutator constituting the BondStereopermutator
looseningMultiplierA loosening factor for the overall model
fixedAngstromPositionsA mapping between atom indices and fixed spatial positions
void Scine::Molassembler::DistanceGeometry::SpatialModel::addDefaultDihedrals_ ( )
private

Adds [0, π] default dihedrals to the model (of connected sequences)

Adds [0, π] default dihedrals to the model. Avoids treating 1-4 pairs as nonbonded.

static void Scine::Molassembler::DistanceGeometry::SpatialModel::checkFixedPositionsPreconditions ( const Molecule molecule,
const Configuration configuration 
)
static

Ensures that the preconditions for fixed positions in Distance Geometry detailed in the Configuration object are met.

Parameters
moleculeThe molecule for which the Configuration is valid
configurationThe Configuration object detailing fixed positions
Exceptions
std::runtime_errorIf the preconditions are not met
static ValueBounds Scine::Molassembler::DistanceGeometry::SpatialModel::clamp ( ValueBounds  bounds,
const ValueBounds clampBounds 
)
static

Analogous to C++17's clamp, reduces ValueBounds to the maximal extent specified by some clamping bounds.

Complexity \(\Theta(1)\)

Parameters
boundsThe bounds to be clamped
clampBoundsThe range in which bounds may exist
Returns
New bounds that are within the clamp bounds
static boost::optional<ValueBounds> Scine::Molassembler::DistanceGeometry::SpatialModel::coneAngle ( const std::vector< AtomIndex > &  baseConstituents,
const ValueBounds coneHeightBounds,
const Graph graph 
)
static

Tries to calculate the cone angle for a possibly haptic ligand site.

Complexity \(O(C)\) where \(C\) is the number of cycles in the molecule

Parameters
baseConstituentsAtom indices constituting the ligand site
coneHeightBoundsBounds on the distance of the ligand site to the central atom
graphThe molecular graph to model
Returns
If calculable, bounds on the cone angle spanned by the ligand site
static std::vector<BondIndex> Scine::Molassembler::DistanceGeometry::SpatialModel::cycleConsistingOfExactly ( const std::vector< AtomIndex > &  atoms,
const PrivateGraph graph 
)
static

Tries to find a cycle that consists of a particular set of atoms.

Complexity \(\Theta(C)\) where \(C\) is the number of cycles in the molecule

Parameters
atomsThe exact set of atoms that constitute the cycle.
graphThe graph context
Returns
A vector of bond indices of exactly the matching cycle if found. An empty vector is returned otherwise.
std::string Scine::Molassembler::DistanceGeometry::SpatialModel::dumpGraphviz ( ) const

Generates a string graphviz representation of the modeled molecule.

The graph contains basic connectivity, stereopermutator information and internal coordinate bounds.

Complexity \(\Theta(N)\)

Returns
A string that can be converted into an image using graphviz of the molecule being modeled.
std::vector<ChiralConstraint> Scine::Molassembler::DistanceGeometry::SpatialModel::getChiralConstraints ( ) const

Yields all collected chiral constraints.

Complexity \(\Theta(1)\)

std::vector<DihedralConstraint> Scine::Molassembler::DistanceGeometry::SpatialModel::getDihedralConstraints ( ) const

Yields all collected dihedral constraints.

Complexity \(\Theta(1)\)

void Scine::Molassembler::DistanceGeometry::SpatialModel::instantiateMissingAtomStereopermutators_ ( Random::Engine engine)
private

It is permissible in some circumstances not to have AtomStereopermutators even on non-terminal atoms. These have to be re-added in order for us to be able to model everywhere.

static ValueBounds Scine::Molassembler::DistanceGeometry::SpatialModel::makeBoundsFromCentralValue ( double  centralValue,
double  absoluteVariance 
)
static

Yields central value plus/minus some absolute variance bounds.

Complexity \(\Theta(1)\)

Parameters
centralValueThe central value
absoluteVarianceThe value to subtract and add to yield new bounds
Returns
bounds with median central value and width 2 * absoluteVariance
static ChiralConstraint Scine::Molassembler::DistanceGeometry::SpatialModel::makeChiralConstraint ( const AtomStereopermutator::MinimalChiralConstraint minimalConstraint,
const AtomStereopermutator permutator,
const double  looseningMultiplier 
)
static

Creates a volume-specific chiral constraint from a list of.

Complexity \(\Theta(1)\)

Parameters
minimalConstraintSite index sequence defining a chiral constraint
permutatorThe AtomStereopermutator from which the minimalConstraint has been generated
looseningMultiplierDescribes the spatial modeling loosening of bounds
Returns
A volume constraint on the site index sequence in a final conformation
BoundsMatrix Scine::Molassembler::DistanceGeometry::SpatialModel::makePairwiseBounds ( ) const

Generate unsmoothed atom-pairwise distance bounds matrix.

Generates a matrix of atom-pairwise distance bounds from the internal coordinate bounds and fixed positions from which this was constructed.

Complexity \(O(P_2 + P_3 + P_4)\) where \(P_i\) is the number of distinct paths of length \(i\) in the graph. That should scale at least linearly in the number of vertices.

Returns
A matrix of atom-pairwise distance bounds. Lower bounds are in the strict lower triangle of the matrix, upper bounds in the strict upper triangle.
static double Scine::Molassembler::DistanceGeometry::SpatialModel::modelDistance ( AtomIndex  i,
AtomIndex  j,
const PrivateGraph graph 
)
static

Models the equilibrium distance between two bonded atoms.

Complexity \(\Theta(1)\)

Precondition
i and j are bonded
static double Scine::Molassembler::DistanceGeometry::SpatialModel::modelDistance ( const BondIndex bond,
const PrivateGraph graph 
)
static

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

bool Scine::Molassembler::DistanceGeometry::SpatialModel::modelPartiallyFixedBond ( const BondStereopermutator permutator,
const std::pair< const AtomStereopermutator &, const AtomStereopermutator & > &  atomStereopermutators,
const std::unordered_map< AtomIndex, Utils::Position > &  fixedAngstromPositions 
)

Sets bond bounds to exact value bounds if unset.

Complexity \(\Theta(1)\)

Parameters
bondIndicesThe atom pair on which to place bond distance bounds
boundsThe distance bounds to enforce
Precondition
bondIndices must be ordered (i.e. front() < back())
static ValueBounds Scine::Molassembler::DistanceGeometry::SpatialModel::modelSiteAngleBounds ( const AtomStereopermutator permutator,
const std::pair< SiteIndex, SiteIndex > &  sites,
const double  looseningMultiplier,
const PrivateGraph inner 
)
static

Models bounds on the angle between AtomStereopermutator sites.

Complexity Varies. For most cases, \(\Omega(1)\)

Parameters
permutatorThe permutator being modeled
sitesThe two sites between which the angle is to be determined
looseningMultiplierLoosening of that particular atom
innerGraph instance being modeled
Returns
Bounds on the angle between two AtomStereopermutator sites
void Scine::Molassembler::DistanceGeometry::SpatialModel::modelSpirocenters_ ( const FixedPositionsMapType fixedAngstromPositions)
private

Adds spirocenter modelling to the data set.

The cross-angle in a spiro-atom is expliticly modeled if the disjoint cycles are both small, i.e. of size \(\le 5\). Only then is the additional variance when calculating cross-angle bounds particularly harmful for the overall 3D structure. We define:

\begin{align} \vec{a} &= R_z\left(\frac{\alpha}{2}\right)\begin{pmatrix} 1\\ 0\\ 0 \end{pmatrix} = \begin{pmatrix} \cos \frac{\alpha}{2}\\ \sin \frac{\alpha}{2}\\ 0 \end{pmatrix}\\ \vec{b} &= R_y\left(\frac{\beta}{2}\right)\begin{pmatrix} -1\\ 0\\ 0 \end{pmatrix} = \begin{pmatrix} - \cos \frac{\beta}{2}\\ 0\\ \sin \frac{\beta}{2} \end{pmatrix} \end{align}

Where \(R_y\) and \(R_z\) are the standard three-dimensional rotation matrices, \(\alpha\) is the cycle internal angle for the first cycle at the spiro-atom, and \(\beta\) is the cycle internal angle for the second cycle at the spiro-atom. The angle between these vectors \(\phi\) is then:

\begin{align} \vec{a}\vec{b} &= |\vec{a}||\vec{b}|\cos\phi\\ \phi &= \arccos \left\{ \left(\cos \frac{\alpha}{2} \right) \left(- \cos \frac{\beta}{2} \right) \right\} \end{align}

void Scine::Molassembler::DistanceGeometry::SpatialModel::setAngleBoundsIfEmpty ( std::array< AtomIndex, 3 >  angleIndices,
ValueBounds  bounds 
)

Adds angle bounds to the model if unset.

Complexity \(\Theta(1)\)

Parameters
angleIndicesThe atom index sequence specifying an angle
boundsThe angle bounds to enforce
Precondition
angleIndices must be ordered (i.e. front() < back())
Note
Passed angle bounds are clamped against the angleClampBounds
void Scine::Molassembler::DistanceGeometry::SpatialModel::setBondBoundsIfEmpty ( std::array< AtomIndex, 2 >  bondIndices,
ValueBounds  bounds 
)

Sets bond bounds to exact value bounds if unset.

Complexity \(\Theta(1)\)

Parameters
bondIndicesThe atom pair on which to place bond distance bounds
boundsThe distance bounds to enforce
Precondition
bondIndices must be ordered (i.e. front() < back())
void Scine::Molassembler::DistanceGeometry::SpatialModel::setDihedralBoundsIfEmpty ( std::array< AtomIndex, 4 >  dihedralIndices,
ValueBounds  bounds 
)

Adds the dihedral bounds to the model if unset.

Complexity \(\Theta(1)\)

Parameters
dihedralIndicesThe atom index sequence specifying a dihedral
boundsThe dihedral bounds to enforce
Precondition
dihedralIndices must be ordered (i.e. front() < back())
static double Scine::Molassembler::DistanceGeometry::SpatialModel::siteCentralAngle ( AtomIndex  placement,
const Shapes::Shape shape,
const RankingInformation ranking,
const AtomStereopermutator::ShapeMap shapeVertexMap,
const std::pair< SiteIndex, SiteIndex > &  sites,
const PrivateGraph inner 
)
static

Determines the central value of the angle between AtomStereopermutator sites.

Complexity Varies. For most cases, \(\Omega(1)\)

Parameters
placementThe atom index of the central index of the angle
shapeThe local shape at placement
rankingThe ranking of substituents at placement
shapeVertexMapThe mapping from site indices to shape vertices
sitesThe two sites between which the angle is to be determined
innerGraph instance being modeled
Returns
The central value of the angle between the sites
static double Scine::Molassembler::DistanceGeometry::SpatialModel::siteCentralAngle ( const AtomStereopermutator permutator,
const std::pair< SiteIndex, SiteIndex > &  sites,
const PrivateGraph inner 
)
static

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

static ValueBounds Scine::Molassembler::DistanceGeometry::SpatialModel::siteDistanceFromCenter ( const std::vector< AtomIndex > &  siteAtomList,
AtomIndex  placement,
const Graph graph 
)
static

Calculates bounds on the distance of a (possibly haptic) ligand site to a central index.

Complexity \(\Theta(A)\) where \(A\) is the number of atoms comprising the site

Parameters
siteAtomListAtom indices constituting the ligand site
placementThe central index to which the ligand is bound
graphThe molecular graph to model
Returns
Bounds on the distance of the ligand site to the central index
static double Scine::Molassembler::DistanceGeometry::SpatialModel::smallestCycleDistortionMultiplier ( const AtomIndex  i,
const Cycles cycles 
)
static

Idealization strictness loosening multiplier due to small cycles.

Returns
A multiplier intended for the absolute angle variance for a particular atom. If that atom is part of a cycle of size < 6, the multiplier is > 1.
static double Scine::Molassembler::DistanceGeometry::SpatialModel::spiroCrossAngle ( double  alpha,
double  beta 
)
static

Calculates the cross angle between opposite cycle atoms in a spirocenter.

Complexity \(\Theta(1)\)

Parameters
alphaThe left cycle-internal angle
betaThe right cycle-internal angle
Returns
The cross angle between cycle atoms in the spirocenter
void Scine::Molassembler::DistanceGeometry::SpatialModel::writeGraphviz ( const std::string &  filename) const

Writes a graphviz representation of a modeled molecule to a file.

The graph contains basic connectivity, stereopermutator information and internal coordinate bounds.

Complexity \(\Theta(N)\)

Parameters
filenameThe filename to which to write the graphviz representation

Field Documentation

constexpr double Scine::Molassembler::DistanceGeometry::SpatialModel::angleRelativeVariance = 0.02
static

Relative angle variance.

0.0x mean x& variance. Must fulfill 0 < x << 1


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