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 for a conformer
47 /*!
48  \param conf Conformer of interest
49  \param axes used to return the principal axes
50  \param moments used to return the principal moments
51  \param ignoreHs If true, ignore hydrogen atoms
52  \param force If true, the calculation will be carried out even if a cached value is present
53  \param weights If present used to weight the atomic coordinates
54 
55  \returns whether or not the calculation was successful
56 */
57 bool computePrincipalAxesAndMoments(
58  const RDKit::Conformer &conf,
59  Eigen::Matrix3d &axes,
60  Eigen::Vector3d &moments,
61  bool ignoreHs = true,
62  bool force = false,
63  const std::vector<double> *weights = NULL);
64 #endif
65 
66 
67 
68 //! Compute the transformation require to orient the conformation
69 //! along the principal axes about the center; i.e. center is made to coincide
70 // with the
71 //! origin, the largest princiapl axis with the x-axis, the next largest with
72 // the y-axis
73 //! and the smallest with the z-axis
74 /*!
75  If center is not specified the the centroid of the conformer will be used
76  \param conf Conformer of interest
77  \param center Center to be used for canonicalization, defaults to
78  the centroid of the
79  conformation
80  \param normalizeCovar Normalize the covariance matrix with the number of
81  atoms
82  \param ignoreHs Optinally ignore hydrogens
83 */
85  const RDKit::Conformer &conf, const RDGeom::Point3D *center = 0,
86  bool normalizeCovar = false, bool ignoreHs = true);
87 
88 //! Transform the conformation using the specified transformation
90  const RDGeom::Transform3D &trans);
91 
92 //! Canonicalize the orientation of a conformer so that its principal axes
93 //! around the specified center point coincide with the x, y, z axes
94 /*!
95  \param conf The conformer of interest
96  \param center Optional center point about which the principal axes are
97  computed
98  if not specified the centroid of the conformer will be
99  used
100  \param normalizeCovar Optionally normalize the covariance matrix by the number
101  of atoms
102  \param ignoreHs If true, ignore hydrogen atoms
103 
104 */
106  const RDGeom::Point3D *center = 0,
107  bool normalizeCovar = false, bool ignoreHs = true);
108 
109 //! Canonicalize all the conformations in a molecule
110 /*!
111  \param mol the molecule of interest
112  \param normalizeCovar Optionally normalize the covariance matrix by the number
113  of atoms
114  \param ignoreHs If true, ignore hydrogens
115 */
116 void canonicalizeMol(RDKit::ROMol &mol, bool normalizeCovar = false,
117  bool ignoreHs = true);
118 
119 //! Get the bond length between the specified atoms i, j
120 double getBondLength(const RDKit::Conformer &conf, unsigned int iAtomId,
121  unsigned int jAtomId);
122 
123 //! Set the bond length between the specified atoms i, j
124 //! (all atoms bonded to atom j are moved)
125 void setBondLength(RDKit::Conformer &conf, unsigned int iAtomId,
126  unsigned int jAtomId, double value);
127 
128 //! Get the angle in radians among the specified atoms i, j, k
129 double getAngleRad(const RDKit::Conformer &conf, unsigned int iAtomId,
130  unsigned int jAtomId, unsigned int kAtomId);
131 
132 //! Get the angle in degrees among the specified atoms i, j, k
133 inline double getAngleDeg(const RDKit::Conformer &conf, unsigned int iAtomId,
134  unsigned int jAtomId, unsigned int kAtomId) {
135  return (180. / M_PI * getAngleRad(conf, iAtomId, jAtomId, kAtomId));
136 }
137 
138 //! Set the angle in radians among the specified atoms i, j, k
139 //! (all atoms bonded to atom k are moved)
140 void setAngleRad(RDKit::Conformer &conf, unsigned int iAtomId,
141  unsigned int jAtomId, unsigned int kAtomId, double value);
142 
143 //! Set the angle in degrees among the specified atoms i, j, k
144 //! (all atoms bonded to atom k are moved)
145 inline void setAngleDeg(RDKit::Conformer &conf, unsigned int iAtomId,
146  unsigned int jAtomId, unsigned int kAtomId,
147  double value) {
148  setAngleRad(conf, iAtomId, jAtomId, kAtomId, value / 180. * M_PI);
149 }
150 
151 //! Get the dihedral angle in radians among the specified atoms i, j, k, l
152 double getDihedralRad(const RDKit::Conformer &conf, unsigned int iAtomId,
153  unsigned int jAtomId, unsigned int kAtomId,
154  unsigned int lAtomId);
155 
156 //! Get the dihedral angle in degrees among the specified atoms i, j, k, l
157 inline double getDihedralDeg(const RDKit::Conformer &conf, unsigned int iAtomId,
158  unsigned int jAtomId, unsigned int kAtomId,
159  unsigned int lAtomId) {
160  return (180. / M_PI *
161  getDihedralRad(conf, iAtomId, jAtomId, kAtomId, lAtomId));
162 }
163 
164 //! Set the dihedral angle in radians among the specified atoms i, j, k, l
165 //! (all atoms bonded to atom l are moved)
166 void setDihedralRad(RDKit::Conformer &conf, unsigned int iAtomId,
167  unsigned int jAtomId, unsigned int kAtomId,
168  unsigned int lAtomId, double value);
169 
170 //! Set the dihedral angle in degrees among the specified atoms i, j, k, l
171 //! (all atoms bonded to atom l are moved)
172 inline void setDihedralDeg(RDKit::Conformer &conf, unsigned int iAtomId,
173  unsigned int jAtomId, unsigned int kAtomId,
174  unsigned int lAtomId, double value) {
175  setDihedralRad(conf, iAtomId, jAtomId, kAtomId, lAtomId, value / 180. * M_PI);
176 }
177 }
178 #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:103
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