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) |
|
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) |
|
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) |
|
OBBond * | GetBond () |
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 |
OBBitVec & | GetFixedBonds () |
double | CalcTorsion (double *coordinates) |
double | CalcBondLength (double *coordinates) |
|
std::vector< double >::iterator | BeginTorIncrement () |
std::vector< double >::iterator | EndTorIncrement () |
|
void | SetDelta (double d) |
double | GetDelta () |
OBBitVec & | GetFixedAtoms () |
void | SetFixedAtoms (OBBitVec &bv) |
OBBitVec & | GetEvalAtoms () |
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
Member Function Documentation
void SetBond |
( |
OBBond * |
bond ) |
[inline] |
void SetIdx |
( |
int |
idx ) |
[inline] |
void SetDihedralAtoms |
( |
std::vector< int > & |
ref ) |
|
void SetDihedralAtoms |
( |
int |
ref[4] ) |
|
Set the dihedral atoms.
- Parameters:
-
ref | The dihedral atom indexes. These indexes start from 1. |
void SetRotAtoms |
( |
std::vector< int > & |
atoms ) |
|
void SetTorsionValues |
( |
std::vector< double > & |
angles ) |
[inline] |
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:
-
coordinates | The coordinates to rotate. |
setang | The 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:
-
coordinates | The coordinates for the molecules as pointer to double. |
next | The index of the new rotor angle. This is an index for the GetTorsionValues() list. |
prev | If 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:
-
coordinates | The 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));
std::vector<double> angles;
angles.push_back(0.0);
angles.push_back(3.1415);
rotor.SetTorsionValues(angles);
rotor.Precompute(coords);
double *coords_1 = new double[numCoords];
for (unsigned int i = 0; i < numCoords; ++i)
coords_1[i] = coords[i];
rotor.Set(coords_1, 0)
double *coords_2 = new double[numCoords];
for (unsigned int i = 0; i < numCoords; ++i)
coords_2[i] = coords[i];
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:
-
coordinates | The coordinates to rotate. |
idx | The 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:
-
conformers | The 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:
-
coordinates | The conformer coordinates. |
conformer | The conformer index in the conformer list given to Precalc(). |
idx | The torsion value index in the GetTorsionValues() list. |
int GetIdx |
( |
void |
) |
const [inline] |
Get the index for this rotor (index in an OBRotorList).
void GetDihedralAtoms |
( |
int |
ref[4] ) |
[inline] |
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.
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:
-
- Returns:
- The torsion angle in radians.
double CalcBondLength |
( |
double * |
coordinates ) |
|
Calculate the bond length for this OBRotor using the specified coordinates.
- Parameters:
-
std::vector<double>::iterator BeginTorIncrement |
( |
) |
[inline] |
std::vector<double>::iterator EndTorIncrement |
( |
) |
[inline] |
void RemoveSymTorsionValues |
( |
int |
fold ) |
|
void SetDelta |
( |
double |
d ) |
[inline] |
double GetDelta |
( |
) |
[inline] |
void SetFixedAtoms |
( |
OBBitVec & |
bv ) |
[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: