Public Member Functions

OBRotor Class Reference

A single rotatable OBBond as part of rotamer searching. More...

#include <openbabel/rotor.h>

List of all members.

Public Member Functions

 OBRotor ()
 ~OBRotor ()
void RemoveSymTorsionValues (int)
Setup
void SetBond (OBBond *bond)
void SetIdx (int idx)
void SetDihedralAtoms (std::vector< int > &ref)
void SetDihedralAtoms (int ref[4])
void SetRotAtoms (std::vector< int > &atoms)
void SetTorsionValues (std::vector< double > &angles)
void SetFixedBonds (OBBitVec &bv)
Performing rotations
void SetToAngle (double *coordinates, double setang)
void SetRotor (double *coordinates, int next, int prev=-1)
void Set (double *coordinates, double sine, double cosine, double translation, double invmag)
void Precompute (double *coordinates)
void Set (double *coordinates, int idx)
void Precalc (std::vector< double * > &conformers)
void Set (double *coordinates, int conformer, int idx)
Methods to retrieve information
OBBondGetBond ()
size_t Size ()
int GetIdx () const
void GetDihedralAtoms (int ref[4])
std::vector< int > & GetDihedralAtoms ()
const std::vector< int > & GetRotAtoms () const
const std::vector< double > & GetTorsionValues () const
OBBitVecGetFixedBonds ()
double CalcTorsion (double *coordinates)
double CalcBondLength (double *coordinates)
Iterator methods
std::vector< double >::iterator BeginTorIncrement ()
std::vector< double >::iterator EndTorIncrement ()
Deprecated
void SetDelta (double d)
double GetDelta ()
OBBitVecGetFixedAtoms ()
void SetFixedAtoms (OBBitVec &bv)
OBBitVecGetEvalAtoms ()
void SetEvalAtoms (OBBitVec &bv)
void * GetRotAtoms ()
std::vector< double > & GetResolution ()
void SetNumCoords (int nc)

Detailed Description

A single rotatable OBBond as part of rotamer searching.


Constructor & Destructor Documentation

OBRotor (  )

Constructor.

~OBRotor (  ) [inline]

Destructor.


Member Function Documentation

void SetBond ( OBBond bond ) [inline]

Set the OBBond associated with this OBRotor.

Referenced by OBRotorList::FindRotors().

void SetIdx ( int  idx ) [inline]

Set the index for this rotor. Used by OBRotorList

Referenced by OBRotorList::FindRotors().

void SetDihedralAtoms ( std::vector< int > &  ref )

Set the dihedral atoms.

Parameters:
refThe dihedral atom indexes. These indexes start from 1.

Referenced by OBRotorList::AssignTorVals(), OBRotorList::SetRotAtoms(), and OBRotorList::SetRotAtomsByFix().

void SetDihedralAtoms ( int  ref[4] )

Set the dihedral atoms.

Parameters:
refThe dihedral atom indexes. These indexes start from 1.
void SetRotAtoms ( std::vector< int > &  atoms )

Set the atom indexes that will be displaced when this rotor changes torsion angle. These indexes start from 0 and are multiplied by 3 for easy coordinate access.

Referenced by OBRotorList::AssignTorVals(), OBRotorList::SetRotAtoms(), and OBRotorList::SetRotAtomsByFix().

void SetTorsionValues ( std::vector< double > &  angles ) [inline]

Set the possible torsion values or angles.

Referenced by OBRotorList::AssignTorVals().

void SetFixedBonds ( OBBitVec bv ) [inline]

Set the bonds that will be fixed.

void SetToAngle ( double *  coordinates,
double  setang 
) [inline]

Rotate the atoms in the specified coordinates to the specified angle.

Parameters:
coordinatesThe coordinates to rotate.
setangThe new torsion angle in radians.

Referenced by OBRotorList::RemoveSymVals().

void SetRotor ( double *  coordinates,
int  next,
int  prev = -1 
)

Rotate the atoms in the specified coordinates. This function does not require any precomputation and will compute all needed information when needed.

Parameters:
coordinatesThe coordinates for the molecules as pointer to double.
nextThe index of the new rotor angle. This is an index for the GetTorsionValues() list.
prevIf specified, the torsion current torsion angle can be looked up and does not have to be calculated again.
void Set ( double *  coordinates,
double  sine,
double  cosine,
double  translation,
double  invmag 
)

Rotate the specified coordinates by using the specified rotation matrix.

void Precompute ( double *  coordinates )

Precompute the reference angle and inverse bond length of this rotor for a single conformer. This function should be used in combination with Set(double *coordinates, int idx).

Parameters:
coordinatesThe coordinates to use in the computation.
 OBMol mol;
 ...

 unsigned int numCoords = mol.NumAtoms() * 3;
 double *coords = mol.GetCoordinates();
 OBRotor rotor;
 rotor.SetBond(mol.GetBond(3));

 // set the possible torsion values
 std::vector<double> angles;
 angles.push_back(0.0);
 angles.push_back(3.1415);
 rotor.SetTorsionValues(angles);

 // precompute inverse bond length (i.e. the bond length of bond with index 3
 // using the specified coordinates) and reference angle (i.e. torsion angle
 //in coords)
 rotor.Precompute(coords);

 // copy coordinates to coords_1
 double *coords_1 = new double[numCoords];
 for (unsigned int i = 0; i < numCoords; ++i)
   coords_1[i] = coords[i];
 // rotate the atoms in coords_1 to angle with index 0 (i.e. 0.0 degrees)
 // note: on input, the coordinates should be the same as the coordinates used
 //       to precompute the inverse bond length and reference angle (in other
 //       words, the inverse magnitude and reference angle in the specfied
 //       coordinates should be the same as the one used for Precompute)
 rotor.Set(coords_1, 0)

 // copy coordinates to coords_2
 double *coords_2 = new double[numCoords];
 for (unsigned int i = 0; i < numCoords; ++i)
   coords_2[i] = coords[i];
 // rotate the atoms in coords_2 to angle with index 1 (i.e. 180.0 degrees)
 rotor.Set(coords_2, 1)

 delete coords_1;
 delete coords_2;
void Set ( double *  coordinates,
int  idx 
)

Rotate the coordinates to set the torsion angle of this rotor to the angle specified by the index idx. Make sure to call Precompute before calling this function.

Parameters:
coordinatesThe coordinates to rotate.
idxThe index of the torsion angle in the GetTorsionValues() list.
void Precalc ( std::vector< double * > &  conformers )

Precompute the inverse bond lengths, rotation matrices for all specified conformers and all possible torsion values. This method is used in combination with Set(double *coordinates, int conformer, int idx).

Parameters:
conformersThe pointers to the conformer coordinates
void Set ( double *  coordinates,
int  conformer,
int  idx 
) [inline]

Rotate the coordinates to set the torsion to the torsion value with the specified index. The coordinates should be the same as the conformer used for calling Precalc (i.e. conformers[conformer] == coordinates). Make sure to call Precalc before calling this method.

Parameters:
coordinatesThe conformer coordinates.
conformerThe conformer index in the conformer list given to Precalc().
idxThe torsion value index in the GetTorsionValues() list.
OBBond* GetBond (  ) [inline]

Get the OBBond object associated with this OBRotor.

Referenced by OBRotorList::AssignTorVals(), OBRotorList::RemoveSymVals(), and OBRotorList::SetEvalAtoms().

size_t Size (  ) [inline]

Get the number of possible torsion angles for this OBRotor. This is the length of the GetTorsionValues() list.

Referenced by OBRotorList::Setup().

int GetIdx ( void   ) const [inline]

Get the index for this rotor (index in an OBRotorList).

void GetDihedralAtoms ( int  ref[4] ) [inline]

Get the dihedral atom indexes. These indexes start from 1.

Referenced by OBRotorList::RemoveSymVals(), OBRotorList::SetRotAtoms(), OBRotorList::SetRotAtomsByFix(), OBRotorList::Setup(), and OBRotamerList::Setup().

std::vector<int>& GetDihedralAtoms (  ) [inline]

Get the dihedral atom indexes. These indexes start from 1.

const std::vector<int>& GetRotAtoms (  ) const [inline]

Get the atom indexes that will be displaced when this rotor changes torsion angle. These indexes start from 1.

const std::vector<double>& GetTorsionValues (  ) const [inline]

Get the possible torsion angles for this OBRotor.

OBBitVec& GetFixedBonds (  ) [inline]

Get an OBBitVec objects with bits set for all bonds that are fixed. Bonds are indexed from 0.

double CalcTorsion ( double *  coordinates )

Calculate the torsion for this OBRotor using the specified coordinates.

Parameters:
coordinatesThe coordinates (e.g. OBMol::GetCoordinates()).
Returns:
The torsion angle in radians.
double CalcBondLength ( double *  coordinates )

Calculate the bond length for this OBRotor using the specified coordinates.

Parameters:
coordinatesThe coordinates (e.g. OBMol::GetCoordinates()).
std::vector<double>::iterator BeginTorIncrement (  ) [inline]
std::vector<double>::iterator EndTorIncrement (  ) [inline]
void RemoveSymTorsionValues ( int  fold )
void SetDelta ( double  d ) [inline]
Deprecated:
Has no effect.

Referenced by OBRotorList::AssignTorVals().

double GetDelta (  ) [inline]
Deprecated:
Has no effect.
OBBitVec& GetFixedAtoms (  ) [inline]
void SetFixedAtoms ( OBBitVec bv ) [inline]
Deprecated:
See SetFixedBonds
OBBitVec& GetEvalAtoms (  ) [inline]
void SetEvalAtoms ( OBBitVec bv ) [inline]
void* GetRotAtoms (  ) [inline]
std::vector<double>& GetResolution (  ) [inline]
void SetNumCoords ( int  nc ) [inline]

The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines