Molassembler  3.0.0
Molecule graph and conformer library
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
SpatialModel.h
Go to the documentation of this file.
1 
18 #ifndef INCLUDE_MOLASSEMBLER_DISTANCE_GEOMETRY_SPATIAL_MODEL_H
19 #define INCLUDE_MOLASSEMBLER_DISTANCE_GEOMETRY_SPATIAL_MODEL_H
20 
22 #include "Molassembler/Molecule.h"
23 #include "Molassembler/Graph.h"
27 
28 #include <unordered_map>
29 
30 namespace Scine {
31 namespace Molassembler {
32 namespace DistanceGeometry {
33 
40 class SpatialModel {
41 public:
45  struct ModelGraphWriter;
46 
48  template<unsigned size>
49  using BoundsMapType = std::unordered_map<
50  std::array<AtomIndex, size>,
53  std::array<AtomIndex, size>
54  >
55  >;
56 
58  using FixedPositionsMapType = std::unordered_map<AtomIndex, Utils::Position>;
59 
67  using BoundsMatrix = Eigen::MatrixXd;
68 
70  inline explicit BoundsMatrixHelper(AtomIndex size) : matrix(size, size) { matrix.setZero(); }
71 
72  void add(AtomIndex i, AtomIndex j, const ValueBounds& bounds);
73  void addMap(const BoundsMapType<2>& boundsMap);
74 
75  ValueBounds get(AtomIndex i, AtomIndex j) const;
76 
77  Eigen::MatrixXd matrix;
78  };
80 
83 
87  static constexpr double bondRelativeVariance = 0.01;
92  static constexpr double angleRelativeVariance = 0.02;
94  static constexpr double dihedralAbsoluteVariance = M_PI / 90; // ~ 2°
95 
97  static constexpr ValueBounds angleClampBounds {0.0, M_PI};
100 
101 /* Static functions */
107  static double modelDistance(AtomIndex i, AtomIndex j, const PrivateGraph& graph);
109  static double modelDistance(const BondIndex& bond, const PrivateGraph& graph);
110 
122  static std::vector<BondIndex> cycleConsistingOfExactly(
123  const std::vector<AtomIndex>& atoms,
124  const PrivateGraph& graph
125  );
126 
139  static boost::optional<ValueBounds> coneAngle(
140  const std::vector<AtomIndex>& baseConstituents,
141  const ValueBounds& coneHeightBounds,
142  const PrivateGraph& graph
143  );
144 
154  static double spiroCrossAngle(double alpha, double beta);
155 
169  const std::vector<AtomIndex>& siteAtomList,
170  AtomIndex placement,
171  const PrivateGraph& graph
172  );
173 
184  double centralValue,
185  double absoluteVariance
186  );
187 
196  static double smallestCycleDistortionMultiplier(
197  AtomIndex i,
198  const Cycles& cycles
199  );
200 
201  /*
202  * @brief Generates a matrix of atom-pairwise distance bounds
203  *
204  * @complexity{@math{O(P_2 + P_3 + P_4)} where @math{P_i} is the number of
205  * distinct paths of length @math{i} in the graph. That should scale at least
206  * linearly in the number of vertices.}
207  *
208  * @param N size of the system
209  * @param fixedPositionbounds Distances to enforce due to fixed positions
210  * @param bondBounds Distances to enforce for bonds
211  * @param angleBounds Angle bounds to enforce on atom index triples
212  * @param dihedralBounds Dihedral bounds to enforce on atom index quadruplet
213  *
214  * @return A matrix of atom-pairwise distance bounds. Lower bounds are in the
215  * strict lower triangle of the matrix, upper bounds in the strict upper
216  * triangle.
217  */
219  unsigned N,
220  const BoundsMapType<2>& fixedPositionBounds,
221  const BoundsMapType<2>& bondBounds,
222  const BoundsMapType<3>& angleBounds,
223  const BoundsMapType<4>& dihedralBounds
224  );
225 
240  static double siteCentralAngle(
241  AtomIndex placement,
242  const Shapes::Shape& shape,
243  const RankingInformation& ranking,
244  const AtomStereopermutator::ShapeMap& shapeVertexMap,
245  const std::pair<SiteIndex, SiteIndex>& sites,
246  const PrivateGraph& inner
247  );
248 
250  static double siteCentralAngle(
251  const AtomStereopermutator& permutator,
252  const std::pair<SiteIndex, SiteIndex>& sites,
253  const PrivateGraph& inner
254  );
255 
268  const AtomStereopermutator& permutator,
269  const Stereopermutators::LocalSpatialModel& localModel,
270  const std::pair<SiteIndex, SiteIndex>& sites,
271  double looseningMultiplier,
272  const PrivateGraph& inner
273  );
274 
289  const AtomStereopermutator::MinimalChiralConstraint& minimalConstraint,
290  const AtomStereopermutator& permutator,
291  const Stereopermutators::LocalSpatialModel& localModel,
292  double looseningMultiplier
293  );
294 
305  static ValueBounds clamp(
306  ValueBounds bounds,
307  const ValueBounds& clampBounds
308  );
309 
320  const Molecule& molecule,
321  const Configuration& configuration
322  );
324 
327 
345  SpatialModel(
346  const Molecule& molecule,
347  const Configuration& configuration
348  );
350 
353 
362  std::array<AtomIndex, 2> bondIndices,
363  ValueBounds bounds
364  );
365 
377  std::array<AtomIndex, 3> angleIndices,
378  ValueBounds bounds
379  );
380 
391  std::array<AtomIndex, 4> dihedralIndices,
392  ValueBounds bounds
393  );
394 
411  const AtomStereopermutator& permutator,
412  const PrivateGraph& graph,
413  double looseningMultiplier,
414  const std::unordered_map<AtomIndex, Utils::Position>& fixedAngstromPositions,
415  bool forceChiralConstraintEmission
416  );
417 
432  const BondStereopermutator& permutator,
433  const AtomStereopermutator& stereopermutatorA,
434  const AtomStereopermutator& stereopermutatorB,
435  double looseningMultiplier,
436  const std::unordered_map<AtomIndex, Utils::Position>& fixedAngstromPositions
437  );
438 
440  const BondStereopermutator& permutator,
441  const std::pair<const AtomStereopermutator&, const AtomStereopermutator&>& atomPermutators,
442  const std::unordered_map<AtomIndex, Utils::Position>& fixedAngstromPositions
443  );
445 
448 
452  std::vector<ChiralConstraint> getChiralConstraints() const;
453 
458  std::vector<DihedralConstraint> getDihedralConstraints() const;
459 
475 
486  std::string dumpGraphviz() const;
487 
497  void writeGraphviz(const std::string& filename) const;
499 
500 private:
501  // Molecule closure
502  const Molecule& molecule_;
503 
505  std::unordered_map<AtomIndex, Stereopermutators::LocalSpatialModel> localModels_;
506 
515 
517  std::vector<ChiralConstraint> chiralConstraints_;
518  std::vector<DihedralConstraint> dihedralConstraints_;
519 
523  void addDefaultAngles_();
524 
531  void addDefaultDihedrals_();
532 
539 
541  void modelBondDistances_(
542  const FixedPositionsMapType& fixedAngstromPositions,
543  double looseningFactor
544  );
545 
547  void modelFlatCycles_(
548  const FixedPositionsMapType& fixedAngstromPositions,
549  double looseningFactor
550  );
551 
593  void modelSpirocenters_(
594  const FixedPositionsMapType& fixedAngstromPositions
595  );
597 };
598 
599 } // namespace DistanceGeometry
600 } // namespace Molassembler
601 } // namespace Scine
602 
603 #endif
std::array< boost::optional< SiteIndex >, 4 > MinimalChiralConstraint
Site index sequence defining a chiral constraint. If a site index is None, then it denotes the positi...
Definition: AtomStereopermutator.h:108
Models a molecule as a graph (connectivity of atoms) and a list of stereopermutators.
Definition: Molecule.h:77
SpatialModel(const Molecule &molecule, const Configuration &configuration)
Model a molecule into internal coordinate bounds stored internally.
Ranking data of substituents around a central vertex.
Definition: RankingInformation.h:23
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.
static constexpr double angleRelativeVariance
Relative angle variance.
Definition: SpatialModel.h:92
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.
Drives a PRNG.
Definition: Prng.h:24
A configuration object for distance geometry runs with sane defaults.
Definition: Conformers.h:73
unsigned size(Shape shape)
Fetch the number of vertices of a shape.
BoundsMapType< 4 > dihedralBounds_
Dihedral value bounds on index quadruplets.
Definition: SpatialModel.h:514
ShapeResult shape(const PositionCollection &normalizedPositions, Shape shape)
Forwarding function to calculate the continuous shape measure.
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.
static constexpr ValueBounds angleClampBounds
Defines clamping bounds on angles.
Definition: SpatialModel.h:97
Class performing spatial modeling of molecules.
Definition: SpatialModel.h:40
BoundsMapType< 2 > constraints_
Constraints by fixed positions.
Definition: SpatialModel.h:508
static double smallestCycleDistortionMultiplier(AtomIndex i, const Cycles &cycles)
Idealization strictness loosening multiplier due to small cycles.
Wrapper class to make working with RDL in C++ more pleasant.
Definition: Cycles.h:44
Decide which stereopermutations are feasible.
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...
Handles the steric permutation of substituents of a non-terminal central atom.
Definition: AtomStereopermutator.h:79
std::unordered_map< AtomIndex, Utils::Position > FixedPositionsMapType
Type used to store fixed positions in angstrom.
Definition: SpatialModel.h:58
static const ValueBounds defaultDihedralBounds
Defines the range (-π,π] using std::nextafter (not constexpr)
Definition: SpatialModel.h:99
static void checkFixedPositionsPreconditions(const Molecule &molecule, const Configuration &configuration)
Ensures that the preconditions for fixed positions in Distance Geometry detailed in the Configuration...
Handle arrangements of substituents at corners of an atom-centered shape.
std::string dumpGraphviz() const
Generates a string graphviz representation of the modeled molecule.
void setDihedralBoundsIfEmpty(std::array< AtomIndex, 4 > dihedralIndices, ValueBounds bounds)
Adds the dihedral bounds to the model if unset.
Molecule class interface.
static constexpr double dihedralAbsoluteVariance
Absolute dihedral angle variance in radians.
Definition: SpatialModel.h:94
Handles specific relative arrangements of two atom stereopermutators joined by a bond.
Definition: BondStereopermutator.h:46
void modelSpirocenters_(const FixedPositionsMapType &fixedAngstromPositions)
Adds spirocenter modelling to the data set.
static double modelDistance(AtomIndex i, AtomIndex j, const PrivateGraph &graph)
Models the equilibrium distance between two bonded atoms.
std::vector< DihedralConstraint > getDihedralConstraints() const
Yields all collected dihedral constraints.
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.
std::size_t AtomIndex
Unsigned integer atom index type. Used to refer to particular atoms.
Definition: Types.h:51
std::vector< std::vector< AtomIndex >> sites(const PrivateGraph &graph, AtomIndex placement, const std::vector< AtomIndex > &excludeAdjacents={})
Differentiate adjacent vertices of a central index into sites.
std::vector< ChiralConstraint > chiralConstraints_
Chiral constraints.
Definition: SpatialModel.h:517
Library internal graph class wrapping BGL types.
Definition: PrivateGraph.h:26
BoundsMapType< 3 > angleBounds_
Angle value bounds on index triples.
Definition: SpatialModel.h:512
std::unordered_map< std::array< AtomIndex, size >, ValueBounds, boost::hash< std::array< AtomIndex, size > > > BoundsMapType
ValueBounds container for various internal coordinates.
Definition: SpatialModel.h:55
static ValueBounds makeBoundsFromCentralValue(double centralValue, double absoluteVariance)
Yields central value plus/minus some absolute variance bounds.
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.
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.
static ValueBounds clamp(ValueBounds bounds, const ValueBounds &clampBounds)
Analogous to C++17&#39;s clamp, reduces ValueBounds to the maximal extent specified by some clamping boun...
BoundsMatrix makePairwiseBounds() const
Generate unsmoothed atom-pairwise distance bounds matrix.
void addDefaultDihedrals_()
Adds [0, π] default dihedrals to the model (of connected sequences)
std::vector< ChiralConstraint > getChiralConstraints() const
Yields all collected chiral constraints.
Type used to refer to particular bonds. Orders first &lt; second.
Definition: Types.h:54
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.
Interface class for the molecular graph.
void setAngleBoundsIfEmpty(std::array< AtomIndex, 3 > angleIndices, ValueBounds bounds)
Adds angle bounds to the model if unset.
BoundsMapType< 2 > bondBounds_
Bond distance value bounds on index pairs.
Definition: SpatialModel.h:510
Data class for bounded values.
Definition: ValueBounds.h:23
WideHashType hash(AtomEnvironmentComponents bitmask, Utils::ElementType elementType, const std::vector< BondInformation > &sortedBonds, const boost::optional< Shapes::Shape > &shapeOptional, const boost::optional< unsigned > &assignedOptional)
Convolutes the atom&#39;s element type and bonds into an unsigned integer.
Shape
Enumeration of all contained symmetry names.
Definition: Shapes.h:28
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 instantiateMissingAtomStereopermutators_(Random::Engine &engine)
static double spiroCrossAngle(double alpha, double beta)
Calculates the cross angle between opposite cycle atoms in a spirocenter.
Data struct representing a chiral constraint.
Definition: DistanceGeometry.h:32
void addDefaultAngles_()
Adds [0, π] default angle bounds for all bonded atom triples.
Owning class storing all stereopermutators in a molecule.
static constexpr double bondRelativeVariance
Relative bond distance variance 0.0x mean x% variance. Must fulfill 0 &lt; x &lt;&lt; 1.
Definition: SpatialModel.h:87
Class storing atom-pairwise distance bounds.
void writeGraphviz(const std::string &filename) const
Writes a graphviz representation of a modeled molecule to a file.
std::unordered_map< AtomIndex, Stereopermutators::LocalSpatialModel > localModels_
Local models of atom stereopermutators.
Definition: SpatialModel.h:505
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...
void setBondBoundsIfEmpty(std::array< AtomIndex, 2 > bondIndices, ValueBounds bounds)
Sets bond bounds to exact value bounds if unset.
Eigen::MatrixXd BoundsMatrix
Type used to represent atom-pairwise distance bounds constructed from bounds on internal coordinates...
Definition: SpatialModel.h:67