Open Babel
3.0
|
#include <openbabel/rotor.h>
Public Member Functions | |
OBRotor () | |
~OBRotor () | |
void | RemoveSymTorsionValues (int) |
Setup | |
void | SetBond (OBBond *bond) |
void | SetRings () |
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 | |
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) |
Iterator methods | |
std::vector< double >::iterator | BeginTorIncrement () |
std::vector< double >::iterator | EndTorIncrement () |
Deprecated | |
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) |
A single rotatable OBBond as part of rotamer searching.
OBRotor | ( | ) |
Constructor.
|
inline |
Destructor.
|
inline |
Set the OBBond associated with this OBRotor.
Referenced by OBRotorList::FindRotors().
void SetRings | ( | ) |
Set the rings associated with this bond (if it's a ring bond)
|
inline |
Set the index for this rotor. Used by OBRotorList
Referenced by OBRotorList::FindRotors().
void SetDihedralAtoms | ( | std::vector< int > & | ref | ) |
Set the dihedral atoms.
ref | The 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.
ref | The 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().
|
inline |
Set the possible torsion values or angles.
Referenced by OBRotorList::AssignTorVals().
|
inline |
Set the bonds that will be fixed.
|
inline |
Rotate the atoms in the specified coordinates
to the specified angle.
coordinates | The coordinates to rotate. |
setang | The new torsion angle in radians. |
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.
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).
coordinates | The coordinates to use in the computation. |
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.
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).
conformers | The pointers to the conformer coordinates |
|
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.
coordinates | The conformer coordinates. |
conformer | The conformer index in the conformer list given to Precalc(). |
idx | The torsion value index in the GetTorsionValues() list. |
|
inline |
Get the OBBond object associated with this OBRotor.
Referenced by OBRotorList::AssignTorVals(), OBRotorList::RemoveSymVals(), OBRotorList::SetEvalAtoms(), OBRotamerList::Setup(), and OBConformerSearch::Setup().
|
inline |
Get the number of possible torsion angles for this OBRotor. This is the length of the GetTorsionValues() list.
Referenced by OBRotorList::RemoveSymVals(), and OBRotorList::Setup().
|
inline |
Get the index for this rotor (index in an OBRotorList).
Referenced by OBRotamerList::Setup().
|
inline |
Get the dihedral atom indexes. These indexes start from 1.
Referenced by OBRotorList::SetRotAtoms(), OBRotorList::SetRotAtomsByFix(), OBRotamerList::Setup(), and OBRotorList::Setup().
|
inline |
Get the dihedral atom indexes. These indexes start from 1.
|
inline |
Get the atom indexes that will be displaced when this rotor changes torsion angle. These indexes start from 1.
|
inline |
Get the possible torsion angles for this OBRotor.
|
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.
coordinates | The coordinates (e.g. OBMol::GetCoordinates()). |
double CalcBondLength | ( | double * | coordinates | ) |
Calculate the bond length for this OBRotor using the specified coordinates.
coordinates | The coordinates (e.g. OBMol::GetCoordinates()). |
|
inline |
|
inline |
void RemoveSymTorsionValues | ( | int | fold | ) |
Remove all torsions angles between 0 and 360/fold.
Referenced by OBRotorList::RemoveSymVals().
|
inline |
Referenced by OBRotorList::AssignTorVals().
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Referenced by OBRotorList::SetEvalAtoms().
|
inline |
|
inline |
Referenced by OBForceField::FastRotorSearch(), OBConformerSearch::GetConformers(), OBForceField::RandomRotorSearchInitialize(), OBRotamerList::Setup(), OBForceField::SystematicRotorSearchInitialize(), OpenBabel::UpdateConformersFromTree(), and OBForceField::WeightedRotorSearch().
|
inline |
Referenced by OBRotorList::FindRotors().