rdkit.Chem.rdMolAlign module¶
Module containing functions to align a molecule to a second molecule
- rdkit.Chem.rdMolAlign.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 one 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 minimizing 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]]]]]])
- rdkit.Chem.rdMolAlign.AlignMolConformers((Mol)mol[, (AtomPairsParameters)atomIds=[][, (AtomPairsParameters)confIds=[][, (AtomPairsParameters)weights=[][, (bool)reflect=False[, (int)maxIters=50[, (AtomPairsParameters)RMSlist=None]]]]]]) None : ¶
Align 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 alignment - 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 minimizing the RMSD
- RMSlist if provided, fills in the RMS values between the reference
conformation and the other aligned conformations
- C++ signature :
void AlignMolConformers(RDKit::ROMol {lvalue} [,boost::python::api::object=[] [,boost::python::api::object=[] [,boost::python::api::object=[] [,bool=False [,unsigned int=50 [,boost::python::api::object=None]]]]]])
- rdkit.Chem.rdMolAlign.CalcRMS((Mol)prbMol, (Mol)refMol[, (int)prbId=-1[, (int)refId=-1[, (AtomPairsParameters)map=None[, (int)maxMatches=1000000[, (bool)symmetrizeConjugatedTerminalGroups=True[, (AtomPairsParameters)weights=[]]]]]]]) float : ¶
- Returns the RMS between two molecules, taking symmetry into account.
In contrast to getBestRMS, the RMS is computed ‘in place’, i.e. probe molecules are not aligned to the reference ahead of the RMS calculation. This is useful, for example, to compute the RMSD between docking poses and the co-crystallized ligand.
Note: This function will attempt to match all permutations of matching atom orders in both molecules, for some molecules it will lead to ‘combinatorial explosion’ especially if hydrogens are present.
- ARGUMENTS
prbMol: the molecule to be aligned to the reference
refMol: the reference molecule
prbCId: (optional) probe conformation to use
refCId: (optional) reference conformation to use
- map: (optional) a list of lists of (probeAtomId, refAtomId)
tuples with the atom-atom mappings of the two molecules. If not provided, these will be generated using a substructure search.
- maxMatches: (optional) if map isn’t specified, this will be
the max number of matches found in a SubstructMatch()
- symmetrizeConjugatedTerminalGroups: (optional) if set, conjugated
terminal functional groups (like nitro or carboxylate) will be considered symmetrically
weights: (optional) weights for mapping
RETURNS The best RMSD found
- C++ signature :
double CalcRMS(RDKit::ROMol {lvalue},RDKit::ROMol {lvalue} [,int=-1 [,int=-1 [,boost::python::api::object=None [,int=1000000 [,bool=True [,boost::python::api::object=[]]]]]]])
- rdkit.Chem.rdMolAlign.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 one 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 minimizing 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]]]]]])
- rdkit.Chem.rdMolAlign.GetAllConformerBestRMS((Mol)mol[, (int)numThreads=1[, (AtomPairsParameters)map=None[, (int)maxMatches=1000000[, (bool)symmetrizeConjugatedTerminalGroups=True[, (AtomPairsParameters)weights=[]]]]]]) tuple : ¶
- Returns the symmetric distance matrix between the conformers of a molecule.
getBestRMS() is used to calculate the inter-conformer distances
- ARGUMENTS
mol: the molecule to be considered
numThreads: (optional) number of threads to use
- map: (optional) a list of lists of (probeAtomId,refAtomId)
tuples with the atom-atom mappings of the two molecules. If not provided, these will be generated using a substructure search.
- maxMatches: (optional) if map isn’t specified, this will be
the max number of matches found in a SubstructMatch()
- symmetrizeConjugatedTerminalGroups: (optional) if set, conjugated
terminal functional groups (like nitro or carboxylate) will be considered symmetrically
weights: (optional) weights for mapping
RETURNS A tuple with the best RMSDS. The ordering is [(1,0),(2,0),(2,1),(3,0),… etc]
- C++ signature :
boost::python::tuple GetAllConformerBestRMS(RDKit::ROMol {lvalue} [,int=1 [,boost::python::api::object=None [,int=1000000 [,bool=True [,boost::python::api::object=[]]]]]])
- rdkit.Chem.rdMolAlign.GetBestAlignmentTransform((Mol)prbMol, (Mol)refMol[, (int)prbCid=-1[, (int)refCid=-1[, (AtomPairsParameters)map=[][, (int)maxMatches=1000000[, (bool)symmetrizeConjugatedTerminalGroups=True[, (AtomPairsParameters)weights=[][, (bool)reflect=False[, (int)maxIters=50[, (int)numThreads=1]]]]]]]]]) object : ¶
- Compute the optimal RMS, transformation and atom map for aligning
two molecules, taking symmetry into account. Molecule coordinates are left unaltered.
This function will attempt to align all permutations of matching atom orders in both molecules, for some molecules it will lead to ‘combinatorial explosion’ especially if hydrogens are present. Use ‘GetAlignmentTransform’ to align molecules without changing the atom order.
- 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)
- map: (optional) a list of lists of (probeAtomId, refAtomId)
tuples with the atom-atom mappings of the two molecules. If not provided, these will be generated using a substructure search.
- maxMatches (optional) if atomMap is empty, this will be the max number of
matches found in a SubstructMatch().
- symmetrizeConjugatedTerminalGroups (optional) if set, conjugated
terminal functional groups (like nitro or carboxylate) will be considered symmetrically.
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 minimizing the RMSD
numThreads (optional) number of threads to use
RETURNS a tuple of (RMSD value, best transform matrix, best atom map)
- C++ signature :
_object* GetBestAlignmentTransform(RDKit::ROMol,RDKit::ROMol [,int=-1 [,int=-1 [,boost::python::api::object=[] [,int=1000000 [,bool=True [,boost::python::api::object=[] [,bool=False [,unsigned int=50 [,int=1]]]]]]]]])
- rdkit.Chem.rdMolAlign.GetBestRMS((Mol)prbMol, (Mol)refMol[, (int)prbId=-1[, (int)refId=-1[, (AtomPairsParameters)map=None[, (int)maxMatches=1000000[, (bool)symmetrizeConjugatedTerminalGroups=True[, (AtomPairsParameters)weights=[][, (int)numThreads=1]]]]]]]) float : ¶
- Returns the optimal RMS for aligning two molecules, taking
symmetry into account. As a side-effect, the probe molecule is left in the aligned state.
Note: This function will attempt to align all permutations of matching atom orders in both molecules, for some molecules it will lead to ‘combinatorial explosion’ especially if hydrogens are present. Use ‘rdkit.Chem.AllChem.AlignMol’ to align molecules without changing the atom order.
- ARGUMENTS
prbMol: the molecule to be aligned to the reference
refMol: the reference molecule
prbId: (optional) probe conformation to use
refId: (optional) reference conformation to use
- map: (optional) a list of lists of (probeAtomId,refAtomId)
tuples with the atom-atom mappings of the two molecules. If not provided, these will be generated using a substructure search.
- maxMatches: (optional) if map isn’t specified, this will be
the max number of matches found in a SubstructMatch()
- symmetrizeConjugatedTerminalGroups: (optional) if set, conjugated
terminal functional groups (like nitro or carboxylate) will be considered symmetrically
weights: (optional) weights for mapping
numThreads: (optional) number of threads to use
RETURNS The best RMSD found
- C++ signature :
double GetBestRMS(RDKit::ROMol {lvalue},RDKit::ROMol {lvalue} [,int=-1 [,int=-1 [,boost::python::api::object=None [,int=1000000 [,bool=True [,boost::python::api::object=[] [,int=1]]]]]]])
- rdkit.Chem.rdMolAlign.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 minimizing 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 The O3A object
- 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=[]]]]]]]]]])
- rdkit.Chem.rdMolAlign.GetCrippenO3AForProbeConfs((Mol)prbMol, (Mol)refMol[, (int)numThreads=1[, (list)prbCrippenContribs=[][, (list)refCrippenContribs=[][, (int)refCid=-1[, (bool)reflect=False[, (int)maxIters=50[, (int)options=0[, (list)constraintMap=[][, (list)constraintWeights=[]]]]]]]]]]) tuple : ¶
- Get a vector of O3A objects for the overlay of all
the probe molecule’s conformations 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
- numThreadsthe number of threads to use, only has an effect if
the RDKit was built with thread support (defaults to 1)
- 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()
- 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 minimizing 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 A vector of O3A objects
- C++ signature :
boost::python::tuple GetCrippenO3AForProbeConfs(RDKit::ROMol {lvalue},RDKit::ROMol {lvalue} [,int=1 [,boost::python::list=[] [,boost::python::list=[] [,int=-1 [,bool=False [,unsigned int=50 [,unsigned int=0 [,boost::python::list=[] [,boost::python::list=[]]]]]]]]]])
- rdkit.Chem.rdMolAlign.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 minimizing 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 The O3A object
- 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=[]]]]]]]]]])
- rdkit.Chem.rdMolAlign.GetO3AForProbeConfs((Mol)prbMol, (Mol)refMol[, (int)numThreads=1[, (AtomPairsParameters)prbPyMMFFMolProperties=None[, (AtomPairsParameters)refPyMMFFMolProperties=None[, (int)refCid=-1[, (bool)reflect=False[, (int)maxIters=50[, (int)options=0[, (list)constraintMap=[][, (list)constraintWeights=[]]]]]]]]]]) tuple : ¶
- Get a vector of O3A objects for the overlay of all
the probe molecule’s conformations 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
- numThreadsthe number of threads to use, only has an effect if
the RDKit was built with thread support (defaults to 1) If set to zero, the max supported by the system will be used.
- prbPyMMFFMolProperties PyMMFFMolProperties object for the probe molecule as returned
by SetupMMFFForceField()
- refPyMMFFMolProperties PyMMFFMolProperties object for the reference molecule as returned
by SetupMMFFForceField()
- 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 minimizing 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 A vector of O3A objects
- C++ signature :
boost::python::tuple GetO3AForProbeConfs(RDKit::ROMol {lvalue},RDKit::ROMol {lvalue} [,int=1 [,boost::python::api::object=None [,boost::python::api::object=None [,int=-1 [,bool=False [,unsigned int=50 [,unsigned int=0 [,boost::python::list=[] [,boost::python::list=[]]]]]]]]]])
- class rdkit.Chem.rdMolAlign.O3A¶
Bases:
instance
Open3DALIGN object
Raises an exception This class cannot be instantiated from Python
- Align((O3A)self) float : ¶
aligns probe molecule onto reference molecule
- C++ signature :
double Align(RDKit::MolAlign::PyO3A {lvalue})
- Matches((O3A)self) list : ¶
returns the AtomMap as found by Open3DALIGN
- C++ signature :
boost::python::list Matches(RDKit::MolAlign::PyO3A {lvalue})
- Score((O3A)self) float : ¶
returns the O3AScore of the alignment
- C++ signature :
double Score(RDKit::MolAlign::PyO3A {lvalue})
- Trans((O3A)self) object : ¶
returns the transformation which aligns probe molecule onto reference molecule
- C++ signature :
_object* Trans(RDKit::MolAlign::PyO3A {lvalue})
- Weights((O3A)self) list : ¶
returns the weight vector as found by Open3DALIGN
- C++ signature :
boost::python::list Weights(RDKit::MolAlign::PyO3A {lvalue})
- rdkit.Chem.rdMolAlign.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]])