00001 /********************************************************************** 00002 align.h - Align two molecules or vectors of vector3 00003 00004 Copyright (C) 2010 by Noel M. O'Boyle 00005 00006 This file is part of the Open Babel project. 00007 For more information, see <http://openbabel.org/> 00008 00009 This program is free software; you can redistribute it and/or modify 00010 it under the terms of the GNU General Public License as published by 00011 the Free Software Foundation version 2 of the License. 00012 00013 This program is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 ***********************************************************************/ 00018 00019 #ifndef OB_ALIGN_H 00020 #define OB_ALIGN_H 00021 00022 #include <openbabel/mol.h> 00023 #include <openbabel/math/vector3.h> 00024 #include <openbabel/math/matrix3x3.h> 00025 #include <openbabel/isomorphism.h> 00026 #include <Eigen/Core> 00027 00028 using namespace std; 00029 00030 namespace OpenBabel 00031 { 00066 class OBAPI OBAlign { 00067 public: 00069 00070 00075 OBAlign(bool includeH=false, bool symmetry=true); 00079 OBAlign(const OBMol &refmol, const OBMol &targetmol, bool includeH=false, bool symmetry=true); 00083 OBAlign(const vector<vector3> &ref, const vector<vector3> &target); 00085 00087 00088 00094 void SetRef(const vector<vector3> &ref); 00099 void SetTarget(const vector<vector3> &target); 00106 void SetRefMol(const OBMol &refmol); 00111 void SetTargetMol(const OBMol &targetmol); 00113 00115 00116 00119 bool Align(); 00121 00123 00124 00129 double GetRMSD(); 00151 matrix3x3 GetRotMatrix(); 00159 vector<vector3> GetAlignment(); 00165 bool UpdateCoords(OBMol* target); 00167 00168 private: 00169 bool _ready; 00170 bool _symmetry; 00171 bool _includeH; 00172 double _rmsd; 00173 OBBitVec _frag_atoms; 00174 Automorphisms _aut; 00175 const OBMol* _prefmol; 00176 const OBMol* _ptargetmol; 00177 Eigen::MatrixXd _rotMatrix; 00178 Eigen::Vector3d _ref_centr, _target_centr; 00179 const vector<vector3> *_pref; 00180 const vector<vector3> *_ptarget; 00181 vector<vector3> _refmol_coords; 00182 vector<vector3> _targetmol_coords; 00183 Eigen::MatrixXd _result; 00184 Eigen::MatrixXd _mref, _mtarget; 00185 void VectorsToMatrix(const vector<vector3> *pcoords, Eigen::MatrixXd &coords); 00186 Eigen::Vector3d MoveToOrigin(Eigen::MatrixXd &coords); 00187 void SimpleAlign(const Eigen::MatrixXd &mtarget); 00188 // Generate a mapping from the permutation map to the index of 00189 // correct column in _mtarget. Need to handle the fact that the 00190 // permutation group contains non-fragment atoms. 00191 // For example, map(213465) will be converted to newidx(102354). 00192 // If the atom with Idx=3 is not in the fragment, it will be 00193 // converted to newidx(10X243) instead. 00194 vector<unsigned int> _newidx; 00195 }; 00196 } 00197 00198 #endif // OB_ALIGN_H 00199
This file is part of the documentation for Open Babel, version 2.3.