Package rdkit :: Package Chem :: Module rdMolAlign
[hide private]
[frames] | no frames]

Module rdMolAlign

Module containing functions to align a molecule to a second molecule

Classes [hide private]
  O3A
Open3DALIGN object
Functions [hide private]
 
AlignMol(...)
AlignMol( (Mol)prbMol, (Mol)refMol [, (int)prbCid=-1 [, (int)refCid=-1 [, (AtomPairsParameters)atomMap=[] [, (AtomPairsParameters)weights=[] [, (bool)reflect=False [, (int)maxIters=50]]]]]]) -> float : Optimally (minimum RMSD) align a molecule to another molecule
 
AlignMolConformers(...)
AlignMolConformers( (Mol)mol [, (AtomPairsParameters)atomIds=[] [, (AtomPairsParameters)confIds=[] [, (AtomPairsParameters)weights=[] [, (bool)reflect=False [, (int)maxIters=50]]]]]) -> None : Alignment conformations in a molecule to each other
 
GetAlignmentTransform(...)
GetAlignmentTransform( (Mol)prbMol, (Mol)refMol [, (int)prbCid=-1 [, (int)refCid=-1 [, (AtomPairsParameters)atomMap=[] [, (AtomPairsParameters)weights=[] [, (bool)reflect=False [, (int)maxIters=50]]]]]]) -> object : Compute the transformation required to align a molecule
 
GetCrippenO3A(...)
GetCrippenO3A( (Mol)prbMol, (Mol)refMol [, (list)prbCrippenContribs=[] [, (list)refCrippenContribs=[] [, (int)prbCid=-1 [, (int)refCid=-1 [, (bool)reflect=False [, (int)maxIters=50 [, (int)options=0 [, (list)constraintMap=[] [, (list)constraintWeights=[]]]]]]]]]]) -> O3A : Get an O3A object with atomMap and weights vectors to overlay the probe molecule onto the reference molecule based on Crippen logP atom contributions
 
GetO3A(...)
GetO3A( (Mol)prbMol, (Mol)refMol [, (AtomPairsParameters)prbPyMMFFMolProperties=None [, (AtomPairsParameters)refPyMMFFMolProperties=None [, (int)prbCid=-1 [, (int)refCid=-1 [, (bool)reflect=False [, (int)maxIters=50 [, (int)options=0 [, (list)constraintMap=[] [, (list)constraintWeights=[]]]]]]]]]]) -> O3A : Get an O3A object with atomMap and weights vectors to overlay the probe molecule onto the reference molecule based on MMFF atom types and charges
 
RandomTransform(...)
RandomTransform( (Mol)mol [, (int)cid=-1 [, (int)seed=-1]]) -> None : Perform a random transformation on a molecule
Variables [hide private]
  __package__ = None
hash(x)
Function Details [hide private]

AlignMol(...)

 

AlignMol( (Mol)prbMol, (Mol)refMol [, (int)prbCid=-1 [, (int)refCid=-1 [, (AtomPairsParameters)atomMap=[] [, (AtomPairsParameters)weights=[] [, (bool)reflect=False [, (int)maxIters=50]]]]]]) -> float :
    Optimally (minimum RMSD) align a molecule to another molecule
         
          The 3D transformation required to align the specied conformation in the probe molecule
          to a specified conformation in the reference molecule is computed so that the root mean
          squared distance between a specified set of atoms is minimized. 
          This transform is then applied to the specified conformation in the probe molecule
         
         ARGUMENTS
          - prbMol    molecule that is to be aligned
          - refMol    molecule used as the reference for the alignment
          - prbCid    ID of the conformation in the probe to be used 
                           for the alignment (defaults to first conformation)
          - refCid    ID of the conformation in the ref molecule to which 
                           the alignment is computed (defaults to first conformation)
          - atomMap   a vector of pairs of atom IDs (probe AtomId, ref AtomId)
                           used to compute the alignments. If this mapping is 
                           not specified an attempt is made to generate on by
                           substructure matching
          - weights   Optionally specify weights for each of the atom pairs
          - reflect   if true reflect the conformation of the probe molecule
          - maxIters  maximum number of iterations used in mimizing the RMSD
           
          RETURNS
          RMSD value
        
    

    C++ signature :
        double AlignMol(RDKit::ROMol {lvalue},RDKit::ROMol [,int=-1 [,int=-1 [,boost::python::api::object=[] [,boost::python::api::object=[] [,bool=False [,unsigned int=50]]]]]])

AlignMolConformers(...)

 

AlignMolConformers( (Mol)mol [, (AtomPairsParameters)atomIds=[] [, (AtomPairsParameters)confIds=[] [, (AtomPairsParameters)weights=[] [, (bool)reflect=False [, (int)maxIters=50]]]]]) -> None :
    Alignment conformations in a molecule to each other
         
          The first conformation in the molecule is used as the reference
         
         ARGUMENTS
          - mol       molecule of interest
          - atomIds   List of atom ids to use a points for alingment - defaults to all atoms
          - confIds   Ids of conformations to align - defaults to all conformers 
          - weights   Optionally specify weights for each of the atom pairs
          - reflect   if true reflect the conformation of the probe molecule
          - maxIters  maximum number of iterations used in mimizing the RMSD
           
          RETURNS
          RMSD value
        
    

    C++ signature :
        void AlignMolConformers(RDKit::ROMol {lvalue} [,boost::python::api::object=[] [,boost::python::api::object=[] [,boost::python::api::object=[] [,bool=False [,unsigned int=50]]]]])

GetAlignmentTransform(...)

 

GetAlignmentTransform( (Mol)prbMol, (Mol)refMol [, (int)prbCid=-1 [, (int)refCid=-1 [, (AtomPairsParameters)atomMap=[] [, (AtomPairsParameters)weights=[] [, (bool)reflect=False [, (int)maxIters=50]]]]]]) -> object :
    Compute the transformation required to align a molecule
         
          The 3D transformation required to align the specied conformation in the probe molecule
          to a specified conformation in the reference molecule is computed so that the root mean
          squared distance between a specified set of atoms is minimized
         
         ARGUMENTS
          - prbMol    molecule that is to be aligned
          - refMol    molecule used as the reference for the alignment
          - prbCid    ID of the conformation in the probe to be used 
                           for the alignment (defaults to first conformation)
          - refCid    ID of the conformation in the ref molecule to which 
                           the alignment is computed (defaults to first conformation)
          - atomMap   a vector of pairs of atom IDs (probe AtomId, ref AtomId)
                           used to compute the alignments. If this mapping is 
                           not specified an attempt is made to generate on by
                           substructure matching
          - weights   Optionally specify weights for each of the atom pairs
          - reflect   if true reflect the conformation of the probe molecule
          - maxIters  maximum number of iterations used in mimizing the RMSD
           
          RETURNS
          a tuple of (RMSD value, transform matrix) 
        
    

    C++ signature :
        _object* GetAlignmentTransform(RDKit::ROMol,RDKit::ROMol [,int=-1 [,int=-1 [,boost::python::api::object=[] [,boost::python::api::object=[] [,bool=False [,unsigned int=50]]]]]])

GetCrippenO3A(...)

 

GetCrippenO3A( (Mol)prbMol, (Mol)refMol [, (list)prbCrippenContribs=[] [, (list)refCrippenContribs=[] [, (int)prbCid=-1 [, (int)refCid=-1 [, (bool)reflect=False [, (int)maxIters=50 [, (int)options=0 [, (list)constraintMap=[] [, (list)constraintWeights=[]]]]]]]]]]) -> O3A :
    Get an O3A object with atomMap and weights vectors to overlay
          the probe molecule onto the reference molecule based on
          Crippen logP atom contributions
         
         ARGUMENTS
          - prbMol                   molecule that is to be aligned
          - refMol                   molecule used as the reference for the alignment
          - prbCrippenContribs       Crippen atom contributions for the probe molecule
                                     as a list of (logp, mr) tuples, as returned
                                     by _CalcCrippenContribs()
          - refCrippenContribs       Crippen atom contributions for the reference molecule
                                     as a list of (logp, mr) tuples, as returned
                                     by _CalcCrippenContribs()
          - prbCid                   ID of the conformation in the probe to be used 
                                     for the alignment (defaults to first conformation)
          - refCid                   ID of the conformation in the ref molecule to which 
                                     the alignment is computed (defaults to first conformation)
          - reflect                  if true reflect the conformation of the probe molecule
                                     (defaults to false)
          - maxIters                 maximum number of iterations used in mimizing the RMSD
                                     (defaults to 50)
          - options                  least 2 significant bits encode accuracy
                                     (0: maximum, 3: minimum; defaults to 0)
                                     bit 3 triggers local optimization of the alignment
                                     (no computation of the cost matrix; defaults: off)
          - constraintMap            a vector of pairs of atom IDs (probe AtomId, ref AtomId)
                                     which shall be used for the alignment (defaults to [])
          - constraintWeights        optionally specify weights for each of the constraints
                                     (weights default to 100.0)
           
          RETURNS
          RMSD value
        
    

    C++ signature :
        RDKit::MolAlign::PyO3A* GetCrippenO3A(RDKit::ROMol {lvalue},RDKit::ROMol {lvalue} [,boost::python::list=[] [,boost::python::list=[] [,int=-1 [,int=-1 [,bool=False [,unsigned int=50 [,unsigned int=0 [,boost::python::list=[] [,boost::python::list=[]]]]]]]]]])

GetO3A(...)

 

GetO3A( (Mol)prbMol, (Mol)refMol [, (AtomPairsParameters)prbPyMMFFMolProperties=None [, (AtomPairsParameters)refPyMMFFMolProperties=None [, (int)prbCid=-1 [, (int)refCid=-1 [, (bool)reflect=False [, (int)maxIters=50 [, (int)options=0 [, (list)constraintMap=[] [, (list)constraintWeights=[]]]]]]]]]]) -> O3A :
    Get an O3A object with atomMap and weights vectors to overlay
          the probe molecule onto the reference molecule based on
          MMFF atom types and charges
         
         ARGUMENTS
          - prbMol                   molecule that is to be aligned
          - refMol                   molecule used as the reference for the alignment
          - prbPyMMFFMolProperties   PyMMFFMolProperties object for the probe molecule as returned
                                     by SetupMMFFForceField()
          - refPyMMFFMolProperties   PyMMFFMolProperties object for the reference molecule as returned
                                     by SetupMMFFForceField()
          - prbCid                   ID of the conformation in the probe to be used 
                                     for the alignment (defaults to first conformation)
          - refCid                   ID of the conformation in the ref molecule to which 
                                     the alignment is computed (defaults to first conformation)
          - reflect                  if true reflect the conformation of the probe molecule
                                     (defaults to false)
          - maxIters                 maximum number of iterations used in mimizing the RMSD
                                     (defaults to 50)
          - options                  least 2 significant bits encode accuracy
                                     (0: maximum, 3: minimum; defaults to 0)
                                     bit 3 triggers local optimization of the alignment
                                     (no computation of the cost matrix; defaults: off)
          - constraintMap            a vector of pairs of atom IDs (probe AtomId, ref AtomId)
                                     which shall be used for the alignment (defaults to [])
          - constraintWeights        optionally specify weights for each of the constraints
                                     (weights default to 100.0)
           
          RETURNS
          RMSD value
        
    

    C++ signature :
        RDKit::MolAlign::PyO3A* GetO3A(RDKit::ROMol {lvalue},RDKit::ROMol {lvalue} [,boost::python::api::object=None [,boost::python::api::object=None [,int=-1 [,int=-1 [,bool=False [,unsigned int=50 [,unsigned int=0 [,boost::python::list=[] [,boost::python::list=[]]]]]]]]]])

RandomTransform(...)

 

RandomTransform( (Mol)mol [, (int)cid=-1 [, (int)seed=-1]]) -> None :
    Perform a random transformation on a molecule
         
         ARGUMENTS
          - mol    molecule that is to be transformed
          - cid    ID of the conformation in the mol to be transformed
                   (defaults to first conformation)
          - seed   seed used to initialize the random generator
                   (defaults to -1, that is no seeding)
           
        
    

    C++ signature :
        void RandomTransform(RDKit::ROMol {lvalue} [,int=-1 [,int=-1]])