RDKit
Open-source cheminformatics and machine learning.
MolTransforms.h
Go to the documentation of this file.
1 //
2 // Copyright (C) 2003-2006 Rational Discovery LLC
3 //
4 // @@ All Rights Reserved @@
5 // This file is part of the RDKit.
6 // The contents are covered by the terms of the BSD license
7 // which is included in the file license.txt, found at the root
8 // of the RDKit source tree.
9 //
10 #ifndef _RD_MOLTRANSFORMS_H_
11 #define _RD_MOLTRANSFORMS_H_
12 
13 #include <Geometry/point.h>
14 #include <Numerics/SymmMatrix.h>
15 
16 #ifdef RDK_HAS_EIGEN3
17 #include <Eigen/Dense>
18 #endif
19 
20 namespace RDKit {
21 class ROMol;
22 class Atom;
23 class Conformer;
24 }
25 
26 namespace RDGeom {
27 class Transform3D;
28 }
29 
30 namespace MolTransforms {
33 
34 //! Compute the centroid of a conformer
35 /*!
36  This is simple the average of the heavy atom locations in the conformer,
37  not attention is paid to hydrogens or the differences in atomic radii
38 
39  \param conf Conformer of interest
40  \param ignoreHs If true, ignore hydrogen atoms
41 */
43  bool ignoreHs = true);
44 
45 #ifdef RDK_HAS_EIGEN3
46 //! Compute principal axes and moments of inertia for a conformer
47 /*!
48  These values are calculated from the inertia tensor:
49  Iij = - sum_{s=1..N}(w_s * r_{si} * r_{sj}) i != j
50  Iii = sum_{s=1..N} sum_{j!=i} (w_s * r_{sj} * r_{sj})
51  where the coordinates are relative to the center of mass.
52 
53 
54  \param conf Conformer of interest
55  \param axes used to return the principal axes
56  \param moments used to return the principal moments
57  \param ignoreHs If true, ignore hydrogen atoms
58  \param force If true, the calculation will be carried out even if a
59  cached value is present
60  \param weights If present used to weight the atomic coordinates
61 
62  \returns whether or not the calculation was successful
63 */
64 bool computePrincipalAxesAndMoments(const RDKit::Conformer &conf,
65  Eigen::Matrix3d &axes,
66  Eigen::Vector3d &moments,
67  bool ignoreHs = false, bool force = false,
68  const std::vector<double> *weights = NULL);
69 //! Compute principal axes and moments from the gyration matrix of a conformer
70 /*!
71 
72  These values are calculated from the gyration matrix/tensor:
73  Iij = sum_{s=1..N}(w_s * r_{si} * r_{sj}) i != j
74  Iii = sum_{s=1..N} sum_{t!=s}(w_s * r_{si} * r_{ti})
75  where the coordinates are relative to the center of mass.
76 
77  \param conf Conformer of interest
78  \param axes used to return the principal axes
79  \param moments used to return the principal moments
80  \param ignoreHs If true, ignore hydrogen atoms
81  \param force If true, the calculation will be carried out even if a
82  cached value is present
83  \param weights If present used to weight the atomic coordinates
84 
85  \returns whether or not the calculation was successful
86 */
87 bool computePrincipalAxesAndMomentsFromGyrationMatrix(
88  const RDKit::Conformer &conf, Eigen::Matrix3d &axes,
89  Eigen::Vector3d &moments, bool ignoreHs = false, bool force = false,
90  const std::vector<double> *weights = NULL);
91 #endif
92 
93 //! Compute the transformation require to orient the conformation
94 //! along the principal axes about the center; i.e. center is made to coincide
95 // with the
96 //! origin, the largest princiapl axis with the x-axis, the next largest with
97 // the y-axis
98 //! and the smallest with the z-axis
99 /*!
100  If center is not specified the the centroid of the conformer will be used
101  \param conf Conformer of interest
102  \param center Center to be used for canonicalization, defaults to
103  the centroid of the
104  conformation
105  \param normalizeCovar Normalize the covariance matrix with the number of
106  atoms
107  \param ignoreHs Optinally ignore hydrogens
108 */
110  const RDKit::Conformer &conf, const RDGeom::Point3D *center = 0,
111  bool normalizeCovar = false, bool ignoreHs = true);
112 
113 //! Transform the conformation using the specified transformation
115  const RDGeom::Transform3D &trans);
116 
117 //! Canonicalize the orientation of a conformer so that its principal axes
118 //! around the specified center point coincide with the x, y, z axes
119 /*!
120  \param conf The conformer of interest
121  \param center Optional center point about which the principal axes are
122  computed
123  if not specified the centroid of the conformer will be
124  used
125  \param normalizeCovar Optionally normalize the covariance matrix by the number
126  of atoms
127  \param ignoreHs If true, ignore hydrogen atoms
128 
129 */
131  const RDGeom::Point3D *center = 0,
132  bool normalizeCovar = false, bool ignoreHs = true);
133 
134 //! Canonicalize all the conformations in a molecule
135 /*!
136  \param mol the molecule of interest
137  \param normalizeCovar Optionally normalize the covariance matrix by the number
138  of atoms
139  \param ignoreHs If true, ignore hydrogens
140 */
141 void canonicalizeMol(RDKit::ROMol &mol, bool normalizeCovar = false,
142  bool ignoreHs = true);
143 
144 //! Get the bond length between the specified atoms i, j
145 double getBondLength(const RDKit::Conformer &conf, unsigned int iAtomId,
146  unsigned int jAtomId);
147 
148 //! Set the bond length between the specified atoms i, j
149 //! (all atoms bonded to atom j are moved)
150 void setBondLength(RDKit::Conformer &conf, unsigned int iAtomId,
151  unsigned int jAtomId, double value);
152 
153 //! Get the angle in radians among the specified atoms i, j, k
154 double getAngleRad(const RDKit::Conformer &conf, unsigned int iAtomId,
155  unsigned int jAtomId, unsigned int kAtomId);
156 
157 //! Get the angle in degrees among the specified atoms i, j, k
158 inline double getAngleDeg(const RDKit::Conformer &conf, unsigned int iAtomId,
159  unsigned int jAtomId, unsigned int kAtomId) {
160  return (180. / M_PI * getAngleRad(conf, iAtomId, jAtomId, kAtomId));
161 }
162 
163 //! Set the angle in radians among the specified atoms i, j, k
164 //! (all atoms bonded to atom k are moved)
165 void setAngleRad(RDKit::Conformer &conf, unsigned int iAtomId,
166  unsigned int jAtomId, unsigned int kAtomId, double value);
167 
168 //! Set the angle in degrees among the specified atoms i, j, k
169 //! (all atoms bonded to atom k are moved)
170 inline void setAngleDeg(RDKit::Conformer &conf, unsigned int iAtomId,
171  unsigned int jAtomId, unsigned int kAtomId,
172  double value) {
173  setAngleRad(conf, iAtomId, jAtomId, kAtomId, value / 180. * M_PI);
174 }
175 
176 //! Get the dihedral angle in radians among the specified atoms i, j, k, l
177 double getDihedralRad(const RDKit::Conformer &conf, unsigned int iAtomId,
178  unsigned int jAtomId, unsigned int kAtomId,
179  unsigned int lAtomId);
180 
181 //! Get the dihedral angle in degrees among the specified atoms i, j, k, l
182 inline double getDihedralDeg(const RDKit::Conformer &conf, unsigned int iAtomId,
183  unsigned int jAtomId, unsigned int kAtomId,
184  unsigned int lAtomId) {
185  return (180. / M_PI *
186  getDihedralRad(conf, iAtomId, jAtomId, kAtomId, lAtomId));
187 }
188 
189 //! Set the dihedral angle in radians among the specified atoms i, j, k, l
190 //! (all atoms bonded to atom l are moved)
191 void setDihedralRad(RDKit::Conformer &conf, unsigned int iAtomId,
192  unsigned int jAtomId, unsigned int kAtomId,
193  unsigned int lAtomId, double value);
194 
195 //! Set the dihedral angle in degrees among the specified atoms i, j, k, l
196 //! (all atoms bonded to atom l are moved)
197 inline void setDihedralDeg(RDKit::Conformer &conf, unsigned int iAtomId,
198  unsigned int jAtomId, unsigned int kAtomId,
199  unsigned int lAtomId, double value) {
200  setDihedralRad(conf, iAtomId, jAtomId, kAtomId, lAtomId, value / 180. * M_PI);
201 }
202 }
203 #endif
double getDihedralDeg(const RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId, unsigned int lAtomId)
Get the dihedral angle in degrees among the specified atoms i, j, k, l.
void canonicalizeConformer(RDKit::Conformer &conf, const RDGeom::Point3D *center=0, bool normalizeCovar=false, bool ignoreHs=true)
void transformAtom(RDKit::Atom *atom, RDGeom::Transform3D &tform)
double getAngleRad(const RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId)
Get the angle in radians among the specified atoms i, j, k.
double getDihedralRad(const RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId, unsigned int lAtomId)
Get the dihedral angle in radians among the specified atoms i, j, k, l.
void canonicalizeMol(RDKit::ROMol &mol, bool normalizeCovar=false, bool ignoreHs=true)
Canonicalize all the conformations in a molecule.
void setDihedralDeg(RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId, unsigned int lAtomId, double value)
ROMol is a molecule class that is intended to have a fixed topology.
Definition: ROMol.h:106
void transformConformer(RDKit::Conformer &conf, const RDGeom::Transform3D &trans)
Transform the conformation using the specified transformation.
double getAngleDeg(const RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId)
Get the angle in degrees among the specified atoms i, j, k.
void setAngleDeg(RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId, double value)
void setAngleRad(RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId, double value)
Includes a bunch of functionality for handling Atom and Bond queries.
Definition: Atom.h:29
RDGeom::Point3D computeCentroid(const RDKit::Conformer &conf, bool ignoreHs=true)
Compute the centroid of a conformer.
double getBondLength(const RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId)
Get the bond length between the specified atoms i, j.
void transformMolsAtoms(RDKit::ROMol *mol, RDGeom::Transform3D &tform)
The class for representing 2D or 3D conformation of a molecule.
Definition: Conformer.h:41
void setBondLength(RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, double value)
void setDihedralRad(RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId, unsigned int lAtomId, double value)
RDGeom::Transform3D * computeCanonicalTransform(const RDKit::Conformer &conf, const RDGeom::Point3D *center=0, bool normalizeCovar=false, bool ignoreHs=true)
origin, the largest princiapl axis with the x-axis, the next largest with
The class for representing atoms.
Definition: Atom.h:68
#define M_PI
Definition: MMFF/Params.h:25