RDKit
Open-source cheminformatics and machine learning.
Loading...
Searching...
No Matches
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#include <RDGeneral/export.h>
11#ifndef _RD_MOLTRANSFORMS_H_
12#define _RD_MOLTRANSFORMS_H_
13
14#include <Geometry/point.h>
15#include <Numerics/SymmMatrix.h>
16
17#ifdef RDK_HAS_EIGEN3
18#include <Eigen/Dense>
19#endif
20
21namespace RDKit {
22class ROMol;
23class Atom;
24class Conformer;
25} // namespace RDKit
26
27namespace RDGeom {
28class Transform3D;
29}
30
31namespace MolTransforms {
33 RDGeom::Transform3D &tform);
35 RDGeom::Transform3D &tform);
36
37//! Compute the centroid of a conformer
38/*!
39 This is simple the average of the heavy atom locations in the conformer,
40 not attention is paid to hydrogens or the differences in atomic radii
41
42 \param conf Conformer of interest
43 \param ignoreHs If true, ignore hydrogen atoms
44 \param weights If present, used to weigh the atomic coordinates
45*/
47 const RDKit::Conformer &conf, bool ignoreHs = true,
48 const std::vector<double> *weights = nullptr);
49
50#ifdef RDK_HAS_EIGEN3
51//! Compute principal axes and moments of inertia for a conformer
52/*!
53 These values are calculated from the inertia tensor:
54 Iij = - sum_{s=1..N}(w_s * r_{si} * r_{sj}) i != j
55 Iii = sum_{s=1..N} sum_{j!=i} (w_s * r_{sj} * r_{sj})
56 where the coordinates are relative to the center of mass.
57
58
59 \param conf Conformer of interest
60 \param axes used to return the principal axes
61 \param moments used to return the principal moments
62 \param ignoreHs If true, ignore hydrogen atoms
63 \param force If true, the calculation will be carried out even if a
64 cached value is present
65 \param weights If present, used to weigh the atomic coordinates
66
67 \returns whether or not the calculation was successful
68*/
69RDKIT_MOLTRANSFORMS_EXPORT bool computePrincipalAxesAndMoments(
70 const RDKit::Conformer &conf, Eigen::Matrix3d &axes,
71 Eigen::Vector3d &moments, bool ignoreHs = false, bool force = false,
72 const std::vector<double> *weights = nullptr);
73//! Compute principal axes and moments from the gyration matrix of a conformer
74/*!
75
76 These values are calculated from the gyration matrix/tensor:
77 Iij = sum_{s=1..N}(w_s * r_{si} * r_{sj}) i != j
78 Iii = sum_{s=1..N} sum_{t!=s}(w_s * r_{si} * r_{ti})
79 where the coordinates are relative to the center of mass.
80
81 \param conf Conformer of interest
82 \param axes used to return the principal axes
83 \param moments used to return the principal moments
84 \param ignoreHs If true, ignore hydrogen atoms
85 \param force If true, the calculation will be carried out even if a
86 cached value is present
87 \param weights If present, used to weigh the atomic coordinates
88
89 \returns whether or not the calculation was successful
90*/
92computePrincipalAxesAndMomentsFromGyrationMatrix(
93 const RDKit::Conformer &conf, Eigen::Matrix3d &axes,
94 Eigen::Vector3d &moments, bool ignoreHs = false, bool force = false,
95 const std::vector<double> *weights = nullptr);
96#endif
97
98//! Compute the transformation required to orient the conformation
99//! along the principal axes about the center; i.e. center is made to coincide
100/// with the
101//! origin, the largest principal axis with the x-axis, the next largest with
102/// the y-axis
103//! and the smallest with the z-axis
104/*!
105 If center is not specified the centroid of the conformer will be used
106 \param conf Conformer of interest
107 \param center Center to be used for canonicalization, defaults to
108 the centroid of the
109 conformation
110 \param normalizeCovar Normalize the covariance matrix with the number of
111 atoms
112 \param ignoreHs Optionally ignore hydrogens
113*/
115 const RDKit::Conformer &conf, const RDGeom::Point3D *center = nullptr,
116 bool normalizeCovar = false, bool ignoreHs = true);
117
118//! Transform the conformation using the specified transformation
120 RDKit::Conformer &conf, const RDGeom::Transform3D &trans);
121
122//! Transforms coordinates in a molecule's substance groups
124 RDKit::ROMol &mol, const RDGeom::Transform3D &trans);
125
126//! Canonicalize the orientation of a conformer so that its principal axes
127//! around the specified center point coincide with the x, y, z axes
128/*!
129 \param conf The conformer of interest
130 \param center Optional center point about which the principal axes are
131 computed
132 if not specified the centroid of the conformer will be
133 used
134 \param normalizeCovar Optionally normalize the covariance matrix by the number
135 of atoms
136 \param ignoreHs If true, ignore hydrogen atoms
137
138*/
140 RDKit::Conformer &conf, const RDGeom::Point3D *center = nullptr,
141 bool normalizeCovar = false, bool ignoreHs = true);
142
143//! Canonicalize all the conformations in a molecule
144/*!
145 \param mol the molecule of interest
146 \param normalizeCovar Optionally normalize the covariance matrix by the number
147 of atoms
148 \param ignoreHs If true, ignore hydrogens
149*/
151 bool normalizeCovar = false,
152 bool ignoreHs = true);
153
154//! Get the bond length between the specified atoms i, j
156 unsigned int iAtomId,
157 unsigned int jAtomId);
158
159//! Set the bond length between the specified atoms i, j
160//! (all atoms bonded to atom j are moved)
162 unsigned int iAtomId,
163 unsigned int jAtomId,
164 double value);
165
166//! Get the angle in radians among the specified atoms i, j, k
168 unsigned int iAtomId,
169 unsigned int jAtomId,
170 unsigned int kAtomId);
171
172//! Get the angle in degrees among the specified atoms i, j, k
173inline double getAngleDeg(const RDKit::Conformer &conf, unsigned int iAtomId,
174 unsigned int jAtomId, unsigned int kAtomId) {
175 return (180. / M_PI * getAngleRad(conf, iAtomId, jAtomId, kAtomId));
176}
177
178//! Set the angle in radians among the specified atoms i, j, k
179//! (all atoms bonded to atom k are moved)
181 unsigned int iAtomId,
182 unsigned int jAtomId,
183 unsigned int kAtomId, double value);
184
185//! Set the angle in degrees among the specified atoms i, j, k
186//! (all atoms bonded to atom k are moved)
187inline void setAngleDeg(RDKit::Conformer &conf, unsigned int iAtomId,
188 unsigned int jAtomId, unsigned int kAtomId,
189 double value) {
190 setAngleRad(conf, iAtomId, jAtomId, kAtomId, value / 180. * M_PI);
191}
192
193//! Get the dihedral angle in radians among the specified atoms i, j, k, l
195 unsigned int iAtomId,
196 unsigned int jAtomId,
197 unsigned int kAtomId,
198 unsigned int lAtomId);
199
200//! Get the dihedral angle in degrees among the specified atoms i, j, k, l
201inline double getDihedralDeg(const RDKit::Conformer &conf, unsigned int iAtomId,
202 unsigned int jAtomId, unsigned int kAtomId,
203 unsigned int lAtomId) {
204 return (180. / M_PI *
205 getDihedralRad(conf, iAtomId, jAtomId, kAtomId, lAtomId));
206}
207
208//! Set the dihedral angle in radians among the specified atoms i, j, k, l
209//! (all atoms bonded to atom l are moved)
211 RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId,
212 unsigned int kAtomId, unsigned int lAtomId, double value);
213
214//! Set the dihedral angle in degrees among the specified atoms i, j, k, l
215//! (all atoms bonded to atom l are moved)
216inline void setDihedralDeg(RDKit::Conformer &conf, unsigned int iAtomId,
217 unsigned int jAtomId, unsigned int kAtomId,
218 unsigned int lAtomId, double value) {
219 setDihedralRad(conf, iAtomId, jAtomId, kAtomId, lAtomId, value / 180. * M_PI);
220}
221} // namespace MolTransforms
222#endif
#define M_PI
Definition MMFF/Params.h:27
The class for representing atoms.
Definition Atom.h:75
The class for representing 2D or 3D conformation of a molecule.
Definition Conformer.h:46
#define RDKIT_MOLTRANSFORMS_EXPORT
Definition export.h:345
RDKIT_MOLTRANSFORMS_EXPORT void transformMolSubstanceGroups(RDKit::ROMol &mol, const RDGeom::Transform3D &trans)
Transforms coordinates in a molecule's substance groups.
RDKIT_MOLTRANSFORMS_EXPORT void canonicalizeConformer(RDKit::Conformer &conf, const RDGeom::Point3D *center=nullptr, bool normalizeCovar=false, bool ignoreHs=true)
RDKIT_MOLTRANSFORMS_EXPORT void setBondLength(RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, double value)
RDKIT_MOLTRANSFORMS_EXPORT void setDihedralRad(RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId, unsigned int lAtomId, double value)
RDKIT_MOLTRANSFORMS_EXPORT void transformMolsAtoms(RDKit::ROMol *mol, RDGeom::Transform3D &tform)
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.
RDKIT_MOLTRANSFORMS_EXPORT void transformConformer(RDKit::Conformer &conf, const RDGeom::Transform3D &trans)
Transform the conformation using the specified transformation.
RDKIT_MOLTRANSFORMS_EXPORT RDGeom::Point3D computeCentroid(const RDKit::Conformer &conf, bool ignoreHs=true, const std::vector< double > *weights=nullptr)
Compute the centroid of a conformer.
RDKIT_MOLTRANSFORMS_EXPORT void canonicalizeMol(RDKit::ROMol &mol, bool normalizeCovar=false, bool ignoreHs=true)
Canonicalize all the conformations in a molecule.
RDKIT_MOLTRANSFORMS_EXPORT 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 setAngleDeg(RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId, double value)
RDKIT_MOLTRANSFORMS_EXPORT void transformAtom(RDKit::Atom *atom, RDGeom::Transform3D &tform)
RDKIT_MOLTRANSFORMS_EXPORT void setAngleRad(RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId, double value)
RDKIT_MOLTRANSFORMS_EXPORT double getBondLength(const RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId)
Get the bond length between the specified atoms i, j.
void setDihedralDeg(RDKit::Conformer &conf, unsigned int iAtomId, unsigned int jAtomId, unsigned int kAtomId, unsigned int lAtomId, double value)
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.
RDKIT_MOLTRANSFORMS_EXPORT 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.
RDKIT_MOLTRANSFORMS_EXPORT RDGeom::Transform3D * computeCanonicalTransform(const RDKit::Conformer &conf, const RDGeom::Point3D *center=nullptr, bool normalizeCovar=false, bool ignoreHs=true)
Std stuff.