00001 /********************************************************************** 00002 rotor.h - Rotate torsional according to rotor rules. 00003 00004 Copyright (C) 1998-2000 by OpenEye Scientific Software, Inc. 00005 Some portions Copyright (C) 2001-2005 by Geoffrey R. Hutchison 00006 00007 This file is part of the Open Babel project. 00008 For more information, see <http://openbabel.org/> 00009 00010 This program is free software; you can redistribute it and/or modify 00011 it under the terms of the GNU General Public License as published by 00012 the Free Software Foundation version 2 of the License. 00013 00014 This program is distributed in the hope that it will be useful, 00015 but WITHOUT ANY WARRANTY; without even the implied warranty of 00016 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00017 GNU General Public License for more details. 00018 ***********************************************************************/ 00019 00020 #ifndef OB_ROTOR_H 00021 #define OB_ROTOR_H 00022 00023 #include <openbabel/parsmart.h> 00024 #include <openbabel/typer.h> 00025 00026 namespace OpenBabel 00027 { 00028 00029 #ifndef SQUARE 00030 #define SQUARE(x) ((x)*(x)) 00031 #endif 00032 00040 class OBAPI OBRotorRule 00041 { 00042 int _ref[4]; 00043 double _delta; 00044 std::string _s; 00045 OBSmartsPattern* _sp; 00046 std::vector<double> _vals; 00047 public: 00048 00049 OBRotorRule(char *buffer,int ref[4],std::vector<double> &vals,double d): 00050 _delta(d), _s(buffer), _vals(vals) 00051 { 00052 _sp = new OBSmartsPattern; 00053 _sp->Init(buffer); 00054 memcpy(_ref,ref,sizeof(int)*4); 00055 } 00056 00057 ~OBRotorRule() 00058 { 00059 if (_sp) 00060 { 00061 delete _sp; 00062 _sp = NULL; 00063 } 00064 } 00065 00067 bool IsValid() { return(_sp->IsValid()); } 00071 void GetReferenceAtoms(int ref[4]) { memcpy(ref,_ref,sizeof(int)*4); } 00073 void SetDelta(double d) { _delta = d; } 00075 double GetDelta() { return(_delta); } 00077 std::vector<double> &GetTorsionVals() { return(_vals); } 00079 std::string &GetSmartsString(){ return(_s); } 00081 OBSmartsPattern *GetSmartsPattern() { return(_sp); } 00082 }; 00083 00090 class OBAPI OBRotorRules : public OBGlobalDataBase 00091 { 00092 bool _quiet; 00093 std::vector<OBRotorRule*> _vr; 00094 std::vector<double> _sp3sp3; 00095 std::vector<double> _sp3sp2; 00096 std::vector<double> _sp2sp2; 00097 public: 00098 OBRotorRules(); 00099 ~OBRotorRules(); 00100 00101 void ParseLine(const char*); 00103 size_t GetSize() { return _vr.size();} 00104 00106 void SetFilename(std::string &s) { _filename = s; } 00107 00114 void GetRotorIncrements(OBMol& mol,OBBond* bond,int refs[4], 00115 std::vector<double> &vals,double &delta); 00117 void Quiet() { _quiet=true; } 00118 }; 00119 00124 class OBAPI OBRotor 00125 { 00126 int _idx; 00127 std::vector<int> _rotatoms; 00128 double _imag, _refang; 00129 OBBond *_bond; 00130 std::vector<int> _ref, _torsion; 00131 OBBitVec _fixedatoms,_fixedbonds, _evalatoms; 00132 std::vector<double> _torsionAngles; 00133 std::vector<double> _invmag; 00134 std::vector<std::vector<double> > _sn,_cs,_t; 00135 public: 00139 OBRotor(); 00143 ~OBRotor() 00144 { 00145 } 00146 00149 00152 void SetBond(OBBond *bond) 00153 { 00154 _bond = bond; 00155 } 00159 void SetIdx(int idx) 00160 { 00161 _idx = idx; 00162 } 00167 void SetDihedralAtoms(std::vector<int> &ref); 00172 void SetDihedralAtoms(int ref[4]); 00178 void SetRotAtoms(std::vector<int> &atoms); 00182 void SetTorsionValues(std::vector<double> &angles) 00183 { 00184 _torsionAngles = angles; 00185 } 00189 void SetFixedBonds(OBBitVec &bv) 00190 { 00191 _fixedbonds = bv; 00192 } 00194 00197 00202 inline void SetToAngle(double *coordinates, double setang) 00203 { 00204 double /*dx,dy,dz,*/ sn,cs,t,ang,mag; 00205 // compute the angle to rotate (radians) 00206 ang = setang - CalcTorsion(coordinates); 00207 // if the angle to rotate is too small, we're done 00208 if (fabs(ang) < 1e-5) 00209 return; 00210 00211 // compute the bond length 00212 mag = CalcBondLength(coordinates); 00213 // compute some rotation matrix elements 00214 sn = sin(ang); 00215 cs = cos(ang); 00216 t = 1 - cs; 00217 00218 // perform rotation 00219 Set(coordinates, sn, cs, t, 1.0 / mag); 00220 } 00231 void SetRotor(double *coordinates, int next, int prev = -1); 00236 void Set(double *coordinates, double sine, double cosine, double translation, double invmag); 00285 void Precompute(double *coordinates); 00293 void Set(double *coordinates, int idx); 00300 void Precalc(std::vector<double*> &conformers); 00310 void Set(double *coordinates, int conformer, int idx) 00311 { 00312 Set(coordinates, _sn[conformer][idx], _cs[conformer][idx], _t[conformer][idx], _invmag[conformer]); 00313 } 00315 00316 00319 00322 OBBond *GetBond() 00323 { 00324 return(_bond); 00325 } 00330 size_t Size() 00331 { 00332 return _torsionAngles.size(); 00333 } 00337 int GetIdx() const 00338 { 00339 return _idx; 00340 } 00344 void GetDihedralAtoms(int ref[4]) 00345 { 00346 for (int i = 0; i < 4; ++i) 00347 ref[i] = _ref[i]; 00348 } 00352 std::vector<int> &GetDihedralAtoms() 00353 { 00354 return _ref; 00355 } 00360 const std::vector<int>& GetRotAtoms() const 00361 { 00362 return _rotatoms; 00363 } 00367 const std::vector<double> &GetTorsionValues() const 00368 { 00369 return _torsionAngles; 00370 } 00375 OBBitVec &GetFixedBonds() 00376 { 00377 return _fixedbonds; 00378 } 00384 double CalcTorsion(double *coordinates); 00389 double CalcBondLength(double *coordinates); 00391 00392 00395 std::vector<double>::iterator BeginTorIncrement() 00396 { 00397 return _torsionAngles.begin(); 00398 } 00399 std::vector<double>::iterator EndTorIncrement() 00400 { 00401 return _torsionAngles.end(); 00402 } 00404 00405 void RemoveSymTorsionValues(int); 00406 00409 00410 void SetDelta(double UNUSED(d)) {} 00412 double GetDelta() { return 10.0; } 00414 OBBitVec &GetFixedAtoms() { return _fixedatoms; } 00416 void SetFixedAtoms(OBBitVec &bv) { _fixedatoms = bv; } 00418 OBBitVec &GetEvalAtoms() { return _evalatoms; } 00420 void SetEvalAtoms(OBBitVec &bv) { _evalatoms = bv; } 00422 void* GetRotAtoms() { return &_rotatoms; } 00424 std::vector<double> &GetResolution() { return _torsionAngles; } 00426 void SetNumCoords(int UNUSED(nc)) {} 00428 00429 }; 00430 00431 00433 typedef std::vector<OBRotor*>::iterator OBRotorIterator; 00434 00439 class OBAPI OBRotorList 00440 { 00441 bool _quiet; 00442 bool _removesym; 00443 OBBitVec _fixedatoms, _fixedbonds; 00444 OBRotorRules _rr; 00445 std::vector<int> _dffv; 00446 std::vector<OBRotor*> _rotor; 00447 00448 std::vector<std::pair<OBSmartsPattern*,std::pair<int,int> > > _vsym2; 00450 std::vector<std::pair<OBSmartsPattern*,std::pair<int,int> > > _vsym3; 00451 public: 00455 OBRotorList(); 00459 ~OBRotorList(); 00460 00464 void Clear(); 00468 size_t Size() 00469 { 00470 return((_rotor.empty()) ? (size_t)0: _rotor.size()); 00471 } 00478 bool IsFixedBond(OBBond*); 00482 bool HasFixedBonds() 00483 { 00484 return !_fixedbonds.Empty(); 00485 } 00489 void RemoveSymVals(OBMol&); 00490 00491 00493 00499 bool Setup(OBMol &mol); 00503 void SetFixedBonds(OBBitVec &fix) 00504 { 00505 _fixedbonds = fix; 00506 _fixedatoms.Clear(); 00507 } 00511 void Init(std::string &fname) 00512 { 00513 _rr.SetFilename(fname); 00514 _rr.Init(); 00515 } 00519 void SetQuiet() { 00520 _quiet=true; 00521 _rr.Quiet(); 00522 } 00527 bool SetRotAtoms(OBMol&); 00539 bool FindRotors(OBMol &mol); 00543 bool SetEvalAtoms(OBMol&); 00552 bool AssignTorVals(OBMol &); 00554 00556 00557 00561 OBRotor *BeginRotor(OBRotorIterator &i) 00562 { 00563 i = _rotor.begin(); 00564 return((i ==_rotor.end()) ? NULL:*i); 00565 } 00570 OBRotor *NextRotor(OBRotorIterator &i) 00571 { 00572 ++i; 00573 return((i ==_rotor.end()) ? NULL:*i); 00574 } 00578 OBRotorIterator BeginRotors() { return(_rotor.begin()); } 00582 OBRotorIterator EndRotors() { return(_rotor.end()); } 00584 00587 // Not declared 00589 bool IdentifyEvalAtoms(OBMol &mol) { return SetEvalAtoms(mol); } 00594 void SetFixAtoms(OBBitVec &fix) 00595 { 00596 _fixedatoms = fix; 00597 _fixedbonds.Clear(); 00598 } 00603 bool HasFixedAtoms() 00604 { 00605 return(!_fixedatoms.Empty()); 00606 } 00609 void IgnoreSymmetryRemoval() { _removesym = false;} 00613 void SetRotAtomsByFix(OBMol&); 00615 00616 }; 00617 00619 class rotor_digit { 00620 public: 00621 rotor_digit(unsigned int rs) 00622 { 00623 resolution_size = rs; 00624 state = 0; 00625 } 00626 00627 rotor_digit() 00628 { 00629 resolution_size = 0; 00630 state = 0; 00631 } 00632 00633 void set_size(unsigned int rs) 00634 { 00635 resolution_size = rs; 00636 state = 0; 00637 } 00638 00639 void set_state(int st) 00640 { 00641 state = st; 00642 } 00643 00644 int get_state() 00645 { 00646 return state; 00647 } 00648 00649 unsigned int size() 00650 { 00651 return resolution_size; 00652 } 00653 00654 bool next() 00655 { 00656 if (state < static_cast<int>(resolution_size - 1)) { 00657 ++state; 00658 return false; 00659 } else 00660 state = 0; 00661 00662 return true; 00663 } 00664 private: 00665 unsigned int resolution_size; 00666 int state; 00667 #ifndef SWIG 00668 } typedef rotor_digit; 00669 #else 00670 }; 00671 #endif 00672 00673 00676 class OBAPI OBRotorKeys 00677 { 00723 public: 00725 OBRotorKeys() 00726 { 00727 _vr.clear(); 00728 } 00729 00731 void Clear(){ 00732 _vr.clear(); 00733 } 00734 00736 unsigned int NumKeys() 00737 { 00738 unsigned int numKeys = 0; 00739 00740 while (Next()) 00741 numKeys++; 00742 00743 return numKeys; 00744 } 00745 00748 void AddRotor(unsigned int size) 00749 { 00750 rotor_digit rd(size); 00751 _vr.push_back(rd); 00752 } 00753 00756 bool Next() 00757 { 00758 if(_vr.size() == 0) 00759 return false; 00760 00761 bool carry = _vr[0].next(); 00762 unsigned int i = 1; 00763 while (carry) { 00764 if(i == _vr.size()) 00765 return false; 00766 00767 carry = _vr[i].next(); 00768 i++; 00769 } 00770 return true; 00771 } 00772 00775 std::vector<int> GetKey() 00776 { 00777 std::vector<int> rt; 00778 rt.clear(); 00779 rt.push_back(0); 00780 for(unsigned int i = 0; i < _vr.size(); i++){ 00781 rt.push_back(_vr[i].get_state()); 00782 } 00783 00784 return rt; 00785 } 00786 00787 private: 00788 std::vector<rotor_digit> _vr; 00789 }; 00790 00791 00792 } // end namespace OpenBabel 00793 00794 #endif // OB_ROTOR_H 00795
This file is part of the documentation for Open Babel, version 2.3.