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 & > &atomPermutators, const std::unordered_map< AtomIndex, Utils::Position > &fixedAngstromPositions) |
| Sets bond bounds to exact value bounds if unset. More... | |
Information | |
| std::vector< ChiralConstraint > | getChiralConstraints () const |
| Yields all collected chiral constraints. More... | |
| std::vector< DihedralConstraint > | getDihedralConstraints () 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 Molecule & | molecule_ |
|
std::unordered_map< AtomIndex, Stereopermutators::LocalSpatialModel > | localModels_ |
| Local models of atom stereopermutators. | |
| 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< ChiralConstraint > | chiralConstraints_ |
| Chiral constraints. | |
| std::vector< DihedralConstraint > | dihedralConstraints_ |
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< BondIndex > | cycleConsistingOfExactly (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 PrivateGraph &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 PrivateGraph &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 (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 Stereopermutators::LocalSpatialModel &localModel, const std::pair< SiteIndex, SiteIndex > &sites, 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 Stereopermutators::LocalSpatialModel &localModel, 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... | |
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.
| using Scine::Molassembler::DistanceGeometry::SpatialModel::BoundsMatrix = Eigen::MatrixXd |
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.
| 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.
| molecule | The molecule that is to be modeled. This may not contain stereopermutators with zero assignments or unassigned stereopermutators. |
| configuration | The Distance Geometry configuration object. Relevant for this stage of the process are the loosening multiplier and fixed positions, if set. |
| 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
| permutator | The AtomStereopermutator to collect information from |
| graph | The graph context |
| looseningMultiplier | A loosening factor for the overall model |
| fixedAngstromPositions | A mapping between atom indices and fixed spatial positions |
| forceChiralConstraintEmission | If 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
| permutator | The BondStereopermutator to collect information from |
| stereopermutatorA | One AtomStereopermutator constituting the BondStereopermutator |
| stereopermutatorB | The other AtomStereopermutator constituting the BondStereopermutator |
| looseningMultiplier | A loosening factor for the overall model |
| fixedAngstromPositions | A mapping between atom indices and fixed spatial positions |
|
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 |
Ensures that the preconditions for fixed positions in Distance Geometry detailed in the Configuration object are met.
| molecule | The molecule for which the Configuration is valid |
| configuration | The Configuration object detailing fixed positions |
| std::runtime_error | If the preconditions are not met |
|
static |
Analogous to C++17's clamp, reduces ValueBounds to the maximal extent specified by some clamping bounds.
Complexity \(\Theta(1)\)
| bounds | The bounds to be clamped |
| clampBounds | The range in which bounds may exist |
|
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
| baseConstituents | Atom indices constituting the ligand site |
| coneHeightBounds | Bounds on the distance of the ligand site to the central atom |
| graph | The molecular graph to model |
|
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
| atoms | The exact set of atoms that constitute the cycle. |
| graph | The graph context |
| 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)\)
| 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)\)
|
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 |
Yields central value plus/minus some absolute variance bounds.
Complexity \(\Theta(1)\)
| centralValue | The central value |
| absoluteVariance | The value to subtract and add to yield new bounds |
|
static |
Creates a volume-specific chiral constraint from a list of.
Complexity \(\Theta(1)\)
| minimalConstraint | Site index sequence defining a chiral constraint |
| permutator | The AtomStereopermutator from which the minimalConstraint has been generated |
| looseningMultiplier | Describes the spatial modeling loosening of bounds |
| 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.
|
static |
Models the equilibrium distance between two bonded atoms.
Complexity \(\Theta(1)\)
i and j are bonded
|
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 & > & | atomPermutators, | ||
| const std::unordered_map< AtomIndex, Utils::Position > & | fixedAngstromPositions | ||
| ) |
Sets bond bounds to exact value bounds if unset.
Complexity \(\Theta(1)\)
| bondIndices | The atom pair on which to place bond distance bounds |
| bounds | The distance bounds to enforce |
bondIndices must be ordered (i.e. front() < back())
|
static |
Models bounds on the angle between AtomStereopermutator sites.
Complexity Varies. For most cases, \(\Omega(1)\)
| permutator | The permutator being modeled |
| sites | The two sites between which the angle is to be determined |
| looseningMultiplier | Loosening of that particular atom |
| inner | Graph instance being modeled |
|
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)\)
| angleIndices | The atom index sequence specifying an angle |
| bounds | The angle bounds to enforce |
angleIndices must be ordered (i.e. front() < back())| 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)\)
| bondIndices | The atom pair on which to place bond distance bounds |
| bounds | The distance bounds to enforce |
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)\)
| dihedralIndices | The atom index sequence specifying a dihedral |
| bounds | The dihedral bounds to enforce |
dihedralIndices must be ordered (i.e. front() < back())
|
static |
Determines the central value of the angle between AtomStereopermutator sites.
Complexity Varies. For most cases, \(\Omega(1)\)
| placement | The atom index of the central index of the angle |
| shape | The local shape at placement |
| ranking | The ranking of substituents at placement |
| shapeVertexMap | The mapping from site indices to shape vertices |
| sites | The two sites between which the angle is to be determined |
| inner | Graph instance being modeled |
|
static |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
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
| siteAtomList | Atom indices constituting the ligand site |
| placement | The central index to which the ligand is bound |
| graph | The molecular graph to model |
|
static |
Idealization strictness loosening multiplier due to small cycles.
|
static |
Calculates the cross angle between opposite cycle atoms in a spirocenter.
Complexity \(\Theta(1)\)
| alpha | The left cycle-internal angle |
| beta | The right cycle-internal angle |
| 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)\)
| filename | The filename to which to write the graphviz representation |
|
static |
Relative angle variance.
0.0x mean x& variance. Must fulfill 0 < x << 1