Open Babel  3.0
rotor.h
Go to the documentation of this file.
1 /**********************************************************************
2 rotor.h - Rotate torsional according to rotor rules.
3 
4 Copyright (C) 1998-2000 by OpenEye Scientific Software, Inc.
5 Some portions Copyright (C) 2001-2005 by Geoffrey R. Hutchison
6 
7 This file is part of the Open Babel project.
8 For more information, see <http://openbabel.org/>
9 
10 This program is free software; you can redistribute it and/or modify
11 it under the terms of the GNU General Public License as published by
12 the Free Software Foundation version 2 of the License.
13 
14 This program is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18 ***********************************************************************/
19 
20 #ifndef OB_ROTOR_H
21 #define OB_ROTOR_H
22 
23 #include <openbabel/parsmart.h>
24 #include <openbabel/typer.h>
25 #include <openbabel/bitvec.h>
26 
27 #ifdef UNUSED
28 #elif (__GNUC__ == 4)
29 # define UNUSED(x) UNUSED_ ## x __attribute__((unused))
30 #elif defined(__LCLINT__)
31 # define UNUSED(x) /*@unused@*/ x
32 #else
33 # define UNUSED(x) x
34 #endif
35 
36 namespace OpenBabel
37 {
38  class OBRing;
39 
40 #ifndef SQUARE
41 #define SQUARE(x) ((x)*(x))
42 #endif
43 
51  class OBAPI OBRotorRule
52  {
53  int _ref[4];
54  double _delta;
55  std::string _s;
56  OBSmartsPattern* _sp;
57  std::vector<double> _vals;
58  public:
59 
60  OBRotorRule(char *buffer,int ref[4],std::vector<double> &vals,double d):
61  _delta(d), _s(buffer), _vals(vals)
62  {
63  _sp = new OBSmartsPattern;
64  _sp->Init(buffer);
65  memcpy(_ref,ref,sizeof(int)*4);
66  }
67 
69  {
70  if (_sp)
71  {
72  delete _sp;
73  _sp = NULL;
74  }
75  }
76 
78  bool IsValid() { return(_sp->IsValid()); }
82  void GetReferenceAtoms(int ref[4]) { memcpy(ref,_ref,sizeof(int)*4); }
84  void SetDelta(double d) { _delta = d; }
86  double GetDelta() { return(_delta); }
88  std::vector<double> &GetTorsionVals() { return(_vals); }
90  std::string &GetSmartsString(){ return(_s); }
92  OBSmartsPattern *GetSmartsPattern() { return(_sp); }
93  };
94 
101  class OBAPI OBRotorRules : public OBGlobalDataBase
102  {
103  bool _quiet;
104  std::vector<OBRotorRule*> _vr;
105  std::vector<double> _sp3sp3;
106  std::vector<double> _sp3sp2;
107  std::vector<double> _sp2sp2;
108  public:
109  OBRotorRules();
110  ~OBRotorRules();
111 
112  void ParseLine(const char*);
114  size_t GetSize() { return _vr.size();}
115 
117  void SetFilename(std::string &s) { _filename = s; }
118 
125  void GetRotorIncrements(OBMol& mol,OBBond* bond,int refs[4],
126  std::vector<double> &vals,double &delta);
128  void Quiet() { _quiet=true; }
129  };
130 
135  class OBAPI OBRotor
136  {
137  int _idx;
138  std::vector<int> _rotatoms;
139  double _imag, _refang;
140  OBBond *_bond;
141  std::vector<int> _ref, _torsion;
142  OBBitVec _fixedatoms,_fixedbonds, _evalatoms;
143  std::vector<double> _torsionAngles;
144  std::vector<double> _invmag;
145  std::vector<std::vector<double> > _sn,_cs,_t;
146  std::vector<OBRing *> _rings;
147  public:
151  OBRotor();
156  {
157  }
158 
161 
164  void SetBond(OBBond *bond)
165  {
166  _bond = bond;
167  SetRings();
168  }
173  void SetRings();
177  void SetIdx(int idx)
178  {
179  _idx = idx;
180  }
185  void SetDihedralAtoms(std::vector<int> &ref);
190  void SetDihedralAtoms(int ref[4]);
196  void SetRotAtoms(std::vector<int> &atoms);
200  void SetTorsionValues(std::vector<double> &angles)
201  {
202  _torsionAngles = angles;
203  }
208  {
209  _fixedbonds = bv;
210  }
212 
215 
220  inline void SetToAngle(double *coordinates, double setang)
221  {
222  double /*dx,dy,dz,*/ sn,cs,t,ang,mag;
223  // compute the angle to rotate (radians)
224  ang = setang - CalcTorsion(coordinates);
225  // if the angle to rotate is too small, we're done
226  if (fabs(ang) < 1e-5)
227  return;
228 
229  // compute the bond length
230  mag = CalcBondLength(coordinates);
231  // compute some rotation matrix elements
232  sn = sin(ang);
233  cs = cos(ang);
234  t = 1 - cs;
235 
236  // perform rotation
237  Set(coordinates, sn, cs, t, 1.0 / mag);
238  }
249  void SetRotor(double *coordinates, int next, int prev = -1);
254  void Set(double *coordinates, double sine, double cosine, double translation, double invmag);
303  void Precompute(double *coordinates);
311  void Set(double *coordinates, int idx);
318  void Precalc(std::vector<double*> &conformers);
328  void Set(double *coordinates, int conformer, int idx)
329  {
330  Set(coordinates, _sn[conformer][idx], _cs[conformer][idx], _t[conformer][idx], _invmag[conformer]);
331  }
333 
334 
337 
341  {
342  return(_bond);
343  }
348  size_t Size()
349  {
350  return _torsionAngles.size();
351  }
355  int GetIdx() const
356  {
357  return _idx;
358  }
362  void GetDihedralAtoms(int ref[4])
363  {
364  for (int i = 0; i < 4; ++i)
365  ref[i] = _ref[i];
366  }
370  std::vector<int> &GetDihedralAtoms()
371  {
372  return _ref;
373  }
378  const std::vector<int>& GetRotAtoms() const
379  {
380  return _rotatoms;
381  }
385  const std::vector<double> &GetTorsionValues() const
386  {
387  return _torsionAngles;
388  }
394  {
395  return _fixedbonds;
396  }
402  double CalcTorsion(double *coordinates);
407  double CalcBondLength(double *coordinates);
409 
410 
413  std::vector<double>::iterator BeginTorIncrement()
414  {
415  return _torsionAngles.begin();
416  }
417  std::vector<double>::iterator EndTorIncrement()
418  {
419  return _torsionAngles.end();
420  }
422 
423  void RemoveSymTorsionValues(int);
424 
427 
428  void SetDelta(double UNUSED(d)) {}
430  double GetDelta() { return 10.0; }
432  OBBitVec &GetFixedAtoms() { return _fixedatoms; }
434  void SetFixedAtoms(OBBitVec &bv) { _fixedatoms = bv; }
436  OBBitVec &GetEvalAtoms() { return _evalatoms; }
438  void SetEvalAtoms(OBBitVec &bv) { _evalatoms = bv; }
440  void* GetRotAtoms() { return &_rotatoms; }
442  std::vector<double> &GetResolution() { return _torsionAngles; }
444  void SetNumCoords(int UNUSED(nc)) {}
446 
447  };
448 
449 
451  typedef std::vector<OBRotor*>::iterator OBRotorIterator;
452 
457  class OBAPI OBRotorList
458  {
459  bool _quiet;
460  bool _removesym;
461  bool _ringRotors;
462  OBBitVec _fixedatoms, _fixedbonds;
463  OBRotorRules _rr;
464  std::vector<int> _dffv;
465  std::vector<OBRotor*> _rotor;
466  std::vector<std::pair<OBSmartsPattern*,std::pair<int,int> > > _vsym2;
469  std::vector<std::pair<OBSmartsPattern*,std::pair<int,int> > > _vsym3;
470  public:
474  OBRotorList();
478  ~OBRotorList();
479 
483  void Clear();
487  size_t Size()
488  {
489  return _rotor.size();
490  }
497  bool IsFixedBond(OBBond*);
502  {
503  return !_fixedbonds.IsEmpty();
504  }
508  void RemoveSymVals(OBMol&);
509 
515  { return _ringRotors; }
516 
518 
525  bool Setup(OBMol &mol, bool sampleRings = false);
530  {
531  _fixedbonds = fix;
532  _fixedatoms.Clear();
533  }
537  void Init(std::string &fname)
538  {
539  _rr.SetFilename(fname);
540  _rr.Init();
541  }
545  void SetQuiet() {
546  _quiet=true;
547  _rr.Quiet();
548  }
553  bool SetRotAtoms(OBMol&);
566  bool FindRotors(OBMol &mol, bool sampleRingBonds = false);
570  bool SetEvalAtoms(OBMol&);
579  bool AssignTorVals(OBMol &);
581 
583 
584 
588  OBRotor *BeginRotor(OBRotorIterator &i)
589  {
590  i = _rotor.begin();
591  return((i ==_rotor.end()) ? NULL:*i);
592  }
597  OBRotor *NextRotor(OBRotorIterator &i)
598  {
599  ++i;
600  return((i ==_rotor.end()) ? NULL:*i);
601  }
605  OBRotorIterator BeginRotors() { return(_rotor.begin()); }
609  OBRotorIterator EndRotors() { return(_rotor.end()); }
611 
614  // Not declared
616  bool IdentifyEvalAtoms(OBMol &mol) { return SetEvalAtoms(mol); }
622  {
623  _fixedatoms = fix;
624  _fixedbonds.Clear();
625  }
631  {
632  return(!_fixedatoms.IsEmpty());
633  }
636  void IgnoreSymmetryRemoval() { _removesym = false;}
640  void SetRotAtomsByFix(OBMol&);
642 
643  };
644 
646  class rotor_digit {
647  public:
648  rotor_digit(unsigned int rs)
649  {
650  resolution_size = rs;
651  state = 0;
652  }
653 
654  rotor_digit()
655  {
656  resolution_size = 0;
657  state = 0;
658  }
659 
660  void set_size(unsigned int rs)
661  {
662  resolution_size = rs;
663  state = 0;
664  }
665 
666  void set_state(int st)
667  {
668  state = st;
669  }
670 
671  int get_state()
672  {
673  return state;
674  }
675 
676  unsigned int size()
677  {
678  return resolution_size;
679  }
680 
681  bool next()
682  {
683  if (state < static_cast<int>(resolution_size - 1)) {
684  ++state;
685  return false;
686  } else
687  state = 0;
688 
689  return true;
690  }
691  private:
692  unsigned int resolution_size;
693  int state;
694 #ifndef SWIG
695  } typedef rotor_digit;
696 #else
697 };
698 #endif
699 
703  class OBAPI OBRotorKeys
704  {
750  public:
753  {
754  _vr.clear();
755  }
756 
758  void Clear(){
759  _vr.clear();
760  }
761 
763  unsigned int NumKeys()
764  {
765  unsigned int numKeys = 0;
766 
767  while (Next())
768  numKeys++;
769 
770  return numKeys;
771  }
772 
775  void AddRotor(unsigned int size)
776  {
777  rotor_digit rd(size);
778  _vr.push_back(rd);
779  }
780 
783  bool Next()
784  {
785  if(_vr.size() == 0)
786  return false;
787 
788  bool carry = _vr[0].next();
789  unsigned int i = 1;
790  while (carry) {
791  if(i == _vr.size())
792  return false;
793 
794  carry = _vr[i].next();
795  i++;
796  }
797  return true;
798  }
799 
802  std::vector<int> GetKey()
803  {
804  std::vector<int> rt;
805  rt.clear();
806  rt.push_back(0);
807  for(unsigned int i = 0; i < _vr.size(); i++){
808  rt.push_back(_vr[i].get_state());
809  }
810 
811  return rt;
812  }
813 
814  private:
815  std::vector<rotor_digit> _vr;
816  };
817 
818 
819 } // end namespace OpenBabel
820 
821 #endif // OB_ROTOR_H
822 
OBSmartsPattern * GetSmartsPattern()
Definition: rotor.h:92
void Clear()
Clear all rotors.
Definition: rotor.h:758
bool HasFixedAtoms()
Definition: rotor.h:630
void SetBond(OBBond *bond)
Definition: rotor.h:164
void GetReferenceAtoms(int ref[4])
Definition: rotor.h:82
bool Next()
Definition: rotor.h:783
void SetFixAtoms(OBBitVec &fix)
Definition: rotor.h:621
void Set(double *coordinates, int conformer, int idx)
Definition: rotor.h:328
size_t GetSize()
Definition: rotor.h:114
size_t Size()
Definition: rotor.h:348
bool IsValid() const
Definition: parsmart.h:229
void SetIdx(int idx)
Definition: rotor.h:177
~OBRotor()
Definition: rotor.h:155
void IgnoreSymmetryRemoval()
Definition: rotor.h:636
OBBitVec & GetFixedBonds()
Definition: rotor.h:393
void SetFixedAtoms(OBBitVec &bv)
Definition: rotor.h:434
bool HasFixedBonds()
Definition: rotor.h:501
Bond class.
Definition: bond.h:58
Given an OBMol, set up a list of possibly rotatable torsions,.
Definition: rotor.h:457
Molecule Class.
Definition: mol.h:118
std::vector< int > & GetDihedralAtoms()
Definition: rotor.h:370
void GetDihedralAtoms(int ref[4])
Definition: rotor.h:362
void SetEvalAtoms(OBBitVec &bv)
Definition: rotor.h:438
bool IsValid()
Definition: rotor.h:78
Database of default hybridization torsional rules and SMARTS-defined OBRotorRule objects.
Definition: rotor.h:101
OBBitVec & GetFixedAtoms()
Definition: rotor.h:432
std::vector< double > & GetResolution()
Definition: rotor.h:442
double GetDelta()
Definition: rotor.h:86
void SetDelta(double d)
Set the resolution (delta) of a torsional step in degrees.
Definition: rotor.h:84
void SetDelta(double d)
Definition: rotor.h:428
OBRotor * NextRotor(OBRotorIterator &i)
Definition: rotor.h:597
A speed-optimized vector of bits.
Definition: bitvec.h:57
OBRotorRule(char *buffer, int ref[4], std::vector< double > &vals, double d)
Definition: rotor.h:60
void SetQuiet()
Definition: rotor.h:545
size_t Size()
Definition: rotor.h:487
OBRotorKeys()
A class to generate all possible rotorKeys.
Definition: rotor.h:752
std::vector< double >::iterator BeginTorIncrement()
Definition: rotor.h:413
void Init(std::string &fname)
Definition: rotor.h:537
SMARTS (SMiles ARbitrary Target Specification) substructure searching.
Definition: parsmart.h:154
OBRotor * BeginRotor(OBRotorIterator &i)
Definition: rotor.h:588
std::vector< int > GetKey()
Definition: rotor.h:802
void AddRotor(unsigned int size)
Definition: rotor.h:775
A class to generate all possible rotorKeys.
Definition: rotor.h:703
bool IdentifyEvalAtoms(OBMol &mol)
Definition: rotor.h:616
void SetFixedBonds(OBBitVec &bv)
Definition: rotor.h:207
void * GetRotAtoms()
Definition: rotor.h:440
unsigned int NumKeys()
Number of rotor keys (= number of possible conformers)
Definition: rotor.h:763
Open Babel atom and aromaticity typer.
void SetTorsionValues(std::vector< double > &angles)
Definition: rotor.h:200
double GetDelta()
Definition: rotor.h:430
std::string & GetSmartsString()
Definition: rotor.h:90
void SetToAngle(double *coordinates, double setang)
Definition: rotor.h:220
A single rotatable OBBond as part of rotamer searching.
Definition: rotor.h:135
void Quiet()
Turn off debugging output from GetRotorIncrements()
Definition: rotor.h:128
A rule for torsional conformer searching, defined by a SMARTS pattern.
Definition: rotor.h:51
void SetNumCoords(int nc)
Definition: rotor.h:444
int GetIdx() const
Definition: rotor.h:355
void Clear()
Set all bits to zero.
Definition: bitvec.cpp:381
bool HasRingRotors()
Definition: rotor.h:514
bool IsEmpty() const
Are there no bits set to 1 in this vector?
Definition: bitvec.cpp:297
const std::vector< double > & GetTorsionValues() const
Definition: rotor.h:385
std::vector< OBRotor * >::iterator OBRotorIterator
A standard iterator over a vector of rotors.
Definition: rotor.h:451
OBBitVec & GetEvalAtoms()
Definition: rotor.h:436
std::vector< double > & GetTorsionVals()
Definition: rotor.h:88
Fast and efficient bitstring class.
const std::vector< int > & GetRotAtoms() const
Definition: rotor.h:378
bool Init(const char *pattern)
Definition: parsmart.cpp:1712
OBRotorIterator BeginRotors()
Definition: rotor.h:605
void SetFixedBonds(OBBitVec &fix)
Definition: rotor.h:529
Daylight SMARTS parser.
std::vector< double >::iterator EndTorIncrement()
Definition: rotor.h:417
void SetFilename(std::string &s)
Set the filename to be used for the database. Default = torlib.txt.
Definition: rotor.h:117
Base data table class, handles reading data files.
Definition: data.h:48
OBRotorIterator EndRotors()
Definition: rotor.h:609
Global namespace for all Open Babel code.
Definition: alias.h:22
void Init()
Read in the data file, falling back as needed.
Definition: data.cpp:666
~OBRotorRule()
Definition: rotor.h:68
OBBond * GetBond()
Definition: rotor.h:340