OBAlign Class Reference

Perform a least-squares alignment of two molecules or two vectors of vector3 objects. More...

#include <openbabel/math/align.h>

List of all members.

Public Member Functions

Constructors
 OBAlign (bool includeH=false, bool symmetry=true)
 OBAlign (const OBMol &refmol, const OBMol &targetmol, bool includeH=false, bool symmetry=true)
 OBAlign (const vector< vector3 > &ref, const vector< vector3 > &target)
Partial Setup
void SetRef (const vector< vector3 > &ref)
void SetTarget (const vector< vector3 > &target)
void SetRefMol (const OBMol &refmol)
void SetTargetMol (const OBMol &targetmol)
Execute the alignment
bool Align ()
Access the result of the alignment
double GetRMSD ()
matrix3x3 GetRotMatrix ()
vector< vector3GetAlignment ()
bool UpdateCoords (OBMol *target)

Detailed Description

Perform a least-squares alignment of two molecules or two vectors of vector3 objects.

This class may be used to perform a least-squares alignment of two OBMol objects or two vector<vector3> objects. The Kabsch algorithm is used for the alignment.

During the alignment, the Target is aligned to the Reference. Note that mutiple alignments to the same Reference will be much faster than multiple alignments to the same Target. When carrying out multiple alignments, a single OBAlign instance should be reused by calling SetTarget() or SetTargetMol() for each additional Target and then calling Align().

When aligning molecules, the atoms of the two molecules must be in the same order for the results to be sensible. By default, hydrogens are not included in the least-squares fitting procedure (set includeH to true if you wish to include them) and so the resulting RMSD is the heavy-atom RMSD (which is usually what you want).

By default, symmetry is taken into account when comparing molecules. For example, if a benzene is flipped by 180 degrees along one of its 2-fold symmetry axes, it will only have an RMSD of 0 (with respect to its original orientation) if symmetry is enabled. To turn off handling of symmetry set symmetry to false (this will speed up the alignment).

Note that neither the Target nor the Reference are modified by the algorithm. As a result, to update a TargetMol with the new coordinates from the alignment, you need to use UpdateCoords().

Since:
version 2.3

Constructor & Destructor Documentation

OBAlign ( bool  includeH = false,
bool  symmetry = true 
)

If this constructor is used, the Target and Reference must be set using SetRef/SetRefMol and SetTarget/SetTargetMol before running the alignment.

OBAlign ( const OBMol refmol,
const OBMol targetmol,
bool  includeH = false,
bool  symmetry = true 
)

Align two molecules.

OBAlign ( const vector< vector3 > &  ref,
const vector< vector3 > &  target 
)

Align two vectors of vector3 objects.


Member Function Documentation

void SetRef ( const vector< vector3 > &  ref )

Set the Reference (to which the Target will be aligned) in terms of a vector of vector3 objects. Note that it is faster to perform multiple alignments to the same Reference, rather than multiple alignments to the same Target.

Referenced by OBRMSDConformerScore::Score().

void SetTarget ( const vector< vector3 > &  target )

Set the Target (which will be aligned to the Reference) in terms of a vector of vector3 objects.

Referenced by OBRMSDConformerScore::Score().

void SetRefMol ( const OBMol refmol )

Set the Reference Molecule (to which the Target Molecule must be aligned). Note that is faster to perform multiple alignments to the same Reference Molecule, rather than multple alignments to the same Target Molecule.

void SetTargetMol ( const OBMol targetmol )

Set the Target Molecule (which will be aligned to the Reference Molecule).

bool Align (  )

Align the Target to the Reference using a least-squares alignment.

Referenced by OBRMSDConformerScore::Score().

double GetRMSD (  )

Return the root mean squared deviation of the target from the reference. This function should only be called after running the alignment (using Align()).

Referenced by OBRMSDConformerScore::Score().

matrix3x3 GetRotMatrix (  )

Return the rotation matrix associated with the least squares alignment. This function should only be called after running the alignment (using Align()).

The following example shows how to use the rotation matrix to rotate all of the atoms in a molecule.

 matrix3x3 rotmatrix = align.GetRotMatrix();
 for (unsigned int i = 1; i <= mol.NumAtoms(); ++i) {
    vector3 tmpvec = mol.GetAtom(i)->GetVector();
    tmpvec *= rotmatrix; //apply the rotation
    mol.GetAtom(i)->SetVector(tmpvec);
 }

Note that if you wish to use the rotation matrix to find the aligned coordinates (that is, the same coordinates returned by GetAlignment()), you should first translate the set of coordinates to the origin by subtracting the centroid, then apply the rotation, and finally add the centroid of the reference coordinates.

vector< vector3 > GetAlignment (  )

Return the actual alignment of the Target to the Reference in terms of a vector of vector3 objects. If you want an OBMol with the aligned coordinates, you should use UpdateCoords() instead. This function should only be called after running the alignment (using Align()).

bool UpdateCoords ( OBMol target )

Set the coordinates of an OBMol to those from the alignment. This function should only be called after running the alignment (using Align()).


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