Package Chem :: Module AllChem
[hide private]
[frames] | no frames]

Module AllChem

source code

Import all RDKit chemistry modules
 



Functions [hide private]
 
TransformMol(mol, tform)
Applies the transformation to a molecule and sets it up with...
source code
 
ComputeMolShape(mol, confId=-1, boxDim=(20, 20, 20), spacing=0.5, **kwargs) source code
 
ComputeMolVolume(mol, confId=-1, gridSpacing=0.1, boxMargin=2.0) source code
 
GenerateDepictionMatching3DStructure(mol, reference, confId=-1, **kwargs) source code
 
GetBestRMS(ref, probe, refConfId=-1, probeConfId=-1, maps=None) source code
 
EnumerateLibraryFromReaction(reaction, sidechainSets)
Returns a generator for the virtual library defined by a reaction and a sequence of sidechain sets >>> import Chem >>> from Chem import AllChem >>> s1=[Chem.MolFromSmiles(x) for x in ('NC','NCC')] >>> s2=[Chem.MolFromSmiles(x) for x in ('OC=O','OC(=O)C')] >>> rxn = AllChem.ReactionFromSmarts('[O:2]=[C:1][OH].[N:3]>>[O:2]=[C:1][N:3]') >>> r = AllChem.EnumerateLibraryFromReaction(rxn,[s2,s1]) >>> [Chem.MolToSmiles(x[0]) for x in list(r)] ['CNC=O', 'CCNC=O', 'CNC(=O)C', 'CCNC(=O)C'] Note that this is all done in a lazy manner, so "infinitely" large libraries can be done without worrying about running out of memory.
source code
 
ConstrainedEmbed(mol, core, useTethers, randomseed=2342) source code
 
_test() source code
 
AddHs(...)
AddHs( (Mol)mol [, (bool)explicitOnly=False [, (bool)addCoords=False]]) -> Mol : Adds hydrogens to the graph of a molecule.
source code
 
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.
source code
 
AlignMolConformers(...)
AlignMolConformers( (Mol)mol [, (AtomPairsParameters)atomIds=[] [, (AtomPairsParameters)confIds=[] [, (AtomPairsParameters)weights=[] [, (bool)reflect=False [, (int)maxIters=50]]]]]) -> None :...
source code
 
AssignAtomChiralCodes(...)
AssignAtomChiralCodes( (Mol)mol [, (bool)cleanIt=False [, (bool)force=False]]) -> None : Does the CIP chirality assignment (R/S) for the molecule's atoms.
source code
 
AssignAtomChiralTagsFromStructure(...)
AssignAtomChiralTagsFromStructure( (Mol)mol [, (int)confId=-1 [, (bool)replaceExistingTags=True]]) -> None : Sets the chiral tags on a molecule's atoms based on a 3D conformation.
source code
 
AssignBondStereoCodes(...)
AssignBondStereoCodes( (Mol)mol [, (bool)cleanIt=False [, (bool)force=False]]) -> None : Does the CIP stereochemistry assignment (Z/E) for the molecule's bonds .
source code
 
BuildFeatureFactory(...)
BuildFeatureFactory( (str)arg1) -> MolChemicalFeatureFactory :...
source code
 
BuildFeatureFactoryFromString(...)
BuildFeatureFactoryFromString( (str)arg1) -> MolChemicalFeatureFactory :...
source code
 
CanonicalizeConformer(...)
CanonicalizeConformer( (Conformer)conf [, (Point3D)center=None [, (bool)normalizeCovar=False [, (bool)ignoreHs=True]]]) -> None :...
source code
 
CanonicalizeMol(...)
CanonicalizeMol( (Mol)mol [, (bool)normalizeCovar=False [, (bool)ignoreHs=True]]) -> None :...
source code
 
Compute2DCoords(...)
Compute2DCoords( (Mol)mol [, (bool)canonOrient=False [, (bool)clearConfs=True [, (dict)coordMap={} [, (int)nFlipsPerSample=0 [, (int)nSample=0 [, (int)sampleSeed=0 [, (bool)permuteDeg4Nodes=False]]]]]]]) -> int : Compute 2D coordinates for a molecule.
source code
 
Compute2DCoordsMimicDistmat(...)
Compute2DCoordsMimicDistmat( (Mol)mol, (AtomPairsParameters)distMat [, (bool)canonOrient=False [, (bool)clearConfs=True [, (float)weightDistMat=0.5 [, (int)nFlipsPerSample=3 [, (int)nSample=100 [, (int)sampleSeed=100 [, (bool)permuteDeg4Nodes=True]]]]]]]) -> int : Compute 2D coordinates for a molecule such that the inter-atom distances mimic those in a user-provided distance matrix.
source code
 
ComputeCanonicalTransform(...)
ComputeCanonicalTransform( (Conformer)conf [, (Point3D)center=None [, (bool)normalizeCovar=False [, (bool)ignoreHs=True]]]) -> object :...
source code
 
ComputeCentroid(...)
ComputeCentroid( (Conformer)conf [, (bool)ignoreHs=True]) -> Point3D :...
source code
 
ComputeConfBox(...)
ComputeConfBox( (Conformer)conf [, (AtomPairsParameters)trans=None [, (float)padding=2.0]]) -> tuple :...
source code
 
ComputeConfDimsAndOffset(...)
ComputeConfDimsAndOffset( (Conformer)conf [, (AtomPairsParameters)trans=None [, (float)padding=2.0]]) -> tuple :...
source code
 
ComputeGasteigerCharges(...)
ComputeGasteigerCharges( (Mol)mol [, (int)nIter=12 [, (bool)throwOnParamFailure=False]]) -> None : Compute Gasteiger partial charges for molecule The charges are computed using an iterative procedure presented in Ref : J.Gasteiger, M.
source code
 
ComputeUnionBox(...)
ComputeUnionBox( (tuple)arg1, (tuple)arg2) -> tuple :...
source code
 
DaylightFingerprint(...)
DaylightFingerprint( (Mol)mol [, (int)minPath=1 [, (int)maxPath=7 [, (int)fpSize=2048 [, (int)nBitsPerHash=4 [, (bool)useHs=True [, (float)tgtDensity=0.0 [, (int)minSize=128]]]]]]]) -> ExplicitBitVect : Returns a "Daylight"-type fingerprint for a molecule Explanation of the algorithm below.
source code
 
DeleteSubstructs(...)
DeleteSubstructs( (Mol)mol, (Mol)query [, (bool)onlyFrags=False]) -> Mol : Removes atoms matching a substructure query from a molecule ARGUMENTS: - mol: the molecule to be modified - query: the molecule to be used as a substructure query - onlyFrags: (optional) if this toggle is set, atoms will only be removed if the entire fragment in which they are found is matched by the query.
source code
 
EmbedMolecule(...)
EmbedMolecule( (Mol)mol [, (int)maxAttempts=0 [, (int)randomSeed=-1 [, (bool)clearConfs=True [, (bool)useRandomCoords=False [, (float)boxSizeMult=2.0 [, (bool)randNegEig=True [, (int)numZeroFail=1 [, (dict)coordMap={}]]]]]]]]) -> int : Use distance geometry to obtain intial coordinates for a molecule ARGUMENTS: - mol : the molecule of interest - maxAttempts : the maximum number of attempts to try embedding - randomSeed : provide a seed for the random number generator so that the same coordinates can be obtained for a molecule on multiple runs.
source code
 
EmbedMultipleConfs(...)
EmbedMultipleConfs( (Mol)mol [, (int)numConfs=10 [, (int)maxAttempts=0 [, (int)randomSeed=-1 [, (bool)clearConfs=True [, (bool)useRandomCoords=False [, (float)boxSizeMult=2.0 [, (bool)randNegEig=True [, (int)numZeroFail=1 [, (float)pruneRmsThresh=-1.0 [, (dict)coordMap={}]]]]]]]]]]) -> _vecti : Use distance geometry to obtain multiple sets of coordinates for a molecule ARGUMENTS: - mol : the molecule of interest - numConfs : the number of conformers to generate - maxAttempts : the maximum number of attempts to try embedding - randomSeed : provide a seed for the random number generator so that the same coordinates can be obtained for a molecule on multiple runs.
source code
 
EncodeShape(...)
EncodeShape( (Mol)mol, (UniformGrid3D_)grid [, (int)confId=-1 [, (AtomPairsParameters)trans=None [, (float)vdwScale=0.80000000000000004 [, (float)stepSize=0.25 [, (int)maxLayers=-1 [, (bool)ignoreHs=True]]]]]]) -> None : Encode the shape of a molecule (one of it conformer) onto a grid ARGUMENTS: - mol : the molecule of interest - grid : grid onto which the encoding is written - confId : id of the conformation of interest on mol (defaults to the first one) - trans : any transformation that needs to be used to encode onto the grid (note the molecule remain unchanged) - vdwScale : Scaling factor for the radius of the atoms to determine the base radius used in the encoding - grid points inside this sphere carry the maximum occupany - setpSize : thickness of the layers outside the base radius, the occupancy value is decreased from layer to layer from the maximum value - maxLayers : the maximum number of layers - defaults to the number allowed the number of bits use per grid point - e.g.
source code
 
FindAllPathsOfLengthN(...)
FindAllPathsOfLengthN( (Mol)mol, (int)length [, (bool)useBonds=True [, (bool)useHs=False]]) -> _listSt6vectorIiSaIiEE : Finds all paths of a particular length in a molecule ARGUMENTS: - mol: the molecule to use - length: an integer with the target length for the paths.
source code
 
FindAllSubgraphsOfLengthN(...)
FindAllSubgraphsOfLengthN( (Mol)mol, (int)length [, (bool)useHs=False]) -> _listSt6vectorIiSaIiEE : Finds all subgraphs of a particular length in a molecule ARGUMENTS: - mol: the molecule to use - length: an integer with the target number of bonds for the subgraphs.
source code
 
FindUniqueSubgraphsOfLengthN(...)
FindUniqueSubgraphsOfLengthN( (Mol)mol, (int)length [, (bool)useHs=False [, (bool)useBO=True]]) -> _listSt6vectorIiSaIiEE : Finds unique subgraphs of a particular length in a molecule ARGUMENTS: - mol: the molecule to use - length: an integer with the target number of bonds for the subgraphs.
source code
 
GetAdjacencyMatrix(...)
GetAdjacencyMatrix( (Mol)mol [, (bool)useBO=False [, (int)emptyVal=0 [, (bool)force=False [, (str)prefix='']]]]) -> object : Returns the molecule's adjacency matrix.
source code
 
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.
source code
 
GetAtomMatch(...)
GetAtomMatch( (AtomPairsParameters)featMatch [, (int)maxAts=1024]) -> object : Returns an empty list if any of the features passed in share an atom.
source code
 
GetDistanceMatrix(...)
GetDistanceMatrix( (Mol)mol [, (bool)useBO=False [, (bool)useAtomWts=False [, (bool)force=False [, (str)prefix='']]]]) -> object : Returns the molecule's distance matrix.
source code
 
GetFormalCharge(...)
GetFormalCharge( (Mol)arg1) -> int : Returns the formal charge for the molecule.
source code
 
GetMolFrags(...)
GetMolFrags( (Mol)mol [, (bool)asMols=False]) -> tuple : Finds the disconnected fragments from a molecule.
source code
 
GetMoleculeBoundsMatrix(...)
GetMoleculeBoundsMatrix( (Mol)mol [, (bool)set15bounds=True [, (bool)scaleVDW=False]]) -> object :...
source code
 
GetPeriodicTable()
Returns the application's PeriodicTable instance.
source code
 
GetSSSR(...)
GetSSSR( (Mol)arg1) -> int : Get the smallest set of simple rings for a molecule.
source code
 
GetSymmSSSR(...)
GetSymmSSSR( (Mol)arg1) -> _vectSt6vectorIiSaIiEE : Get a symmetrized SSSR for a molecule.
source code
 
Kekulize(...)
Kekulize( (Mol)mol [, (bool)clearAromaticFlags=False]) -> None : Kekulizes the molecule ARGUMENTS: - mol: the molecule to use - clearAromaticFlags: (optional) if this toggle is set, all atoms and bonds in the molecule will be marked non-aromatic following the kekulization.
source code
 
MolFromMolBlock(...)
MolFromMolBlock( (str)molBlock [, (bool)sanitize=True [, (bool)removeHs=True]]) -> Mol : Construct a molecule from a Mol block.
source code
 
MolFromMolFile(...)
MolFromMolFile( (str)molFileName [, (bool)sanitize=True [, (bool)removeHs=True]]) -> Mol : Construct a molecule from a Mol file.
source code
 
MolFromSmarts(...)
MolFromSmarts( (str)SMARTS [, (bool)mergeHs=False]) -> Mol : Construct a molecule from a SMARTS string.
source code
 
MolFromSmiles(...)
MolFromSmiles( (str)SMILES [, (bool)sanitize=True]) -> Mol : Construct a molecule from a SMILES string.
source code
 
MolFromTPLBlock(...)
MolFromTPLBlock( (str)tplBlock [, (bool)sanitize=True [, (bool)skipFirstConf=False]]) -> Mol : Construct a molecule from a TPL block.
source code
 
MolFromTPLFile(...)
MolFromTPLFile( (str)fileName [, (bool)sanitize=True [, (bool)skipFirstConf=False]]) -> Mol : Construct a molecule from a TPL file.
source code
 
MolToMolBlock(...)
MolToMolBlock( (Mol)mol [, (bool)includeStereo=False [, (int)confId=-1]]) -> str :...
source code
 
MolToSmarts(...)
MolToSmarts( (Mol)mol [, (bool)isomericSmiles=False]) -> str : Returns a SMARTS string for a molecule ARGUMENTS: - mol: the molecule - isomericSmarts: (optional) include information about stereochemistry in the SMARTS.
source code
 
MolToSmiles(...)
MolToSmiles( (Mol)mol [, (bool)isomericSmiles=False [, (bool)kekuleSmiles=False [, (int)rootedAtAtom=-1]]]) -> str : Returns the canonical SMILES string for a molecule ARGUMENTS: - mol: the molecule - isomericSmiles: (optional) include information about stereochemistry in the SMILES.
source code
 
MolToTPLBlock(...)
MolToTPLBlock( (Mol)mol [, (str)partialChargeProp='_GasteigerCharge' [, (bool)writeFirstConfTwice=False]]) -> str : Returns the Tpl block for a molecule.
source code
 
MolToTPLFile(...)
MolToTPLFile( (Mol)mol, (str)fileName [, (str)partialChargeProp='_GasteigerCharge' [, (bool)writeFirstConfTwice=False]]) -> None : Writes a molecule to a TPL file.
source code
 
RDKFingerprint(...)
RDKFingerprint( (Mol)mol [, (int)minPath=1 [, (int)maxPath=7 [, (int)fpSize=2048 [, (int)nBitsPerHash=4 [, (bool)useHs=True [, (float)tgtDensity=0.0 [, (int)minSize=128]]]]]]]) -> ExplicitBitVect : Returns an RDKit topological fingerprint for a molecule Explanation of the algorithm below.
source code
 
ReactionFromRxnBlock(...)
ReactionFromRxnBlock( (str)arg1) -> ChemicalReaction :...
source code
 
ReactionFromRxnFile(...)
ReactionFromRxnFile( (str)arg1) -> ChemicalReaction :...
source code
 
ReactionFromSmarts(...)
ReactionFromSmarts( (str)arg1) -> ChemicalReaction :...
source code
 
RemoveHs(...)
RemoveHs( (Mol)mol [, (bool)implicitOnly=False]) -> Mol : Removes any hydrogens from the graph of a molecule.
source code
 
ReplaceCore(...)
ReplaceCore( (Mol)mol, (Mol)coreQuery [, (bool)replaceDummies=True]) -> Mol : Removes the core of a molecule and labels the sidechains with dummy atoms.
source code
 
ReplaceSidechains(...)
ReplaceSidechains( (Mol)mol, (Mol)coreQuery) -> Mol : Replaces sidechains in a molecule with dummy atoms for their attachment points.
source code
 
ReplaceSubstructs(...)
ReplaceSubstructs( (Mol)mol, (Mol)query, (Mol)replacement [, (bool)replaceAll=False]) -> object : Replaces atoms matching a substructure query in a molecule ARGUMENTS: - mol: the molecule to be modified - query: the molecule to be used as a substructure query - replacement: the molecule to be used as the replacement - replaceAll: (optional) if this toggle is set, all substructures matching the query will be replaced in a single result, otherwise each result will contain a separate replacement.
source code
 
SanitizeMol(...)
SanitizeMol( (Mol)arg1) -> None : Kekulize, check valencies, set aromaticity, conjugation and hybridization - The molecule is modified in place.
source code
 
ShapeProtrudeDist(...)
ShapeProtrudeDist( (Mol)mol1, (Mol)mol2 [, (int)confId1=-1 [, (int)confId2=-1 [, (float)gridSpacing=0.5 [, (DiscreteValueType)bitsPerPoint=DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE [, (float)vdwScale=0.80000000000000004 [, (float)stepSize=0.25 [, (int)maxLayers=-1 [, (bool)ignoreHs=True [, (bool)allowReordering=True]]]]]]]]]) -> float : Compute the shape protrude distance between two molecule based on a predefined alignment ARGUMENTS: - mol1 : The first molecule of interest - mol2 : The second molecule of interest - confId1 : Conformer in the first molecule (defaults to first conformer) - confId2 : Conformer in the second molecule (defaults to first conformer) - gridSpacing : resolution of the grid used to encode the molecular shapes - bitsPerPoint : number of bit used to encode the occupancy at each grid point defaults to two bits per grid point - vdwScale : Scaling factor for the radius of the atoms to determine the base radius used in the encoding - grid points inside this sphere carry the maximum occupan - stepSize : thickness of the each layer outside the base radius, the occupancy value is decreased from layer to layer from the maximum value - maxLayers : the maximum number of layers - defaults to the number allowed the number of bits use per grid point - e.g.
source code
 
ShapeTanimotoDist(...)
ShapeTanimotoDist( (Mol)mol1, (Mol)mol2 [, (int)confId1=-1 [, (int)confId2=-1 [, (float)gridSpacing=0.5 [, (DiscreteValueType)bitsPerPoint=DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE [, (float)vdwScale=0.80000000000000004 [, (float)stepSize=0.25 [, (int)maxLayers=-1 [, (bool)ignoreHs=True]]]]]]]]) -> float : Compute the shape tanimoto distance between two molecule based on a predefined alignment ARGUMENTS: - mol1 : The first molecule of interest - mol2 : The second molecule of interest - confId1 : Conformer in the first molecule (defaults to first conformer) - confId2 : Conformer in the second molecule (defaults to first conformer) - gridSpacing : resolution of the grid used to encode the molecular shapes - bitsPerPoint : number of bit used to encode the occupancy at each grid point defaults to two bits per grid point - vdwScale : Scaling factor for the radius of the atoms to determine the base radius used in the encoding - grid points inside this sphere carry the maximum occupan - stepSize : thickness of the each layer outside the base radius, the occupancy value is decreased from layer to layer from the maximum value - maxLayers : the maximum number of layers - defaults to the number allowed the number of bits use per grid point - e.g.
source code
 
SmilesMolSupplierFromText(...)
SmilesMolSupplierFromText( (str)text [, (str)delimiter=' ' [, (int)smilesColumn=0 [, (int)nameColumn=1 [, (bool)titleLine=True [, (bool)sanitize=True]]]]]) -> SmilesMolSupplier :...
source code
 
TransformConformer(...)
TransformConformer( (Conformer)arg1, (AtomPairsParameters)arg2) -> None :...
source code
 
UFFGetMoleculeForceField(...)
UFFGetMoleculeForceField( (Mol)mol [, (float)vdwThresh=10.0 [, (int)confId=-1]]) -> ForceField :...
source code
 
UFFOptimizeMolecule(...)
UFFOptimizeMolecule( (Mol)self [, (int)maxIters=200 [, (float)vdwThresh=10.0 [, (int)confId=-1]]]) -> int :...
source code
 
WedgeMolBonds(...)
WedgeMolBonds( (Mol)arg1, (Conformer)arg2) -> None : Set the wedging on single bonds in a molecule.
source code
 
tossit()
C++ signature :...
source code
Function Details [hide private]

TransformMol(mol, tform)

source code 
Applies the transformation to a molecule and sets it up with
a single conformer

EnumerateLibraryFromReaction(reaction, sidechainSets)

source code 
Returns a generator for the virtual library defined by
 a reaction and a sequence of sidechain sets

>>> import Chem
>>> from Chem import AllChem
>>> s1=[Chem.MolFromSmiles(x) for x in ('NC','NCC')]
>>> s2=[Chem.MolFromSmiles(x) for x in ('OC=O','OC(=O)C')]
>>> rxn = AllChem.ReactionFromSmarts('[O:2]=[C:1][OH].[N:3]>>[O:2]=[C:1][N:3]')
>>> r = AllChem.EnumerateLibraryFromReaction(rxn,[s2,s1])
>>> [Chem.MolToSmiles(x[0]) for x in list(r)]
['CNC=O', 'CCNC=O', 'CNC(=O)C', 'CCNC(=O)C']

Note that this is all done in a lazy manner, so "infinitely" large libraries can
be done without worrying about running out of memory. Your patience will run out first:

Define a set of 10000 amines:
>>> amines = (Chem.MolFromSmiles('N'+'C'*x) for x in range(10000))

... a set of 10000 acids
>>> acids = (Chem.MolFromSmiles('OC(=O)'+'C'*x) for x in range(10000))

... now the virtual library (1e8 compounds in principle):
>>> r = AllChem.EnumerateLibraryFromReaction(rxn,[acids,amines])

... look at the first 4 compounds:
>>> [Chem.MolToSmiles(r.next()[0]) for x in range(4)]
['NC=O', 'CNC=O', 'CCNC=O', 'CCCNC=O']

AddHs(...)

source code 

AddHs( (Mol)mol [, (bool)explicitOnly=False [, (bool)addCoords=False]]) -> Mol :
    Adds hydrogens to the graph of a molecule.
    
      ARGUMENTS:
    
        - mol: the molecule to be modified
    
        - explicitOnly: (optional) if this toggle is set, only explicit Hs will
          be added to the molecule.  Default value is 0 (add implicit and explicit Hs).
    
        - addCoords: (optional) if this toggle is set, The Hs will have 3D coordinates
          set.  Default value is 0 (no 3D coords).
    
      RETURNS: a new molecule with added Hs
    
      NOTES:
    
        - The original molecule is *not* modified.
    
        - Much of the code assumes that Hs are not included in the molecular
          topology, so be *very* careful with the molecule that comes back from
          this function.
    
    

    C++ signature :
        RDKit::ROMol* AddHs(RDKit::ROMol [,bool=False [,bool=False]])

AlignMol(...)

source code 

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 iteration 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(...)

source code 

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 iteration 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]]]]])

AssignAtomChiralCodes(...)

source code 

AssignAtomChiralCodes( (Mol)mol [, (bool)cleanIt=False [, (bool)force=False]]) -> None :
    Does the CIP chirality assignment (R/S) 
      for the molecule's atoms.
      Chiral atoms will have a property '_CIPCode' indicating
      their chiral code.
    
      ARGUMENTS:
    
        - mol: the molecule to use
        - cleanIt: (optional) if provided, atoms with a chiral specifier that aren't
          actually chiral (e.g. atoms with duplicate substituents or only 2 substituents,
          etc.) will have their chiral code set to CHI_UNSPECIFIED
        - force: (optional) causes the calculation to be repeated, even if it has already
          been done
    
    

    C++ signature :
        void AssignAtomChiralCodes(RDKit::ROMol {lvalue} [,bool=False [,bool=False]])

AssignAtomChiralTagsFromStructure(...)

source code 

AssignAtomChiralTagsFromStructure( (Mol)mol [, (int)confId=-1 [, (bool)replaceExistingTags=True]]) -> None :
    Sets the chiral tags on a molecule's atoms based on 
      a 3D conformation.
    
      ARGUMENTS:
    
        - mol: the molecule to use
        - cleanIt: (optional) if provided, atoms with a chiral specifier that aren't
          actually chiral (e.g. atoms with duplicate substituents or only 2 substituents,
          etc.) will have their chiral code set to CHI_UNSPECIFIED
        - force: (optional) causes the calculation to be repeated, even if it has already
          been done
    
    

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

AssignBondStereoCodes(...)

source code 

AssignBondStereoCodes( (Mol)mol [, (bool)cleanIt=False [, (bool)force=False]]) -> None :
    Does the CIP stereochemistry assignment (Z/E)
       for the molecule's bonds .
      Qualifying bonds will have a property '_CIPCode' indicating
      their stereochemistry.
    
      ARGUMENTS:
    
        - mol: the molecule to use
        - cleanIt: (optional) ignored
        - force: (optional) causes the calculation to be repeated, even if it has already
          been done
    
    

    C++ signature :
        void AssignBondStereoCodes(RDKit::ROMol {lvalue} [,bool=False [,bool=False]])

BuildFeatureFactory(...)

source code 

BuildFeatureFactory( (str)arg1) -> MolChemicalFeatureFactory :
    Construct a feature factory given a feature definition in a file

    C++ signature :
        RDKit::MolChemicalFeatureFactory* BuildFeatureFactory(std::string)

BuildFeatureFactoryFromString(...)

source code 

BuildFeatureFactoryFromString( (str)arg1) -> MolChemicalFeatureFactory :
    Construct a feature factory given a feature definition block

    C++ signature :
        RDKit::MolChemicalFeatureFactory* BuildFeatureFactoryFromString(std::string)

CanonicalizeConformer(...)

source code 

CanonicalizeConformer( (Conformer)conf [, (Point3D)center=None [, (bool)normalizeCovar=False [, (bool)ignoreHs=True]]]) -> None :
    Canonicalize the orientation of a conformer so that its principal axes
                   around the specified center point coincide with the x, y, z axes
      
      ARGUMENTS:
        - conf : conformer of interest 
        - center : optionally center point about which the principal axes are computed 
                              if not specified the centroid of the conformer will be used
        - normalizeCovar : Optionally normalize the covariance matrix by the number of atoms
    

    C++ signature :
        void CanonicalizeConformer(RDKit::Conformer {lvalue} [,RDGeom::Point3D const*=None [,bool=False [,bool=True]]])

CanonicalizeMol(...)

source code 

CanonicalizeMol( (Mol)mol [, (bool)normalizeCovar=False [, (bool)ignoreHs=True]]) -> None :
    Loop over the conformers in a molecule and canonicalize their orientation

    C++ signature :
        void CanonicalizeMol(RDKit::ROMol {lvalue} [,bool=False [,bool=True]])

Compute2DCoords(...)

source code 

Compute2DCoords( (Mol)mol [, (bool)canonOrient=False [, (bool)clearConfs=True [, (dict)coordMap={} [, (int)nFlipsPerSample=0 [, (int)nSample=0 [, (int)sampleSeed=0 [, (bool)permuteDeg4Nodes=False]]]]]]]) -> int :
    Compute 2D coordinates for a molecule. 
      The resulting coordinates are stored on each atom of the molecule 
    
      ARGUMENTS: 
    
         mol - the molecule of interest
         canonOrient - orient the molecule in a canonical way
         clearConfs - if true, all existing conformations on the molecule
                 will be cleared
         coordMap - a dictionary mapping atom Ids -> Point2D objects 
                    with starting coordinates for atoms that should
                    have their positions locked.
         nFlipsPerSample - number of rotatable bonds that are
                    flipped at random at a time.
         nSample - Number of random samplings of rotatable bonds.
         sampleSeed - seed for the random sampling process.
         permuteDeg4Nodes - allow permutation of bonds at a degree 4
                     node during the sampling process 
    
      RETURNS: 
    
         ID of the conformation added to the molecule
    

    C++ signature :
        unsigned int Compute2DCoords(RDKit::ROMol {lvalue} [,bool=False [,bool=True [,boost::python::dict {lvalue}={} [,unsigned int=0 [,unsigned int=0 [,int=0 [,bool=False]]]]]]])

Compute2DCoordsMimicDistmat(...)

source code 

Compute2DCoordsMimicDistmat( (Mol)mol, (AtomPairsParameters)distMat [, (bool)canonOrient=False [, (bool)clearConfs=True [, (float)weightDistMat=0.5 [, (int)nFlipsPerSample=3 [, (int)nSample=100 [, (int)sampleSeed=100 [, (bool)permuteDeg4Nodes=True]]]]]]]) -> int :
    Compute 2D coordinates for a molecule such 
      that the inter-atom distances mimic those in a user-provided
      distance matrix. 
      The resulting coordinates are stored on each atom of the molecule 
    
      ARGUMENTS: 
    
         mol - the molecule of interest
         distMat - distance matrix that we want the 2D structure to mimic
         canonOrient - orient the molecule in a canonical way
         clearConfs - if true, all existing conformations on the molecule
                 will be cleared
         weightDistMat - weight assigned in the cost function to mimicing
                         the distance matrix.
                         This must be between (0.0,1.0). (1.0-weightDistMat)
                         is then the weight assigned to improving 
                         the density of the 2D structure i.e. try to
                         make it spread out
         nFlipsPerSample - number of rotatable bonds that are
                    flipped at random at a time.
         nSample - Number of random samplings of rotatable bonds.
         sampleSeed - seed for the random sampling process.
         permuteDeg4Nodes - allow permutation of bonds at a degree 4
                     node during the sampling process 
    
      RETURNS: 
    
         ID of the conformation added to the molecule
    

    C++ signature :
        unsigned int Compute2DCoordsMimicDistmat(RDKit::ROMol {lvalue},boost::python::api::object [,bool=False [,bool=True [,double=0.5 [,unsigned int=3 [,unsigned int=100 [,int=100 [,bool=True]]]]]]])

ComputeCanonicalTransform(...)

source code 

ComputeCanonicalTransform( (Conformer)conf [, (Point3D)center=None [, (bool)normalizeCovar=False [, (bool)ignoreHs=True]]]) -> object :
    Compute the transformation required aligna conformer so that
                   the the principal axes align up with the x,y, z axes
                   The conformer itself is left unchanged
      ARGUMENTS:
        - conf : the conformer of interest
        - center : optional center point to compute the principal axes around (defaults to the centroid)
        - normalizeCovar : optionally normalize the covariance matrix by the number of atoms
    

    C++ signature :
        _object* ComputeCanonicalTransform(RDKit::Conformer [,RDGeom::Point3D const*=None [,bool=False [,bool=True]]])

ComputeCentroid(...)

source code 

ComputeCentroid( (Conformer)conf [, (bool)ignoreHs=True]) -> Point3D :
    Compute the centroid of the conformation - hydrogens are ignored and no attention
                               if paid to the difference in sizes of the heavy atoms
    

    C++ signature :
        RDGeom::Point3D ComputeCentroid(RDKit::Conformer [,bool=True])

ComputeConfBox(...)

source code 

ComputeConfBox( (Conformer)conf [, (AtomPairsParameters)trans=None [, (float)padding=2.0]]) -> tuple :
    Compute the lower and upper corners of a cuboid that will fit the conformer

    C++ signature :
        boost::python::tuple ComputeConfBox(RDKit::Conformer [,boost::python::api::object=None [,double=2.0]])

ComputeConfDimsAndOffset(...)

source code 

ComputeConfDimsAndOffset( (Conformer)conf [, (AtomPairsParameters)trans=None [, (float)padding=2.0]]) -> tuple :
    Compute the size of the box that can fit the conformations, and offset 
       of the box from the origin
    

    C++ signature :
        boost::python::tuple ComputeConfDimsAndOffset(RDKit::Conformer [,boost::python::api::object=None [,double=2.0]])

ComputeGasteigerCharges(...)

source code 

ComputeGasteigerCharges( (Mol)mol [, (int)nIter=12 [, (bool)throwOnParamFailure=False]]) -> None :
    Compute Gasteiger partial charges for molecule
    
     The charges are computed using an iterative procedure presented in 
     
     Ref : J.Gasteiger, M. Marseli, Iterative Equalization of Oribital Electronegatiity 
     A Rapid Access to Atomic Charges, Tetrahedron Vol 36 p3219 1980
     
     The computed charges are stored on each atom are stored a computed property ( under the name 
     _GasteigerCharge). In addition, each atom also stored the total charge for the implicit hydrogens 
     on the atom (under the property name _GasteigerHCharge)
     
     ARGUMENTS:
    
        - mol : the molecule of interrest
        - nIter : number of iteration (defaults to 12)
        - throwOnParamFailure : toggles whether or not an exception should be raised if parameters
          for an atom cannot be found.  If this is false (the default), all parameters for unknown
          atoms will be set to zero.  This has the effect of removing that atom from the iteration.
    
    

    C++ signature :
        void ComputeGasteigerCharges(RDKit::ROMol const* [,int=12 [,bool=False]])

ComputeUnionBox(...)

source code 

ComputeUnionBox( (tuple)arg1, (tuple)arg2) -> tuple :
    Compute the union of two boxes, so that all the points in both boxes are 
        contain in the new box

    C++ signature :
        boost::python::tuple ComputeUnionBox(boost::python::tuple,boost::python::tuple)

DaylightFingerprint(...)

source code 

DaylightFingerprint( (Mol)mol [, (int)minPath=1 [, (int)maxPath=7 [, (int)fpSize=2048 [, (int)nBitsPerHash=4 [, (bool)useHs=True [, (float)tgtDensity=0.0 [, (int)minSize=128]]]]]]]) -> ExplicitBitVect :
    Returns a "Daylight"-type fingerprint for a molecule
    
      Explanation of the algorithm below.
    
      ARGUMENTS:
    
        - mol: the molecule to use
    
        - minPath: (optional) minimum number of bonds to include in the subgraphs
          Defaults to 1.
    
        - maxPath: (optional) maximum number of bonds to include in the subgraphs
          Defaults to 7.
    
        - fpSize: (optional) number of bits in the fingerprint
          Defaults to 2048.
    
        - nBitsPerPath: (optional) number of bits to set per path
          Defaults to 4.
    
        - useHs: (optional) include information about number of Hs on each
          atom when calculating path hashes.
          Defaults to 1.
    
        - tgtDensity: (optional) fold the fingerprint until this minimum density has
          been reached
          Defaults to 0.
    
        - minSize: (optional) the minimum size the fingerprint will be folded to when
          trying to reach tgtDensity
          Defaults to 128.
    
      RETURNS: a DataStructs.ExplicitBitVect with _fpSize_ bits
    
      ALGORITHM:
    
       This algorithm functions by find all paths between minPath and maxPath in
        length.  For each path:
    
         1) The Balaban J value is calculated.
    
         2) The 32 bit Balaban J value is used to seed a random-number generator
    
         3) _nBitsPerPath_ random numbers are generated and used to set the corresponding
            bits in the fingerprint
    
    
    

    C++ signature :
        ExplicitBitVect* DaylightFingerprint(RDKit::ROMol [,unsigned int=1 [,unsigned int=7 [,unsigned int=2048 [,unsigned int=4 [,bool=True [,double=0.0 [,unsigned int=128]]]]]]])

DeleteSubstructs(...)

source code 

DeleteSubstructs( (Mol)mol, (Mol)query [, (bool)onlyFrags=False]) -> Mol :
    Removes atoms matching a substructure query from a molecule
    
      ARGUMENTS:
    
        - mol: the molecule to be modified
    
        - query: the molecule to be used as a substructure query
    
        - onlyFrags: (optional) if this toggle is set, atoms will only be removed if
          the entire fragment in which they are found is matched by the query.
          See below for examples.
          Default value is 0 (remove the atoms whether or not the entire fragment matches)
    
      RETURNS: a new molecule with the substructure removed
    
      NOTES:
    
        - The original molecule is *not* modified.
    
      EXAMPLES:
    
       The following examples substitute SMILES/SMARTS strings for molecules, you'd have
       to actually use molecules:
    
        - DeleteSubstructs('CCOC','OC') -> 'CC'
    
        - DeleteSubstructs('CCOC','OC',1) -> 'CCOC'
    
        - DeleteSubstructs('CCOCCl.Cl','Cl',1) -> 'CCOCCl'
    
        - DeleteSubstructs('CCOCCl.Cl','Cl') -> 'CCOC'
    
    

    C++ signature :
        RDKit::ROMol* DeleteSubstructs(RDKit::ROMol,RDKit::ROMol [,bool=False])

EmbedMolecule(...)

source code 

EmbedMolecule( (Mol)mol [, (int)maxAttempts=0 [, (int)randomSeed=-1 [, (bool)clearConfs=True [, (bool)useRandomCoords=False [, (float)boxSizeMult=2.0 [, (bool)randNegEig=True [, (int)numZeroFail=1 [, (dict)coordMap={}]]]]]]]]) -> int :
    Use distance geometry to obtain intial 
     coordinates for a molecule
    
     
     ARGUMENTS:
    
        - mol : the molecule of interest
        - maxAttempts : the maximum number of attempts to try embedding 
        - randomSeed : provide a seed for the random number generator 
                       so that the same coordinates can be obtained 
                       for a molecule on multiple runs. The default 
                       (-1) uses a random seed 
        - clearConfs : clear all existing conformations on the molecule
        - useRandomCoords : Start the embedding from random coordinates instead of
                            using eigenvalues of the distance matrix.
        - boxSizeMult    Determines the size of the box that is used for
                         random coordinates. If this is a positive number, the 
                         side length will equal the largest element of the distance
                         matrix times boxSizeMult. If this is a negative number,
                         the side length will equal -boxSizeMult (i.e. independent
                         of the elements of the distance matrix).
        - randNegEig : If the embedding yields a negative eigenvalue, 
                       pick coordinates that correspond 
                       to this component at random 
        - numZeroFail : fail embedding is we have this more zero eigenvalues 
    
     RETURNS:
    
        ID of the new conformation added to the molecule 
    
    

    C++ signature :
        int EmbedMolecule(RDKit::ROMol {lvalue} [,unsigned int=0 [,int=-1 [,bool=True [,bool=False [,double=2.0 [,bool=True [,unsigned int=1 [,boost::python::dict {lvalue}={}]]]]]]]])

EmbedMultipleConfs(...)

source code 

EmbedMultipleConfs( (Mol)mol [, (int)numConfs=10 [, (int)maxAttempts=0 [, (int)randomSeed=-1 [, (bool)clearConfs=True [, (bool)useRandomCoords=False [, (float)boxSizeMult=2.0 [, (bool)randNegEig=True [, (int)numZeroFail=1 [, (float)pruneRmsThresh=-1.0 [, (dict)coordMap={}]]]]]]]]]]) -> _vecti :
    Use distance geometry to obtain multiple sets of 
     coordinates for a molecule
     
     ARGUMENTS:
    
      - mol : the molecule of interest
      - numConfs : the number of conformers to generate 
      - maxAttempts : the maximum number of attempts to try embedding 
      - randomSeed : provide a seed for the random number generator 
                     so that the same coordinates can be obtained 
                     for a molecule on multiple runs. The default 
                     (-1) uses a random seed 
      - clearConfs : clear all existing conformations on the molecule
      - useRandomCoords : Start the embedding from random coordinates instead of
                          using eigenvalues of the distance matrix.
      - boxSizeMult    Determines the size of the box that is used for
                       random coordinates. If this is a positive number, the 
                       side length will equal the largest element of the distance
                       matrix times boxSizeMult. If this is a negative number,
                       the side length will equal -boxSizeMult (i.e. independent
                       of the elements of the distance matrix).
      - randNegEig : If the embedding yields a negative eigenvalue, 
                     pick coordinates that correspond 
                     to this component at random 
      - numZeroFail : fail embedding is we have this more zero eigenvalues 
      - pruneRmsThresh : Retain only the conformations out of 'numConfs' 
                        after embedding that are at least 
                        this far apart from each other. 
              RMSD is computed on the heavy atoms. 
              Pruning is greedy; i.e. the first embedded conformation
              is retained and from then on only those that are at
              least pruneRmsThresh away from all retained conformations
              are kept. The pruning is done after embedding and 
              bounds violation minimization. No pruning by default.
     RETURNS:
    
        List of new conformation IDs 
    
    

    C++ signature :
        std::vector<int, std::allocator<int> > EmbedMultipleConfs(RDKit::ROMol {lvalue} [,unsigned int=10 [,unsigned int=0 [,int=-1 [,bool=True [,bool=False [,double=2.0 [,bool=True [,unsigned int=1 [,double=-1.0 [,boost::python::dict {lvalue}={}]]]]]]]]]])

EncodeShape(...)

source code 

EncodeShape( (Mol)mol, (UniformGrid3D_)grid [, (int)confId=-1 [, (AtomPairsParameters)trans=None [, (float)vdwScale=0.80000000000000004 [, (float)stepSize=0.25 [, (int)maxLayers=-1 [, (bool)ignoreHs=True]]]]]]) -> None :
    Encode the shape of a molecule (one of it conformer) onto a grid
    
     
     ARGUMENTS:
    
        - mol : the molecule of interest
        - grid : grid onto which the encoding is written 
        - confId : id of the conformation of interest on mol (defaults to the first one) 
        - trans : any transformation that needs to be used to encode onto the grid (note the molecule remain unchanged) 
        - vdwScale : Scaling factor for the radius of the atoms to determine the base radius 
                     used in the encoding - grid points inside this sphere carry the maximum occupany 
        - setpSize : thickness of the layers outside the base radius, the occupancy value is decreased 
                     from layer to layer from the maximum value 
        - maxLayers : the maximum number of layers - defaults to the number allowed the number of bits 
                      use per grid point - e.g. two bits per grid point will allow 3 layers
        - ignoreHs : when set, the contribution of Hs to the shape will be ignored
    

    C++ signature :
        void EncodeShape(RDKit::ROMol,RDGeom::UniformGrid3D {lvalue} [,int=-1 [,boost::python::api::object=None [,double=0.80000000000000004 [,double=0.25 [,int=-1 [,bool=True]]]]]])

FindAllPathsOfLengthN(...)

source code 

FindAllPathsOfLengthN( (Mol)mol, (int)length [, (bool)useBonds=True [, (bool)useHs=False]]) -> _listSt6vectorIiSaIiEE :
    Finds all paths of a particular length in a molecule
    
      ARGUMENTS:
    
        - mol: the molecule to use
    
        - length: an integer with the target length for the paths.
    
        - useBonds: (optional) toggles the use of bond indices in the paths.
          Otherwise atom indices are used.  *Note* this behavior is different
          from that for subgraphs.
          Defaults to 1.
    
      RETURNS: a tuple of tuples with IDs for the bonds.
    
      NOTES: 
    
       - Difference between _subgraphs_ and _paths_ :: 
    
           Subgraphs are potentially branched, whereas paths (in our 
           terminology at least) cannot be.  So, the following graph: 
    
                C--0--C--1--C--3--C
                      |
                      2
                      |
                      C
    
           has 3 _subgraphs_ of length 3: (0,1,2),(0,1,3),(2,1,3)
           but only 2 _paths_ of length 3: (0,1,3),(2,1,3)
    
    

    C++ signature :
        std::list<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > FindAllPathsOfLengthN(RDKit::ROMol,unsigned int [,bool=True [,bool=False]])

FindAllSubgraphsOfLengthN(...)

source code 

FindAllSubgraphsOfLengthN( (Mol)mol, (int)length [, (bool)useHs=False]) -> _listSt6vectorIiSaIiEE :
    Finds all subgraphs of a particular length in a molecule
    
      ARGUMENTS:
    
        - mol: the molecule to use
    
        - length: an integer with the target number of bonds for the subgraphs.
    
        - useHs: (optional) toggles whether or not bonds to Hs that are part of the graph
          should be included in the results.
          Defaults to 0.
    
        - verbose: (optional, internal use) toggles verbosity in the search algorithm.
          Defaults to 0.
    
      RETURNS: a tuple of 2-tuples with bond IDs
    
      NOTES: 
    
       - Difference between _subgraphs_ and _paths_ :: 
    
           Subgraphs are potentially branched, whereas paths (in our 
           terminology at least) cannot be.  So, the following graph: 
    
                C--0--C--1--C--3--C
                      |
                      2
                      |
                      C
      has 3 _subgraphs_ of length 3: (0,1,2),(0,1,3),(2,1,3)
      but only 2 _paths_ of length 3: (0,1,3),(2,1,3)
    
    

    C++ signature :
        std::list<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > FindAllSubgraphsOfLengthN(RDKit::ROMol,unsigned int [,bool=False])

FindUniqueSubgraphsOfLengthN(...)

source code 

FindUniqueSubgraphsOfLengthN( (Mol)mol, (int)length [, (bool)useHs=False [, (bool)useBO=True]]) -> _listSt6vectorIiSaIiEE :
    Finds unique subgraphs of a particular length in a molecule
    
      ARGUMENTS:
    
        - mol: the molecule to use
    
        - length: an integer with the target number of bonds for the subgraphs.
    
        - useHs: (optional) toggles whether or not bonds to Hs that are part of the graph
          should be included in the results.
          Defaults to 0.
    
        - useBO: (optional) Toggles use of bond orders in distinguishing one subgraph from
          another.
          Defaults to 1.
    
      RETURNS: a tuple of tuples with bond IDs
    
    
    

    C++ signature :
        std::list<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > FindUniqueSubgraphsOfLengthN(RDKit::ROMol,unsigned int [,bool=False [,bool=True]])

GetAdjacencyMatrix(...)

source code 

GetAdjacencyMatrix( (Mol)mol [, (bool)useBO=False [, (int)emptyVal=0 [, (bool)force=False [, (str)prefix='']]]]) -> object :
    Returns the molecule's adjacency matrix.
    
      ARGUMENTS:
    
        - mol: the molecule to use
    
        - useBO: (optional) toggles use of bond orders in calculating the matrix.
          Default value is 0.
    
        - emptyVal: (optional) sets the elements of the matrix between non-adjacent atoms
          Default value is 0.
    
        - force: (optional) forces the calculation to proceed, even if there is a cached value.
          Default value is 0.
    
        - prefix: (optional, internal use) sets the prefix used in the property cache
          Default value is .
    
      RETURNS: a Numeric array of floats containing the adjacency matrix
    
    

    C++ signature :
        _object* GetAdjacencyMatrix(RDKit::ROMol {lvalue} [,bool=False [,int=0 [,bool=False [,char const*='']]]])

GetAlignmentTransform(...)

source code 

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 iteration 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]]]]]])

GetAtomMatch(...)

source code 

GetAtomMatch( (AtomPairsParameters)featMatch [, (int)maxAts=1024]) -> object :
    Returns an empty list if any of the features passed in share an atom.
     Otherwise a list of lists of atom indices is returned.
    

    C++ signature :
        boost::python::api::object GetAtomMatch(boost::python::api::object [,int=1024])

GetDistanceMatrix(...)

source code 

GetDistanceMatrix( (Mol)mol [, (bool)useBO=False [, (bool)useAtomWts=False [, (bool)force=False [, (str)prefix='']]]]) -> object :
    Returns the molecule's distance matrix.
    
      ARGUMENTS:
    
        - mol: the molecule to use
    
        - useBO: (optional) toggles use of bond orders in calculating the distance matrix.
          Default value is 0.
    
        - useAtomWts: (optional) toggles using atom weights for the diagonal elements of the
          matrix (to return a "Balaban" distance matrix).
          Default value is 0.
    
        - force: (optional) forces the calculation to proceed, even if there is a cached value.
          Default value is 0.
    
        - prefix: (optional, internal use) sets the prefix used in the property cache
          Default value is .
    
      RETURNS: a Numeric array of floats with the distance matrix
    
    

    C++ signature :
        _object* GetDistanceMatrix(RDKit::ROMol {lvalue} [,bool=False [,bool=False [,bool=False [,char const*='']]]])

GetFormalCharge(...)

source code 

GetFormalCharge( (Mol)arg1) -> int :
    Returns the formal charge for the molecule.
    
      ARGUMENTS:
    
        - mol: the molecule to use
    
    

    C++ signature :
        int GetFormalCharge(RDKit::ROMol)

GetMolFrags(...)

source code 

GetMolFrags( (Mol)mol [, (bool)asMols=False]) -> tuple :
    Finds the disconnected fragments from a molecule.
    
      For example, for the molecule 'CC(=O)[O-].[NH3+]C' GetMolFrags() returns
      ((0, 1, 2, 3), (4, 5))
    
      ARGUMENTS:
    
        - mol: the molecule to use
        - asMols: (optional) if this is provided and true, the fragments
          will be returned as molecules instead of atom ids.
    
      RETURNS: a tuple of tuples with IDs for the atoms in each fragment
               or a tuple of molecules.
    
    

    C++ signature :
        boost::python::tuple GetMolFrags(RDKit::ROMol [,bool=False])

GetMoleculeBoundsMatrix(...)

source code 

GetMoleculeBoundsMatrix( (Mol)mol [, (bool)set15bounds=True [, (bool)scaleVDW=False]]) -> object :
    Returns the distance bounds matrix for a molecule
     
     ARGUMENTS:
    
        - mol : the molecule of interest
        - set15bounds : set bounds for 1-5 atom distances based on 
                        topology (otherwise stop at 1-4s)
        - scaleVDW : scale down the sum of VDW radii when setting the 
                     lower bounds for atoms less than 5 bonds apart 
     RETURNS:
    
        the bounds matrix as a Numeric array with lower bounds in 
        the lower triangle and upper bounds in the upper triangle
    
    

    C++ signature :
        _object* GetMoleculeBoundsMatrix(RDKit::ROMol {lvalue} [,bool=True [,bool=False]])

GetPeriodicTable()

source code 
    Returns the application's PeriodicTable instance.
    
    

    C++ signature :
        RDKit::PeriodicTable* GetPeriodicTable()

Returns:
PeriodicTable :

GetSSSR(...)

source code 

GetSSSR( (Mol)arg1) -> int :
    Get the smallest set of simple rings for a molecule.
    
      ARGUMENTS:
    
        - mol: the molecule to use.
    
      RETURNS: the number of rings found
             This will be equal to NumBonds-NumAtoms+1 for single-fragment molecules.
    
    

    C++ signature :
        int GetSSSR(RDKit::ROMol {lvalue})

GetSymmSSSR(...)

source code 

GetSymmSSSR( (Mol)arg1) -> _vectSt6vectorIiSaIiEE :
    Get a symmetrized SSSR for a molecule.
    
      The symmetrized SSSR is at least as large as the SSSR for a molecule.
      In certain highly-symmetric cases (e.g. cubane), the symmetrized SSSR can be
      a bit larger (i.e. the number of symmetrized rings is >= NumBonds-NumAtoms+1).
    
      ARGUMENTS:
    
        - mol: the molecule to use.
    
      RETURNS: the number of rings found
    
    

    C++ signature :
        std::vector<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > GetSymmSSSR(RDKit::ROMol {lvalue})

Kekulize(...)

source code 

Kekulize( (Mol)mol [, (bool)clearAromaticFlags=False]) -> None :
    Kekulizes the molecule
    
      ARGUMENTS:
    
        - mol: the molecule to use
    
        - clearAromaticFlags: (optional) if this toggle is set, all atoms and bonds in the 
          molecule will be marked non-aromatic following the kekulization.
          Default value is 0.
    
      NOTES:
    
        - The molecule is modified in place.
    
    

    C++ signature :
        void Kekulize(RDKit::ROMol {lvalue} [,bool=False])

MolFromMolBlock(...)

source code 

MolFromMolBlock( (str)molBlock [, (bool)sanitize=True [, (bool)removeHs=True]]) -> Mol :
    Construct a molecule from a Mol block.
    
      ARGUMENTS:
    
        - molBlock: string containing the Mol block
    
        - sanitize: (optional) toggles sanitization of the molecule.
          Defaults to 1.
    
        - removeHs: (optional) toggles removing hydrogens from the molecule.
          This only make sense when sanitization is done.
          Defaults to true.
    
      RETURNS:
    
        a Mol object, None on failure.
    
    

    C++ signature :
        RDKit::ROMol* MolFromMolBlock(std::string [,bool=True [,bool=True]])

MolFromMolFile(...)

source code 

MolFromMolFile( (str)molFileName [, (bool)sanitize=True [, (bool)removeHs=True]]) -> Mol :
    Construct a molecule from a Mol file.
    
      ARGUMENTS:
    
        - fileName: name of the file to read
    
        - sanitize: (optional) toggles sanitization of the molecule.
          Defaults to true.
    
        - removeHs: (optional) toggles removing hydrogens from the molecule.
          This only make sense when sanitization is done.
          Defaults to true.
    
      RETURNS:
    
        a Mol object, None on failure.
    
    

    C++ signature :
        RDKit::ROMol* MolFromMolFile(char const* [,bool=True [,bool=True]])

MolFromSmarts(...)

source code 

MolFromSmarts( (str)SMARTS [, (bool)mergeHs=False]) -> Mol :
    Construct a molecule from a SMARTS string.
    
      ARGUMENTS:
    
        - SMARTS: the smarts string
    
        - mergeHs: (optional) toggles the merging of explicit Hs in the query into the attached
          atoms.  So, for example, 'C[H]' becomes '[C;!H0]'.
          Defaults to 0.
    
      RETURNS:
    
        a Mol object, None on failure.
    
    

    C++ signature :
        RDKit::ROMol* MolFromSmarts(char const* [,bool=False])

MolFromSmiles(...)

source code 

MolFromSmiles( (str)SMILES [, (bool)sanitize=True]) -> Mol :
    Construct a molecule from a SMILES string.
    
      ARGUMENTS:
    
        - SMILES: the smiles string
    
        - sanitize: (optional) toggles sanitization of the molecule.
          Defaults to 1.
    
      RETURNS:
    
        a Mol object, None on failure.
    
    

    C++ signature :
        RDKit::ROMol* MolFromSmiles(std::string [,bool=True])

MolFromTPLBlock(...)

source code 

MolFromTPLBlock( (str)tplBlock [, (bool)sanitize=True [, (bool)skipFirstConf=False]]) -> Mol :
    Construct a molecule from a TPL block.
    
      ARGUMENTS:
    
        - fileName: name of the file to read
    
        - sanitize: (optional) toggles sanitization of the molecule.
          Defaults to True.
    
        - skipFirstConf: (optional) skips reading the first conformer.
          Defaults to False.
          This should be set to True when reading TPLs written by 
          the CombiCode.
    
      RETURNS:
    
        a Mol object, None on failure.
    
    

    C++ signature :
        RDKit::ROMol* MolFromTPLBlock(std::string [,bool=True [,bool=False]])

MolFromTPLFile(...)

source code 

MolFromTPLFile( (str)fileName [, (bool)sanitize=True [, (bool)skipFirstConf=False]]) -> Mol :
    Construct a molecule from a TPL file.
    
      ARGUMENTS:
    
        - fileName: name of the file to read
    
        - sanitize: (optional) toggles sanitization of the molecule.
          Defaults to True.
    
        - skipFirstConf: (optional) skips reading the first conformer.
          Defaults to False.
          This should be set to True when reading TPLs written by 
          the CombiCode.
    
      RETURNS:
    
        a Mol object, None on failure.
    
    

    C++ signature :
        RDKit::ROMol* MolFromTPLFile(char const* [,bool=True [,bool=False]])

MolToMolBlock(...)

source code 

MolToMolBlock( (Mol)mol [, (bool)includeStereo=False [, (int)confId=-1]]) -> str :
    Returns the a Mol block for a molecule
      ARGUMENTS:
    
        - mol: the molecule
        - includeStereo: (optional) toggles inclusion of stereochemical
                         information in the output
        - confId: (optional) selects which conformation to output (-1 = default)
    
      RETURNS:
    
        a string
    
    

    C++ signature :
        std::string MolToMolBlock(RDKit::ROMol [,bool=False [,int=-1]])

MolToSmarts(...)

source code 

MolToSmarts( (Mol)mol [, (bool)isomericSmiles=False]) -> str :
    Returns a SMARTS string for a molecule
      ARGUMENTS:
    
        - mol: the molecule
        - isomericSmarts: (optional) include information about stereochemistry in
          the SMARTS.  Defaults to false.
    
      RETURNS:
    
        a string
    
    

    C++ signature :
        std::string MolToSmarts(RDKit::ROMol {lvalue} [,bool=False])

MolToSmiles(...)

source code 

MolToSmiles( (Mol)mol [, (bool)isomericSmiles=False [, (bool)kekuleSmiles=False [, (int)rootedAtAtom=-1]]]) -> str :
    Returns the canonical SMILES string for a molecule
      ARGUMENTS:
    
        - mol: the molecule
        - isomericSmiles: (optional) include information about stereochemistry in
          the SMILES.  Defaults to false.
        - kekuleSmiles: (optional) use the Kekule form (no aromatic bonds) in
          the SMILES.  Defaults to false.
        - rootedAtAtom: (optional) if non-negative, this forces the SMILES 
          to start at a particular atom. Defaults to -1.
    
      RETURNS:
    
        a string
    
    

    C++ signature :
        std::string MolToSmiles(RDKit::ROMol {lvalue} [,bool=False [,bool=False [,int=-1]]])

MolToTPLBlock(...)

source code 

MolToTPLBlock( (Mol)mol [, (str)partialChargeProp='_GasteigerCharge' [, (bool)writeFirstConfTwice=False]]) -> str :
    Returns the Tpl block for a molecule.
    
      ARGUMENTS:
    
        - mol: the molecule
        - partialChargeProp: name of the property to use for partial charges
          Defaults to '_GasteigerCharge'.
        - writeFirstConfTwice: Defaults to False.
          This should be set to True when writing TPLs to be read by 
          the CombiCode.
    
      RETURNS:
    
        a string
    
    

    C++ signature :
        std::string MolToTPLBlock(RDKit::ROMol [,std::string='_GasteigerCharge' [,bool=False]])

MolToTPLFile(...)

source code 

MolToTPLFile( (Mol)mol, (str)fileName [, (str)partialChargeProp='_GasteigerCharge' [, (bool)writeFirstConfTwice=False]]) -> None :
    Writes a molecule to a TPL file.
    
      ARGUMENTS:
    
        - mol: the molecule
        - fileName: name of the file to write
        - partialChargeProp: name of the property to use for partial charges
          Defaults to '_GasteigerCharge'.
        - writeFirstConfTwice: Defaults to False.
          This should be set to True when writing TPLs to be read by 
          the CombiCode.
    
    

    C++ signature :
        void MolToTPLFile(RDKit::ROMol,std::string [,std::string='_GasteigerCharge' [,bool=False]])

RDKFingerprint(...)

source code 

RDKFingerprint( (Mol)mol [, (int)minPath=1 [, (int)maxPath=7 [, (int)fpSize=2048 [, (int)nBitsPerHash=4 [, (bool)useHs=True [, (float)tgtDensity=0.0 [, (int)minSize=128]]]]]]]) -> ExplicitBitVect :
    Returns an RDKit topological fingerprint for a molecule
    
      Explanation of the algorithm below.
    
      ARGUMENTS:
    
        - mol: the molecule to use
    
        - minPath: (optional) minimum number of bonds to include in the subgraphs
          Defaults to 1.
    
        - maxPath: (optional) maximum number of bonds to include in the subgraphs
          Defaults to 7.
    
        - fpSize: (optional) number of bits in the fingerprint
          Defaults to 2048.
    
        - nBitsPerPath: (optional) number of bits to set per path
          Defaults to 4.
    
        - useHs: (optional) include information about number of Hs on each
          atom when calculating path hashes.
          Defaults to 1.
    
        - tgtDensity: (optional) fold the fingerprint until this minimum density has
          been reached
          Defaults to 0.
    
        - minSize: (optional) the minimum size the fingerprint will be folded to when
          trying to reach tgtDensity
          Defaults to 128.
    
      RETURNS: a DataStructs.ExplicitBitVect with _fpSize_ bits
    
      ALGORITHM:
    
       This algorithm functions by find all paths between minPath and maxPath in
        length.  For each path:
    
         1) A hash is calculated.
    
         2) The hash is used to seed a random-number generator
    
         3) _nBitsPerPath_ random numbers are generated and used to set the corresponding
            bits in the fingerprint
    
    
    

    C++ signature :
        ExplicitBitVect* RDKFingerprint(RDKit::ROMol [,unsigned int=1 [,unsigned int=7 [,unsigned int=2048 [,unsigned int=4 [,bool=True [,double=0.0 [,unsigned int=128]]]]]]])

ReactionFromRxnBlock(...)

source code 

ReactionFromRxnBlock( (str)arg1) -> ChemicalReaction :
    construct a ChemicalReaction from an string in MDL rxn format

    C++ signature :
        RDKit::ChemicalReaction* ReactionFromRxnBlock(std::string)

ReactionFromRxnFile(...)

source code 

ReactionFromRxnFile( (str)arg1) -> ChemicalReaction :
    construct a ChemicalReaction from an MDL rxn file

    C++ signature :
        RDKit::ChemicalReaction* ReactionFromRxnFile(std::string)

ReactionFromSmarts(...)

source code 

ReactionFromSmarts( (str)arg1) -> ChemicalReaction :
    construct a ChemicalReaction from a reaction SMARTS string

    C++ signature :
        RDKit::ChemicalReaction* ReactionFromSmarts(std::string)

RemoveHs(...)

source code 

RemoveHs( (Mol)mol [, (bool)implicitOnly=False]) -> Mol :
    Removes any hydrogens from the graph of a molecule.
    
      ARGUMENTS:
    
        - mol: the molecule to be modified
    
        - implicitOnly: (optional) if this toggle is set, only implicit Hs will
          be removed from the graph.  Default value is 0 (remove implicit and explicit Hs).
    
      RETURNS: a new molecule with the Hs removed
    
      NOTES:
    
        - The original molecule is *not* modified.
    
    

    C++ signature :
        RDKit::ROMol* RemoveHs(RDKit::ROMol [,bool=False])

ReplaceCore(...)

source code 

ReplaceCore( (Mol)mol, (Mol)coreQuery [, (bool)replaceDummies=True]) -> Mol :
    Removes the core of a molecule and labels the sidechains with dummy atoms.
    
      ARGUMENTS:
    
        - mol: the molecule to be modified
    
        - coreQuery: the molecule to be used as a substructure query for recognizing the core
    
        - replaceDummies: toggles replacement of atoms that match dummies in the query
    
      RETURNS: a new molecule with the core removed
    
      NOTES:
    
        - The original molecule is *not* modified.
    
      EXAMPLES:
    
       The following examples substitute SMILES/SMARTS strings for molecules, you'd have
       to actually use molecules:
    
        - ReplaceCore('CCC1CCC1','C1CCC1') -> 'CC[Xa]'
    
        - ReplaceCore('CCC1CC1','C1CCC1') -> ''
    
        - ReplaceCore('C1CC2C1CCC2','C1CCC1') -> '[Xa]C1CCC1[Xb]'
    
        - ReplaceCore('C1CNCC1','N') -> '[Xa]CCCC[Xb]'
    
        - ReplaceCore('C1CCC1CN','C1CCC1[*]',False) -> '[Xa]CN'
    
    

    C++ signature :
        RDKit::ROMol* ReplaceCore(RDKit::ROMol,RDKit::ROMol [,bool=True])

ReplaceSidechains(...)

source code 

ReplaceSidechains( (Mol)mol, (Mol)coreQuery) -> Mol :
    Replaces sidechains in a molecule with dummy atoms for their attachment points.
    
      ARGUMENTS:
    
        - mol: the molecule to be modified
    
        - coreQuery: the molecule to be used as a substructure query for recognizing the core
    
      RETURNS: a new molecule with the sidechains removed
    
      NOTES:
    
        - The original molecule is *not* modified.
    
      EXAMPLES:
    
       The following examples substitute SMILES/SMARTS strings for molecules, you'd have
       to actually use molecules:
    
        - ReplaceSidechains('CCC1CCC1','C1CCC1') -> '[Xa]C1CCC1'
    
        - ReplaceSidechains('CCC1CC1','C1CCC1') -> ''
    
        - ReplaceSidechains('C1CC2C1CCC2','C1CCC1') -> '[Xa]C1CCC1[Xb]'
    
    

    C++ signature :
        RDKit::ROMol* ReplaceSidechains(RDKit::ROMol,RDKit::ROMol)

ReplaceSubstructs(...)

source code 

ReplaceSubstructs( (Mol)mol, (Mol)query, (Mol)replacement [, (bool)replaceAll=False]) -> object :
    Replaces atoms matching a substructure query in a molecule
    
      ARGUMENTS:
    
        - mol: the molecule to be modified
    
        - query: the molecule to be used as a substructure query
    
        - replacement: the molecule to be used as the replacement
    
        - replaceAll: (optional) if this toggle is set, all substructures matching
          the query will be replaced in a single result, otherwise each result will
          contain a separate replacement.
          Default value is False (return multiple replacements)
    
      RETURNS: a tuple of new molecules with the substructures replaced removed
    
      NOTES:
    
        - The original molecule is *not* modified.
    
      EXAMPLES:
    
       The following examples substitute SMILES/SMARTS strings for molecules, you'd have
       to actually use molecules:
    
        - ReplaceSubstructs('CCOC','OC','NC') -> ('CCNC',)
    
        - ReplaceSubstructs('COCCOC','OC','NC') -> ('COCCNC','CNCCOC')
    
        - ReplaceSubstructs('COCCOC','OC','NC',True) -> ('CNCCNC',)
    
    

    C++ signature :
        _object* ReplaceSubstructs(RDKit::ROMol,RDKit::ROMol,RDKit::ROMol [,bool=False])

SanitizeMol(...)

source code 

SanitizeMol( (Mol)arg1) -> None :
    Kekulize, check valencies, set aromaticity, conjugation and hybridization
    
        - The molecule is modified in place.
    
        - If sanitization fails, an exception will be thrown
    
      ARGUMENTS:
    
        - mol: the molecule to be modified
    
      NOTES:
    
    

    C++ signature :
        void SanitizeMol(RDKit::ROMol {lvalue})

ShapeProtrudeDist(...)

source code 

ShapeProtrudeDist( (Mol)mol1, (Mol)mol2 [, (int)confId1=-1 [, (int)confId2=-1 [, (float)gridSpacing=0.5 [, (DiscreteValueType)bitsPerPoint=DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE [, (float)vdwScale=0.80000000000000004 [, (float)stepSize=0.25 [, (int)maxLayers=-1 [, (bool)ignoreHs=True [, (bool)allowReordering=True]]]]]]]]]) -> float :
    Compute the shape protrude distance between two molecule based on a predefined alignment
      
      ARGUMENTS:
        - mol1 : The first molecule of interest 
        - mol2 : The second molecule of interest 
        - confId1 : Conformer in the first molecule (defaults to first conformer) 
        - confId2 : Conformer in the second molecule (defaults to first conformer) 
        - gridSpacing : resolution of the grid used to encode the molecular shapes 
        - bitsPerPoint : number of bit used to encode the occupancy at each grid point 
                              defaults to two bits per grid point 
        - vdwScale : Scaling factor for the radius of the atoms to determine the base radius 
                    used in the encoding - grid points inside this sphere carry the maximum occupan 
        - stepSize : thickness of the each layer outside the base radius, the occupancy value is decreased 
                     from layer to layer from the maximum value 
        - maxLayers : the maximum number of layers - defaults to the number allowed the number of bits 
                      use per grid point - e.g. two bits per grid point will allow 3 layers 
        - ignoreHs : when set, the contribution of Hs to the shape will be ignored
        - allowReordering : when set, the order will be automatically updated so that the value calculated
                            is the protrusion of the smaller shape from the larger one.
    

    C++ signature :
        double ShapeProtrudeDist(RDKit::ROMol,RDKit::ROMol [,int=-1 [,int=-1 [,double=0.5 [,RDKit::DiscreteValueVect::DiscreteValueType=DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE [,double=0.80000000000000004 [,double=0.25 [,int=-1 [,bool=True [,bool=True]]]]]]]]])

ShapeTanimotoDist(...)

source code 

ShapeTanimotoDist( (Mol)mol1, (Mol)mol2 [, (int)confId1=-1 [, (int)confId2=-1 [, (float)gridSpacing=0.5 [, (DiscreteValueType)bitsPerPoint=DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE [, (float)vdwScale=0.80000000000000004 [, (float)stepSize=0.25 [, (int)maxLayers=-1 [, (bool)ignoreHs=True]]]]]]]]) -> float :
    Compute the shape tanimoto distance between two molecule based on a predefined alignment
      
      ARGUMENTS:
        - mol1 : The first molecule of interest 
        - mol2 : The second molecule of interest 
        - confId1 : Conformer in the first molecule (defaults to first conformer) 
        - confId2 : Conformer in the second molecule (defaults to first conformer) 
        - gridSpacing : resolution of the grid used to encode the molecular shapes 
        - bitsPerPoint : number of bit used to encode the occupancy at each grid point 
                              defaults to two bits per grid point 
        - vdwScale : Scaling factor for the radius of the atoms to determine the base radius 
                    used in the encoding - grid points inside this sphere carry the maximum occupan 
        - stepSize : thickness of the each layer outside the base radius, the occupancy value is decreased 
                     from layer to layer from the maximum value 
        - maxLayers : the maximum number of layers - defaults to the number allowed the number of bits 
                      use per grid point - e.g. two bits per grid point will allow 3 layers 
        - ignoreHs : when set, the contribution of Hs to the shape will be ignored
    

    C++ signature :
        double ShapeTanimotoDist(RDKit::ROMol,RDKit::ROMol [,int=-1 [,int=-1 [,double=0.5 [,RDKit::DiscreteValueVect::DiscreteValueType=DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE [,double=0.80000000000000004 [,double=0.25 [,int=-1 [,bool=True]]]]]]]])

SmilesMolSupplierFromText(...)

source code 

SmilesMolSupplierFromText( (str)text [, (str)delimiter=' ' [, (int)smilesColumn=0 [, (int)nameColumn=1 [, (bool)titleLine=True [, (bool)sanitize=True]]]]]) -> SmilesMolSupplier :

    C++ signature :
        RDKit::SmilesMolSupplier* SmilesMolSupplierFromText(std::string [,std::string=' ' [,int=0 [,int=1 [,bool=True [,bool=True]]]]])

TransformConformer(...)

source code 

TransformConformer( (Conformer)arg1, (AtomPairsParameters)arg2) -> None :
    Transform the coordinates of a conformer

    C++ signature :
        void TransformConformer(RDKit::Conformer {lvalue},boost::python::api::object)

UFFGetMoleculeForceField(...)

source code 

UFFGetMoleculeForceField( (Mol)mol [, (float)vdwThresh=10.0 [, (int)confId=-1]]) -> ForceField :
    returns a UFF force field for a molecule
    
     
     ARGUMENTS:
    
        - mol : the molecule of interrest
        - vdwThresh : used to exclude long-range van der Waals interactions
                      (defaults to 10.0)
        - confId : indicates which conformer to optimize
    
    

    C++ signature :
        ForceFields::PyForceField* UFFGetMoleculeForceField(RDKit::ROMol {lvalue} [,double=10.0 [,int=-1]])

UFFOptimizeMolecule(...)

source code 

UFFOptimizeMolecule( (Mol)self [, (int)maxIters=200 [, (float)vdwThresh=10.0 [, (int)confId=-1]]]) -> int :
    Use UFF to optimize a molecule's structure
    
     
     ARGUMENTS:
    
        - mol : the molecule of interrest
        - maxIters : the maximum number of iterations (defaults to 100)
        - vdwThresh : used to exclude long-range van der Waals interactions
                      (defaults to 10.0)
        - confId : indicates which conformer to optimize
    
    

    C++ signature :
        int UFFOptimizeMolecule(RDKit::ROMol {lvalue} [,int=200 [,double=10.0 [,int=-1]]])

WedgeMolBonds(...)

source code 

WedgeMolBonds( (Mol)arg1, (Conformer)arg2) -> None :
    Set the wedging on single bonds in a molecule.
        The wedging scheme used is that from Mol files.
     
      ARGUMENTS:
    
        - molecule: the molecule to update
     
    
    

    C++ signature :
        void WedgeMolBonds(RDKit::ROMol {lvalue},RDKit::Conformer const*)

tossit()

source code 
    C++ signature :
        void tossit()

Returns:
None :