Open Babel  3.0
align.h
Go to the documentation of this file.
1 /**********************************************************************
2 align.h - Align two molecules or vectors of vector3
3 
4 Copyright (C) 2010 by Noel M. O'Boyle
5 
6 This file is part of the Open Babel project.
7 For more information, see <http://openbabel.org/>
8 
9 This program is free software; you can redistribute it and/or modify
10 it under the terms of the GNU General Public License as published by
11 the Free Software Foundation version 2 of the License.
12 
13 This program is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17 ***********************************************************************/
18 
19 #ifndef OB_ALIGN_H
20 #define OB_ALIGN_H
21 
22 #include <openbabel/mol.h>
23 #include <openbabel/math/vector3.h>
25 #include <openbabel/isomorphism.h>
26 #include <Eigen/Core>
27 
28 namespace OpenBabel
29 {
64  class OBAPI OBAlign {
65  public:
67 
68 
73  OBAlign(bool includeH=false, bool symmetry=true);
77  OBAlign(const OBMol &refmol, const OBMol &targetmol, bool includeH=false, bool symmetry=true);
81  OBAlign(const std::vector<vector3> &ref, const std::vector<vector3> &target);
83 
85 
86 
92  void SetRef(const std::vector<vector3> &ref);
97  void SetTarget(const std::vector<vector3> &target);
104  void SetRefMol(const OBMol &refmol);
109  void SetTargetMol(const OBMol &targetmol);
111 
113 
114 
117  bool Align();
118 
119  enum AlignMethod {
120  Kabsch = 0, // Returns matrix and RMSD
121  QCP = 1 // Returns just RMSD (fast)
122  };
123 
124  void SetMethod(enum AlignMethod method);
126 
128 
129 
134  double GetRMSD();
156  matrix3x3 GetRotMatrix();
164  std::vector<vector3> GetAlignment();
170  bool UpdateCoords(OBMol* target);
172 
173  private:
174  bool _ready;
175  bool _fail;
176  bool _symmetry;
177  bool _includeH;
178  enum AlignMethod _method;
179  double _rmsd;
180  OBBitVec _frag_atoms;
181  Automorphisms _aut;
182  const OBMol* _prefmol;
183  const OBMol* _ptargetmol;
184  Eigen::MatrixXd _rotMatrix;
185  Eigen::Vector3d _ref_centr, _target_centr;
186  const std::vector<vector3> *_pref;
187  const std::vector<vector3> *_ptarget;
188  std::vector<vector3> _refmol_coords;
189  std::vector<vector3> _targetmol_coords;
190  Eigen::MatrixXd _result;
191  Eigen::MatrixXd _mref, _mtarget;
192  void VectorsToMatrix(const std::vector<vector3> *pcoords, Eigen::MatrixXd &coords);
193  Eigen::Vector3d MoveToOrigin(Eigen::MatrixXd &coords);
194  void SimpleAlign(const Eigen::MatrixXd &mtarget);
195  void TheobaldAlign(const Eigen::MatrixXd &mtarget);
196  // Generate a mapping from the permutation map to the index of
197  // correct column in _mtarget. Need to handle the fact that the
198  // permutation group contains non-fragment atoms.
199  // For example, map(213465) will be converted to newidx(102354).
200  // If the atom with Idx=3 is not in the fragment, it will be
201  // converted to newidx(10X243) instead.
202  std::vector<unsigned int> _newidx;
203  };
204 }
205 
206 #endif // OB_ALIGN_H
207 
OBIsomorphismMapper class for finding isomorphisms.
Handle 3D coordinates.
AlignMethod
Definition: align.h:119
OBIsomorphismMapper::Mappings Automorphisms
A group of automorphic permutations.
Definition: isomorphism.h:206
Molecule Class.
Definition: mol.h:118
Handle molecules. Declarations of OBMol, OBAtom, OBBond, OBResidue. (the main header for Open Babel) ...
A speed-optimized vector of bits.
Definition: bitvec.h:57
Perform a least-squares alignment of two molecules or two vectors of vector3 objects.
Definition: align.h:64
Represents a real 3x3 matrix.
Definition: matrix3x3.h:42
Handle 3D Rotation matrix.
Global namespace for all Open Babel code.
Definition: alias.h:22