rdkit.Chem.rdmolops module

Module containing RDKit functionality for manipulating molecules.

rdkit.Chem.rdmolops.AddHs((Mol)mol[, (bool)explicitOnly=False[, (bool)addCoords=False[, (AtomPairsParameters)onlyOnAtoms=None[, (bool)addResidueInfo=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).

  • onlyOnAtoms: (optional) if this sequence is provided, only these atoms will be considered to have Hs added to them

  • addResidueInfo: (optional) if this is true, add residue info to hydrogen atoms (useful for PDB files).

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 [,boost::python::api::object=None [,bool=False]]]])

rdkit.Chem.rdmolops.AddRecursiveQuery((Mol)mol, (Mol)query, (int)atomIdx[, (bool)preserveExistingQuery=True]) None :

Adds a recursive query to an atom

ARGUMENTS:

  • mol: the molecule to be modified

  • query: the molecule to be used as the recursive query (this will be copied)

  • atomIdx: the atom to modify

  • preserveExistingQuery: (optional) if this is set, existing query information on the atom will be preserved

RETURNS: None

C++ signature :

void AddRecursiveQuery(RDKit::ROMol {lvalue},RDKit::ROMol,unsigned int [,bool=True])

rdkit.Chem.rdmolops.AddStereoAnnotations((Mol)mol[, (str)absLabel='abs ({cip})'[, (str)orLabel='or{id}'[, (str)andLabel='and{id}'[, (str)cipLabel='({cip})'[, (str)bondLabel='({cip})']]]]]) None :

add R/S, relative stereo, and E/Z annotations to atoms and bonds

Arguments:
  • mol: molecule to modify

  • absLabel: label for atoms in an ABS stereo group

  • orLabel: label for atoms in an OR stereo group

  • andLabel: label for atoms in an AND stereo group

  • cipLabel: label for chiral atoms that aren’t in a stereo group.

  • bondLabel: label for CIP stereochemistry on bonds

If any label is empty, the corresponding annotations will not be added.

The labels can contain the following placeholders:
  • {id} - the stereo group’s index

  • {cip} - the atom or bond’s CIP stereochemistry

Note that CIP labels will only be added if CIP stereochemistry has been assigned to the molecule.

C++ signature :

void AddStereoAnnotations(RDKit::ROMol {lvalue} [,std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> >=’abs ({cip})’ [,std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> >=’or{id}’ [,std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> >=’and{id}’ [,std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> >=’({cip})’ [,std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> >=’({cip})’]]]]])

rdkit.Chem.rdmolops.AddWavyBondsForStereoAny((Mol)mol[, (bool)clearDoubleBondFlags=True[, (int)addWhenImpossible=1000]]) None :
set wavy bonds around double bonds with STEREOANY stereo
ARGUMENTS :
  • molecule : the molecule to updaten -

  • conformer : the conformer to use to determine wedge direction

C++ signature :

void AddWavyBondsForStereoAny(RDKit::ROMol {lvalue} [,bool=True [,unsigned int=1000]])

class rdkit.Chem.rdmolops.AdjustQueryParameters((object)arg1)

Bases: instance

Parameters controlling which components of the query atoms/bonds are adjusted.

Note that some of the options here are either directly contradictory or make

no sense when combined with each other. We generally assume that client code is doing something sensible and don’t attempt to detect possible conflicts or problems.

A note on the flags controlling which atoms/bonds are modified:

These generally limit the set of atoms/bonds to be modified. For example:

  • ADJUST_IGNORERINGS atoms/bonds in rings will not be modified.

  • ADJUST_IGNORENONE causes all atoms/bonds to be modified

  • ADJUST_IGNOREALL no atoms/bonds will be modified

Some of the options obviously make no sense for bonds

C++ signature :

void __init__(_object*)

static NoAdjustments() AdjustQueryParameters :

Returns an AdjustQueryParameters object with all parameters set to false

C++ signature :

RDKit::MolOps::AdjustQueryParameters NoAdjustments()

property adjustConjugatedFiveRings

set bond queries in conjugated five-rings to SINGLE|DOUBLE|AROMATIC

property adjustDegree

add degree queries

property adjustDegreeFlags

controls which atoms have their degree queries changed

property adjustHeavyDegree

adjust the heavy-atom degree

property adjustHeavyDegreeFlags

controls which atoms have their heavy-atom degree queries changed

property adjustRingChain

add ring-chain queries to atoms

property adjustRingChainFlags

controls which atoms have ring-chain queries added

property adjustRingCount

add ring-count queries

property adjustRingCountFlags

controls which atoms have ring-count queries added

property adjustSingleBondsBetweenAromaticAtoms

sets non-ring single bonds between two aromatic or conjugated atoms to SINGLE|AROMATIC

property adjustSingleBondsToDegreeOneNeighbors

set single bonds bewteen aromatic or conjugated atoms and degree-one neighbors to SINGLE|AROMATIC

property aromatizeIfPossible

perceive and set aromaticity

property makeAtomsGeneric

convert atoms to generic queries (any atoms)

property makeAtomsGenericFlags

controls which atoms are converted to generic queries

property makeBondsGeneric

converts bonds to generic queries (any bonds)

property makeBondsGenericFlags

controls which bonds are converted to generic queries

property makeDummiesQueries

convert dummy atoms without isotope labels to any-atom queries

property setMDLFiveRingAromaticity

uses the 5-ring aromaticity behavior of the (former) MDL software as documented in the Chemical Representation Guide

property useStereoCareForBonds

if this is set sterochemistry information will be removed from double bonds that do not have the stereoCare property set

rdkit.Chem.rdmolops.AdjustQueryProperties((Mol)mol[, (AtomPairsParameters)params=None]) Mol :

Returns a new molecule where the query properties of atoms have been modified.

C++ signature :

RDKit::ROMol* AdjustQueryProperties(RDKit::ROMol [,boost::python::api::object=None])

rdkit.Chem.rdmolops.AdjustQueryPropertiesWithGenericGroups((Mol)mol[, (AtomPairsParameters)params=None]) Mol :

Returns a new molecule where the query properties of atoms have been modified and generic group queries have been prepared.

C++ signature :

RDKit::ROMol* AdjustQueryPropertiesWithGenericGroups(RDKit::ROMol [,boost::python::api::object=None])

class rdkit.Chem.rdmolops.AdjustQueryWhichFlags

Bases: enum

ADJUST_IGNOREALL = rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNOREALL
ADJUST_IGNORECHAINS = rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORECHAINS
ADJUST_IGNOREDUMMIES = rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNOREDUMMIES
ADJUST_IGNOREMAPPED = rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNOREMAPPED
ADJUST_IGNORENONDUMMIES = rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORENONDUMMIES
ADJUST_IGNORENONE = rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORENONE
ADJUST_IGNORERINGS = rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORERINGS
names = {'ADJUST_IGNOREALL': rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNOREALL, 'ADJUST_IGNORECHAINS': rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORECHAINS, 'ADJUST_IGNOREDUMMIES': rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNOREDUMMIES, 'ADJUST_IGNOREMAPPED': rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNOREMAPPED, 'ADJUST_IGNORENONDUMMIES': rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORENONDUMMIES, 'ADJUST_IGNORENONE': rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORENONE, 'ADJUST_IGNORERINGS': rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORERINGS}
values = {0: rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORENONE, 1: rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORECHAINS, 2: rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNOREDUMMIES, 4: rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORERINGS, 8: rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNORENONDUMMIES, 16: rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNOREMAPPED, 268435455: rdkit.Chem.rdmolops.AdjustQueryWhichFlags.ADJUST_IGNOREALL}
class rdkit.Chem.rdmolops.AromaticityModel

Bases: enum

AROMATICITY_CUSTOM = rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_CUSTOM
AROMATICITY_DEFAULT = rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_DEFAULT
AROMATICITY_MDL = rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_MDL
AROMATICITY_RDKIT = rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_RDKIT
AROMATICITY_SIMPLE = rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_SIMPLE
names = {'AROMATICITY_CUSTOM': rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_CUSTOM, 'AROMATICITY_DEFAULT': rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_DEFAULT, 'AROMATICITY_MDL': rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_MDL, 'AROMATICITY_RDKIT': rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_RDKIT, 'AROMATICITY_SIMPLE': rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_SIMPLE}
values = {0: rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_DEFAULT, 1: rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_RDKIT, 2: rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_SIMPLE, 4: rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_MDL, 268435455: rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_CUSTOM}
rdkit.Chem.rdmolops.AssignAtomChiralTagsFromMolParity((Mol)mol[, (bool)replaceExistingTags=True]) None :
Sets the chiral tags on a molecule’s atoms based on

the molParity atom property.

ARGUMENTS:

  • mol: the molecule to use

  • replaceExistingTags: if True, existing stereochemistry information will be cleared

before running the calculation.

C++ signature :

void AssignAtomChiralTagsFromMolParity(RDKit::ROMol {lvalue} [,bool=True])

rdkit.Chem.rdmolops.AssignAtomChiralTagsFromStructure((Mol)mol[, (int)confId=-1[, (bool)replaceExistingTags=True]]) None :
Sets the chiral tags on a molecule’s atoms based on

a 3D conformation. NOTE that this does not check to see if atoms are chiral centers (i.e. all substituents are different), it merely sets the chiral type flags based on the coordinates and atom ordering. Use AssignStereochemistryFrom3D() if you want chiral flags only on actual stereocenters.

ARGUMENTS:

  • mol: the molecule to use

  • confId: the conformer id to use, -1 for the default

  • replaceExistingTags: if True, existing stereochemistry information will be cleared

before running the calculation.

C++ signature :

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

rdkit.Chem.rdmolops.AssignChiralTypesFromBondDirs((Mol)mol[, (int)confId=-1[, (bool)replaceExistingTags=True]]) None :

Uses bond directions to assign ChiralTypes to a molecule’s atoms.

ARGUMENTS:

  • mol: the molecule to use

  • confId: (optional) the conformation to use

  • replaceExistingTags: (optional) replace any existing information about stereochemistry

C++ signature :

void AssignChiralTypesFromBondDirs(RDKit::ROMol {lvalue} [,int=-1 [,bool=True]])

rdkit.Chem.rdmolops.AssignRadicals((Mol)mol) None :

Assigns radical counts to atoms

ARGUMENTS:

  • mol: the molecule to use

NOTES:

  • The molecule is modified in place.

C++ signature :

void AssignRadicals(RDKit::ROMol {lvalue})

rdkit.Chem.rdmolops.AssignStereochemistry((Mol)mol[, (bool)cleanIt=False[, (bool)force=False[, (bool)flagPossibleStereoCenters=False]]]) None :
Assign stereochemistry tags to atoms and bonds.

If useLegacyStereoPerception is true, it also does the CIP stereochemistry assignment for the molecule’s atoms (R/S) and double bonds (Z/E). This assignment is based on legacy code which is fast, but is known to incorrectly assign CIP labels in some cases. instead, to assign CIP labels based on an accurate, though slower, implementation of the CIP rules, call CIPLabeler::assignCIPLabels(). Chiral atoms will have a property ‘_CIPCode’ indicating their chiral code.

ARGUMENTS:

  • mol: the molecule to use

  • cleanIt: (optional) if provided, any existing values of the property _CIPCode

    will be cleared, 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. Bonds with STEREOCIS/STEREOTRANS specified that have duplicate substituents based upon the CIP atom ranks will be marked STEREONONE.

  • force: (optional) causes the calculation to be repeated, even if it has already been done

  • flagPossibleStereoCenters (optional) set the _ChiralityPossible property on atoms that are possible stereocenters

C++ signature :

void AssignStereochemistry(RDKit::ROMol {lvalue} [,bool=False [,bool=False [,bool=False]]])

rdkit.Chem.rdmolops.AssignStereochemistryFrom3D((Mol)mol[, (int)confId=-1[, (bool)replaceExistingTags=True]]) None :
Uses a conformer (should be 3D) to assign ChiralTypes to a molecule’s atoms

and stereo flags to its bonds

ARGUMENTS:

  • mol: the molecule to use

  • confId: (optional) the conformation to use

  • replaceExistingTags: (optional) replace any existing information about stereochemistry

C++ signature :

void AssignStereochemistryFrom3D(RDKit::ROMol {lvalue} [,int=-1 [,bool=True]])

class rdkit.Chem.rdmolops.BondWedgingParameters((object)arg1)

Bases: instance

Parameters controlling how bond wedging is done.

C++ signature :

void __init__(_object*)

property wedgeTwoBondsIfPossible

If this is enabled then two bonds will be wedged at chiral centers subject to the following constraints:

  1. ring bonds will not be wedged

  2. bonds to chiral centers will not be wedged

  3. bonds separated by more than 120 degrees will not be

    wedged

rdkit.Chem.rdmolops.Cleanup((Mol)mol) None :

cleans up certain common bad functionalities in the molecule

ARGUMENTS:

  • mol: the molecule to use

NOTES:

  • The molecule is modified in place.

C++ signature :

void Cleanup(RDKit::ROMol {lvalue})

rdkit.Chem.rdmolops.CleanupOrganometallics((Mol)mol) None :

cleans up certain common bad functionalities in the organometallic molecule

Note that this function is experimental and may either change in behavior or be replaced with something else in future releases.

ARGUMENTS :

  • mol : the molecule to use

NOTES :

  • The molecule is modified in place.

C++ signature :

void CleanupOrganometallics(RDKit::ROMol {lvalue})

rdkit.Chem.rdmolops.CollapseAttachmentPoints((Mol)mol[, (bool)markedOnly=True]) None :

dummy atoms in the graph are removed and replaced with attachment point annotations on the attached atoms

Arguments:
  • mol: molecule to be modified

  • markedOnly: if true, only dummy atoms with the _fromAttachPoint property will be collapsed

In order for a dummy atom to be considered for collapsing it must have:
  • degree 1 with a single or unspecified bond

  • the bond to it can not be wedged

  • either no query or be an AtomNullQuery

C++ signature :

void CollapseAttachmentPoints(RDKit::ROMol {lvalue} [,bool=True])

rdkit.Chem.rdmolops.CombineMols((Mol)mol1, (Mol)mol2[, (Point3D)offset=<rdkit.Geometry.rdGeometry.Point3D object at 0x746192768d40>]) Mol :

Combine the atoms from two molecules to produce a third

C++ signature :

RDKit::ROMol* CombineMols(RDKit::ROMol,RDKit::ROMol [,RDGeom::Point3D=<rdkit.Geometry.rdGeometry.Point3D object at 0x746192768d40>])

rdkit.Chem.rdmolops.ConvertGenericQueriesToSubstanceGroups((Mol)mol) None :

documentation

C++ signature :

void ConvertGenericQueriesToSubstanceGroups(RDKit::ROMol {lvalue})

rdkit.Chem.rdmolops.DativeBondsToHaptic((Mol)mol) Mol :

Does the reverse of hapticBondsToDative. If there are multiple contiguous atoms attached by dative bonds to an atom (probably a metal atom), the dative bonds will be replaced by a dummy atom in their centre attached to the (metal) atom by a dative bond, which is labelled with ENDPTS of the atoms that had the original dative bonds.

ARGUMENTS:

  • mol: the molecule to use

RETURNS:

a modified copy of the molecule

C++ signature :

RDKit::ROMol* DativeBondsToHaptic(RDKit::ROMol)

rdkit.Chem.rdmolops.DeleteSubstructs((Mol)mol, (Mol)query[, (bool)onlyFrags=False[, (bool)useChirality=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)

  • useChirality: (optional) match the substructure query using chirality

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 [,bool=False]])

rdkit.Chem.rdmolops.DetectBondStereoChemistry((Mol)mol, (Conformer)conformer) None :
Assign stereochemistry to bonds based on coordinates and a conformer.

DEPRECATED

ARGUMENTS:

  • mol: the molecule to be modified

  • conformer: Conformer providing the coordinates

C++ signature :

void DetectBondStereoChemistry(RDKit::ROMol {lvalue},RDKit::Conformer const*)

rdkit.Chem.rdmolops.DetectBondStereochemistry((Mol)mol[, (int)confId=-1]) None :
  • mol: the molecule to be modified

  • confId: Conformer to use for the coordinates

C++ signature :

void DetectBondStereochemistry(RDKit::ROMol {lvalue} [,int=-1])

rdkit.Chem.rdmolops.DetectChemistryProblems((Mol)mol[, (int)sanitizeOps=rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_ALL]) tuple :

checks for chemistry problems

C++ signature :

boost::python::tuple DetectChemistryProblems(RDKit::ROMol [,unsigned int=rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_ALL])

rdkit.Chem.rdmolops.ExpandAttachmentPoints((Mol)mol[, (bool)addAsQueries=True[, (bool)addCoords=True]]) None :

attachment points encoded as attachPt properties are added to the graph as dummy atoms

Arguments:
  • mol: molecule to be modified

  • addAsQueries: if true, the dummy atoms will be added as null queries

    (i.e. they will match any atom in a substructure search)

  • addCoords: if true and the molecule has one or more conformers,

    positions for the attachment points will be added to the conformer(s)

C++ signature :

void ExpandAttachmentPoints(RDKit::ROMol {lvalue} [,bool=True [,bool=True]])

rdkit.Chem.rdmolops.FastFindRings((Mol)mol) None :

Does a non-SSSR ring finding for a molecule.

ARGUMENTS:

  • mol: the molecule to use.

RETURNS: Nothing

C++ signature :

void FastFindRings(RDKit::ROMol)

rdkit.Chem.rdmolops.FindAllPathsOfLengthN((Mol)mol, (int)length[, (bool)useBonds=True[, (bool)useHs=False[, (int)rootedAtAtom=-1[, (bool)onlyShortestPaths=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.

  • rootedAtAtom: (optional) if nonzero, only paths from the specified atom will be returned.

  • onlyShortestPaths: (optional) if set then only paths which are <= the shortest path between the begin and end atoms will be included in the results

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::__cxx11::list<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > FindAllPathsOfLengthN(RDKit::ROMol,unsigned int [,bool=True [,bool=False [,int=-1 [,bool=False]]]])

rdkit.Chem.rdmolops.FindAllSubgraphsOfLengthMToN((Mol)mol, (int)min, (int)max[, (bool)useHs=False[, (int)rootedAtAtom=-1]]) object :
Finds all subgraphs of a particular length in a molecule

See documentation for FindAllSubgraphsOfLengthN for definitions

C++ signature :

boost::python::api::object FindAllSubgraphsOfLengthMToN(RDKit::ROMol,unsigned int,unsigned int [,bool=False [,int=-1]])

rdkit.Chem.rdmolops.FindAllSubgraphsOfLengthN((Mol)mol, (int)length[, (bool)useHs=False[, (int)rootedAtAtom=-1]]) _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.

  • rootedAtAtom: (optional) if nonzero, only subgraphs from the specified atom will be returned.

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::__cxx11::list<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > FindAllSubgraphsOfLengthN(RDKit::ROMol,unsigned int [,bool=False [,int=-1]])

rdkit.Chem.rdmolops.FindAtomEnvironmentOfRadiusN((Mol)mol, (int)radius, (int)rootedAtAtom[, (bool)useHs=False[, (bool)enforceSize=True[, (AtomPairsParameters)atomMap=None]]]) _vecti :
Find bonds of a particular radius around an atom.

Return empty result if there is no bond at the requested radius.

ARGUMENTS:

  • mol: the molecule to use

  • radius: an integer with the target radius for the environment.

  • rootedAtAtom: the atom to consider

  • useHs: (optional) toggles whether or not bonds to Hs that are part of the graph should be included in the results. Defaults to 0.

  • enforceSize (optional) If set to False, all bonds within the requested radius is collected. Defaults to 1.

  • atomMap: (optional) If provided, it will measure the minimum distance of the atom from the rooted atom (start with 0 from the rooted atom). The result is a pair of the atom ID and the distance.

RETURNS: a vector of bond IDs

C++ signature :

std::vector<int, std::allocator<int> > FindAtomEnvironmentOfRadiusN(RDKit::ROMol,unsigned int,unsigned int [,bool=False [,bool=True [,boost::python::api::object=None]]])

rdkit.Chem.rdmolops.FindPotentialStereo((Mol)mol[, (bool)cleanIt=False[, (bool)flagPossible=True]]) _vectN5RDKit9Chirality10StereoInfoE :

find potential stereo elements in a molecule and returns them as StereoInfo objects Note that this function is still somewhat experimental and the API and results may change in a future release.

C++ signature :

std::vector<RDKit::Chirality::StereoInfo, std::allocator<RDKit::Chirality::StereoInfo> > FindPotentialStereo(RDKit::ROMol {lvalue} [,bool=False [,bool=True]])

rdkit.Chem.rdmolops.FindPotentialStereoBonds((Mol)mol[, (bool)cleanIt=False]) None :
Find bonds than can be cis/trans in a molecule and mark them as ‘any’.

This function finds any double bonds that can potentially be part of a cis/trans system. No attempt is made here to mark them cis or trans

ARGUMENTS:

  • mol: the molecule to use

  • cleanIt: (optional) if this option is set to true, any previous marking of _CIPCode

    on the bond is cleared - otherwise it is left untouched

C++ signature :

void FindPotentialStereoBonds(RDKit::ROMol {lvalue} [,bool=False])

rdkit.Chem.rdmolops.FindRingFamilies((Mol)mol) None :

generate Unique Ring Families

C++ signature :

void FindRingFamilies(RDKit::ROMol)

rdkit.Chem.rdmolops.FindUniqueSubgraphsOfLengthN((Mol)mol, (int)length[, (bool)useHs=False[, (bool)useBO=True[, (int)rootedAtAtom=-1]]]) _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.

  • rootedAtAtom: (optional) if nonzero, only subgraphs from the specified atom will be returned.

RETURNS: a tuple of tuples with bond IDs

C++ signature :

std::__cxx11::list<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > FindUniqueSubgraphsOfLengthN(RDKit::ROMol,unsigned int [,bool=False [,bool=True [,int=-1]]])

rdkit.Chem.rdmolops.FragmentOnBRICSBonds((Mol)mol) Mol :

Return a new molecule with all BRICS bonds broken

C++ signature :

RDKit::ROMol* FragmentOnBRICSBonds(RDKit::ROMol)

rdkit.Chem.rdmolops.FragmentOnBonds((Mol)mol, (AtomPairsParameters)bondIndices[, (bool)addDummies=True[, (AtomPairsParameters)dummyLabels=None[, (AtomPairsParameters)bondTypes=None[, (list)cutsPerAtom=[]]]]]) Mol :

Return a new molecule with all specified bonds broken

ARGUMENTS:

  • mol - the molecule to be modified

  • bondIndices - indices of the bonds to be broken

  • addDummies - toggles addition of dummy atoms to indicate where bonds were broken

  • dummyLabels - used to provide the labels to be used for the dummies. the first element in each pair is the label for the dummy that replaces the bond’s beginAtom, the second is for the dummy that replaces the bond’s endAtom. If not provided, the dummies are labeled with atom indices.

  • bondTypes - used to provide the bond type to use between the fragments and the dummy atoms. If not provided, defaults to single.

  • cutsPerAtom - used to return the number of cuts made at each atom.

RETURNS:

a new Mol with the modifications

C++ signature :

RDKit::ROMol* FragmentOnBonds(RDKit::ROMol,boost::python::api::object [,bool=True [,boost::python::api::object=None [,boost::python::api::object=None [,boost::python::list=[]]]]])

rdkit.Chem.rdmolops.FragmentOnSomeBonds((Mol)mol, (AtomPairsParameters)bondIndices[, (int)numToBreak=1[, (bool)addDummies=True[, (AtomPairsParameters)dummyLabels=None[, (AtomPairsParameters)bondTypes=None[, (bool)returnCutsPerAtom=False]]]]]) tuple :

fragment on some bonds

C++ signature :

boost::python::tuple FragmentOnSomeBonds(RDKit::ROMol,boost::python::api::object [,unsigned int=1 [,bool=True [,boost::python::api::object=None [,boost::python::api::object=None [,bool=False]]]]])

rdkit.Chem.rdmolops.Get3DDistanceMatrix((Mol)mol[, (int)confId=-1[, (bool)useAtomWts=False[, (bool)force=False[, (str)prefix='']]]]) object :

Returns the molecule’s 3D distance matrix.

ARGUMENTS:

  • mol: the molecule to use

  • confId: (optional) chooses the conformer Id to use Default value is -1.

  • 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* Get3DDistanceMatrix(RDKit::ROMol {lvalue} [,int=-1 [,bool=False [,bool=False [,char const*=’’]]]])

rdkit.Chem.rdmolops.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*=’’]]]])

rdkit.Chem.rdmolops.GetAllowNontetrahedralChirality() bool :

returns whether or not recognition of non-tetrahedral chirality from 3D structures is enabled

C++ signature :

bool GetAllowNontetrahedralChirality()

rdkit.Chem.rdmolops.GetDistanceMatrix((Mol)mol[, (bool)useBO=False[, (bool)useAtomWts=False[, (bool)force=False[, (str)prefix='']]]]) object :

Returns the molecule’s topological 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*=’’]]]])

rdkit.Chem.rdmolops.GetFormalCharge((Mol)mol) int :

Returns the formal charge for the molecule.

ARGUMENTS:

  • mol: the molecule to use

C++ signature :

int GetFormalCharge(RDKit::ROMol)

rdkit.Chem.rdmolops.GetMolFrags((Mol)mol[, (bool)asMols=False[, (bool)sanitizeFrags=True[, (AtomPairsParameters)frags=None[, (AtomPairsParameters)fragsMolAtomMapping=None]]]]) 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.

  • sanitizeFrags: (optional) if this is provided and true, the fragments molecules will be sanitized before returning them.

  • frags: (optional, defaults to None) if asMols is true and this is provided

    as an empty list, the result will be mol.GetNumAtoms() long on return and will contain the fragment assignment for each Atom

  • fragsMolAtomMapping: (optional, defaults to None) if asMols is true and this is provided as an empty list, the result will be numFrags long on return, and each entry will contain the indices of the Atoms in that fragment: [(0, 1, 2, 3), (4, 5)]

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 [,bool=True [,boost::python::api::object=None [,boost::python::api::object=None]]]])

rdkit.Chem.rdmolops.GetMostSubstitutedCoreMatch((Mol)mol, (Mol)core, (AtomPairsParameters)matches) object :

Postprocesses the results of a mol.GetSubstructMatches(core) call where mol has explicit Hs and core bears terminal dummy atoms (i.e., R groups). It returns the match with the largest number of non-hydrogen matches to the terminal dummy atoms.

ARGUMENTS:

  • mol: the molecule GetSubstructMatches was run on

  • core: the molecule used as a substructure query

  • matches: the result returned by GetSubstructMatches

RETURNS: the tuple where terminal dummy atoms in the core match the largest

number of non-hydrogen atoms in mol

C++ signature :

_object* GetMostSubstitutedCoreMatch(RDKit::ROMol,RDKit::ROMol,boost::python::api::object)

rdkit.Chem.rdmolops.GetSSSR((Mol)mol[, (bool)includeDativeBonds=False]) _vectSt6vectorIiSaIiEE :

Get the smallest set of simple rings for a molecule.

ARGUMENTS:

  • mol: the molecule to use.

  • includeDativeBonds: whether or not dative bonds should be included in the ring finding.

RETURNS: a sequence of sequences containing the rings found as atom ids

The length of this will be equal to NumBonds-NumAtoms+1 for single-fragment molecules.

C++ signature :

std::vector<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > GetSSSR(RDKit::ROMol {lvalue} [,bool=False])

rdkit.Chem.rdmolops.GetShortestPath((Mol)mol, (int)aid1, (int)aid2) tuple :

Find the shortest path between two atoms using the Bellman-Ford algorithm.

ARGUMENTS:

  • mol: the molecule to use

  • idx1: index of the first atom

  • idx2: index of the second atom

C++ signature :

boost::python::tuple GetShortestPath(RDKit::ROMol,int,int)

rdkit.Chem.rdmolops.GetSymmSSSR((Mol)mol[, (bool)includeDativeBonds=False]) _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.

  • includeDativeBonds: whether or not dative bonds should be included in the ring finding.

RETURNS: a sequence of sequences containing the rings found as atom ids

C++ signature :

std::vector<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > GetSymmSSSR(RDKit::ROMol {lvalue} [,bool=False])

rdkit.Chem.rdmolops.GetUseLegacyStereoPerception() bool :

returns whether or not the legacy stereo perception code is being used

C++ signature :

bool GetUseLegacyStereoPerception()

rdkit.Chem.rdmolops.HapticBondsToDative((Mol)mol) Mol :

One way of showing haptic bonds (such as cyclopentadiene to iron in ferrocene) is to use a dummy atom with a dative bond to the iron atom with the bond labelled with the atoms involved in the organic end of the bond. Another way is to have explicit dative bonds from the atoms of the haptic group to the metal atom. This function converts the former representation to the latter.

ARGUMENTS:

  • mol: the molecule to use

RETURNS:

a modified copy of the molecule

C++ signature :

RDKit::ROMol* HapticBondsToDative(RDKit::ROMol)

rdkit.Chem.rdmolops.HasQueryHs((Mol)mol) tuple :

Check to see if the molecule has query Hs, this is normally used on query molecules such as thos returned from MolFromSmarts Example:

(hasQueryHs, hasUnmergableQueryHs) = HasQueryHs(mol)

if hasUnmergableQueryHs, these query hs cannot be removed by calling MergeQueryHs

C++ signature :

boost::python::tuple HasQueryHs(RDKit::ROMol)

rdkit.Chem.rdmolops.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 False.

NOTES:

  • The molecule is modified in place.

  • this does not modify query bonds which have bond type queries (like those which come from SMARTS) or rings containing them.

  • even if clearAromaticFlags is False the BondType for all modified aromatic bonds will be changed from AROMATIC to SINGLE or DOUBLE Kekulization.

C++ signature :

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

rdkit.Chem.rdmolops.KekulizeIfPossible((Mol)mol[, (bool)clearAromaticFlags=False]) None :

Kekulizes the molecule if possible. Otherwise the molecule is not modified

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 if the kekulization succeds. Default value is False.

NOTES:

  • The molecule is modified in place.

C++ signature :

void KekulizeIfPossible(RDKit::ROMol {lvalue} [,bool=False])

rdkit.Chem.rdmolops.LayeredFingerprint((Mol)mol[, (int)layerFlags=4294967295[, (int)minPath=1[, (int)maxPath=7[, (int)fpSize=2048[, (list)atomCounts=[][, (ExplicitBitVect)setOnlyBits=None[, (bool)branchedPaths=True[, (AtomPairsParameters)fromAtoms=0]]]]]]]]) ExplicitBitVect :

Returns a layered fingerprint for a molecule

NOTE: This function is experimental. The API or results may change from

release to release.

Explanation of the algorithm below.

ARGUMENTS:

  • mol: the molecule to use

  • layerFlags: (optional) which layers to include in the fingerprint See below for definitions. Defaults to all.

  • 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.

  • atomCounts: (optional) if provided, this should be a list at least as long as the number of atoms in the molecule. It will be used to provide the count of the number of paths that set bits each atom is involved in. NOTE: the list is not zeroed out here.

  • setOnlyBits: (optional) if provided, only bits that are set in this bit vector will be set in the result. This is essentially the same as doing:

    res &= setOnlyBits

    but also has an impact on the atomCounts (if being used)

  • branchedPaths: (optional) if set both branched and unbranched paths will be used in the fingerprint. Defaults to True.

  • fromAtoms: (optional) a sequence of atom indices. If provided, only paths/subgraphs starting from these atoms will be used. Defaults to empty.

RETURNS: a DataStructs.ExplicitBitVect with _fpSize_ bits

Layer definitions:
  • 0x01: pure topology

  • 0x02: bond order

  • 0x04: atom types

  • 0x08: presence of rings

  • 0x10: ring sizes

  • 0x20: aromaticity

C++ signature :

ExplicitBitVect* LayeredFingerprint(RDKit::ROMol [,unsigned int=4294967295 [,unsigned int=1 [,unsigned int=7 [,unsigned int=2048 [,boost::python::list=[] [,ExplicitBitVect*=None [,bool=True [,boost::python::api::object=0]]]]]]]])

rdkit.Chem.rdmolops.MergeQueryHs((Mol)mol[, (bool)mergeUnmappedOnly=False[, (bool)mergeIsotopes=False]]) Mol :

merges hydrogens into their neighboring atoms as queries

C++ signature :

RDKit::ROMol* MergeQueryHs(RDKit::ROMol [,bool=False [,bool=False]])

rdkit.Chem.rdmolops.MolAddRecursiveQueries((Mol)mol, (dict)queries, (str)propName) None :

Adds named recursive queries to atoms

C++ signature :

void MolAddRecursiveQueries(RDKit::ROMol {lvalue},boost::python::dict,std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> >)

class rdkit.Chem.rdmolops.MolzipLabel

Bases: enum

AtomMapNumber = rdkit.Chem.rdmolops.MolzipLabel.AtomMapNumber
AtomType = rdkit.Chem.rdmolops.MolzipLabel.AtomType
FragmentOnBonds = rdkit.Chem.rdmolops.MolzipLabel.FragmentOnBonds
Isotope = rdkit.Chem.rdmolops.MolzipLabel.Isotope
names = {'AtomMapNumber': rdkit.Chem.rdmolops.MolzipLabel.AtomMapNumber, 'AtomType': rdkit.Chem.rdmolops.MolzipLabel.AtomType, 'FragmentOnBonds': rdkit.Chem.rdmolops.MolzipLabel.FragmentOnBonds, 'Isotope': rdkit.Chem.rdmolops.MolzipLabel.Isotope}
values = {0: rdkit.Chem.rdmolops.MolzipLabel.AtomMapNumber, 1: rdkit.Chem.rdmolops.MolzipLabel.Isotope, 2: rdkit.Chem.rdmolops.MolzipLabel.FragmentOnBonds, 3: rdkit.Chem.rdmolops.MolzipLabel.AtomType}
class rdkit.Chem.rdmolops.MolzipParams((object)self)

Bases: instance

Parameters controllnig how to zip molecules together

OPTIONS:

label : set the MolzipLabel option [default MolzipLabel.AtomMapNumber]

MolzipLabel.AtomMapNumber: atom maps are on dummy atoms, zip together the corresponding

attaced atoms, i.e. zip ‘C[:1]’ ‘N[:1]’ results in ‘CN’

MolzipLabel.Isotope: isotope labels are on dummy atoms, zip together the corresponding

attaced atoms, i.e. zip ‘C[1*]’ ‘N[1*]’ results in ‘CN’

MolzipLabel.FragmentOnBonds: zip together molecules generated by fragment on bonds.

Note the atom indices cannot change or be reorderd from the output of fragmentOnBonds

MolzipLabel.AtomTypes: choose the atom types to act as matching dummy atoms.

i.e. ‘C[V]’ and ‘N[Xe]’ with atoms pairs [(‘V’, ‘Xe’)] results in ‘CN’

C++ signature :

void __init__(_object*)

property enforceValenceRules

If true (default) enforce valences after zipping Setting this to false allows assembling chemically incorrect fragments.

property generateCoordinates

If true will add depiction coordinates to input molecules and zipped molecule (for molzipFragments only)

property label

Set the atom labelling system to zip together

setAtomSymbols((MolzipParams)self, (AtomPairsParameters)symbols) None :

Set the atom symbols used to zip mols together when using AtomType labeling

C++ signature :

void setAtomSymbols(RDKit::MolzipParams {lvalue},boost::python::api::object)

rdkit.Chem.rdmolops.MurckoDecompose((Mol)mol) Mol :

Do a Murcko decomposition and return the scaffold

C++ signature :

RDKit::ROMol* MurckoDecompose(RDKit::ROMol)

rdkit.Chem.rdmolops.ParseMolQueryDefFile((AtomPairsParameters)fileobj[, (bool)standardize=True[, (str)delimiter='\t'[, (str)comment='//'[, (int)nameColumn=0[, (int)smartsColumn=1]]]]]) dict :

reads query definitions from a simply formatted file

C++ signature :

boost::python::dict ParseMolQueryDefFile(boost::python::api::object {lvalue} [,bool=True [,std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> >=’t’ [,std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> >=’//’ [,unsigned int=0 [,unsigned int=1]]]]])

rdkit.Chem.rdmolops.PathToSubmol((Mol)mol, (AtomPairsParameters)path[, (bool)useQuery=False[, (AtomPairsParameters)atomMap=None]]) Mol :
C++ signature :

RDKit::ROMol* PathToSubmol(RDKit::ROMol,boost::python::api::object {lvalue} [,bool=False [,boost::python::api::object=None]])

rdkit.Chem.rdmolops.PatternFingerprint((Mol)mol[, (int)fpSize=2048[, (list)atomCounts=[][, (ExplicitBitVect)setOnlyBits=None[, (bool)tautomerFingerprints=False]]]]) ExplicitBitVect :

A fingerprint using SMARTS patterns

NOTE: This function is experimental. The API or results may change from

release to release.

C++ signature :

ExplicitBitVect* PatternFingerprint(RDKit::ROMol [,unsigned int=2048 [,boost::python::list=[] [,ExplicitBitVect*=None [,bool=False]]]])

PatternFingerprint( (MolBundle)mol [, (int)fpSize=2048 [, (ExplicitBitVect)setOnlyBits=None [, (bool)tautomerFingerprints=False]]]) -> ExplicitBitVect :

A fingerprint using SMARTS patterns

NOTE: This function is experimental. The API or results may change from

release to release.

C++ signature :

ExplicitBitVect* PatternFingerprint(RDKit::MolBundle [,unsigned int=2048 [,ExplicitBitVect*=None [,bool=False]]])

rdkit.Chem.rdmolops.RDKFingerprint((Mol)mol[, (int)minPath=1[, (int)maxPath=7[, (int)fpSize=2048[, (int)nBitsPerHash=2[, (bool)useHs=True[, (float)tgtDensity=0.0[, (int)minSize=128[, (bool)branchedPaths=True[, (bool)useBondOrder=True[, (AtomPairsParameters)atomInvariants=0[, (AtomPairsParameters)fromAtoms=0[, (AtomPairsParameters)atomBits=None[, (AtomPairsParameters)bitInfo=None]]]]]]]]]]]]]) 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.

  • nBitsPerHash: (optional) number of bits to set per path Defaults to 2.

  • useHs: (optional) include paths involving Hs in the fingerprint if the molecule has explicit Hs. Defaults to True.

  • 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.

  • branchedPaths: (optional) if set both branched and unbranched paths will be used in the fingerprint. Defaults to True.

  • useBondOrder: (optional) if set both bond orders will be used in the path hashes Defaults to True.

  • atomInvariants: (optional) a sequence of atom invariants to use in the path hashes Defaults to empty.

  • fromAtoms: (optional) a sequence of atom indices. If provided, only paths/subgraphs starting from these atoms will be used. Defaults to empty.

  • atomBits: (optional) an empty list. If provided, the result will contain a list containing the bits each atom sets. Defaults to empty.

  • bitInfo: (optional) an empty dict. If provided, the result will contain a dict with bits as keys and corresponding bond paths as values. Defaults to empty.

RETURNS: a DataStructs.ExplicitBitVect with _fpSize_ bits

ALGORITHM:

This algorithm functions by find all subgraphs between minPath and maxPath in length. For each subgraph:

  1. A hash is calculated.

  2. The hash is used to seed a random-number generator

  3. _nBitsPerHash_ 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=2 [,bool=True [,double=0.0 [,unsigned int=128 [,bool=True [,bool=True [,boost::python::api::object=0 [,boost::python::api::object=0 [,boost::python::api::object=None [,boost::python::api::object=None]]]]]]]]]]]]])

rdkit.Chem.rdmolops.ReapplyMolBlockWedging((Mol)mol) None :
Set the wedging to that which was read from the original

MolBlock, over-riding anything that was originally there.

ARGUMENTS:

  • molecule: the molecule to update

C++ signature :

void ReapplyMolBlockWedging(RDKit::ROMol {lvalue})

rdkit.Chem.rdmolops.RemoveAllHs((Mol)mol[, (bool)sanitize=True]) Mol :

Returns a copy of the molecule with all Hs removed.

C++ signature :

RDKit::ROMol* RemoveAllHs(RDKit::ROMol [,bool=True])

rdkit.Chem.rdmolops.RemoveHs((Mol)mol[, (bool)implicitOnly=False[, (bool)updateExplicitCount=False[, (bool)sanitize=True]]]) 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).

  • updateExplicitCount: (optional) if this toggle is set, the explicit H count on atoms with Hs will be updated. Default value is 0 (do not update explicit H count).

  • sanitize: (optional) if this toggle is set, the molecule will be sanitized after the Hs are removed. Default value is 1 (do sanitize).

RETURNS: a new molecule with the Hs removed

NOTES:

  • The original molecule is not modified.

  • Hydrogens which aren’t connected to a heavy atom will not be removed. This prevents molecules like [H][H] from having all atoms removed.

  • Labelled hydrogen (e.g. atoms with atomic number=1, but isotope > 1), will not be removed.

  • two coordinate Hs, like the central H in C[H-]C, will not be removed

  • Hs connected to dummy atoms will not be removed

  • Hs that are part of the definition of double bond Stereochemistry will not be removed

  • Hs that are not connected to anything else will not be removed

C++ signature :

RDKit::ROMol* RemoveHs(RDKit::ROMol [,bool=False [,bool=False [,bool=True]]])

RemoveHs( (Mol)mol, (RemoveHsParameters)params [, (bool)sanitize=True]) -> Mol :

Returns a copy of the molecule with Hs removed. Which Hs are removed is controlled by the params argument

C++ signature :

RDKit::ROMol* RemoveHs(RDKit::ROMol,RDKit::MolOps::RemoveHsParameters [,bool=True])

class rdkit.Chem.rdmolops.RemoveHsParameters((object)arg1)

Bases: instance

Parameters controlling which Hs are removed.

C++ signature :

void __init__(_object*)

property removeAndTrackIsotopes

hydrogens with non-default isotopes and store them in the _isotopicHs atom property such that AddHs() can add the same isotope at a later stage

property removeDefiningBondStereo

hydrogens defining bond stereochemistry

property removeDegreeZero

hydrogens that have no bonds

property removeDummyNeighbors

hydrogens with at least one dummy-atom neighbor

property removeHigherDegrees

hydrogens with two (or more) bonds

property removeHydrides

hydrogens with formal charge -1

property removeInSGroups

hydrogens involved in SubstanceGroups

property removeIsotopes

hydrogens with non-default isotopes

property removeMapped

mapped hydrogens

property removeNonimplicit

DEPRECATED

property removeNontetrahedralNeighbors

hydrogens with neighbors that have non-tetrahedral stereochemistry

property removeOnlyHNeighbors

hydrogens with bonds only to other hydrogens

property removeWithQuery

hydrogens with queries defined

property removeWithWedgedBond

hydrogens with wedged bonds to them

property showWarnings

display warning messages for some classes of removed Hs

property updateExplicitCount

DEPRECATED

rdkit.Chem.rdmolops.RemoveNonExplicit3DChirality((Mol)mol) None :
Remove chiral markings that were derived from a 3D mol but were not

explicity marked in the mol block. (wedge bond or CFG indication

ARGUMENTS:

  • molecule: the molecule to update

C++ signature :

void RemoveNonExplicit3DChirality(RDKit::ROMol {lvalue})

rdkit.Chem.rdmolops.RemoveStereochemistry((Mol)mol) None :

Removes all stereochemistry info from the molecule.

C++ signature :

void RemoveStereochemistry(RDKit::ROMol {lvalue})

rdkit.Chem.rdmolops.RenumberAtoms((Mol)mol, (AtomPairsParameters)newOrder) Mol :

Returns a copy of a molecule with renumbered atoms

ARGUMENTS:

  • mol: the molecule to be modified

  • newOrder: the new ordering the atoms (should be numAtoms long) for example: if newOrder is [3,2,0,1], then atom 3 in the original molecule will be atom 0 in the new one

C++ signature :

RDKit::ROMol* RenumberAtoms(RDKit::ROMol,boost::python::api::object {lvalue})

rdkit.Chem.rdmolops.ReplaceCore((Mol)mol, (Mol)core, (AtomPairsParameters)matches[, (bool)replaceDummies=True[, (bool)labelByIndex=False[, (bool)requireDummyMatch=False]]]) Mol :
Removes the core of a molecule and labels the sidechains with dummy atoms based on

The matches indices given in the matching vector matches. Calling:

ReplaceCore(mol,core,mol.GetSubstructMatch(core))

ARGUMENTS:

  • mol: the molecule to be modified

  • coreQuery: the molecule to be used as a substructure query for recognizing the core

  • matches: a matching vector of the type returned by mol.GetSubstructMatch(…)

  • replaceDummies: toggles replacement of atoms that match dummies in the query

  • labelByIndex: toggles labeling the attachment point dummy atoms with the index of the core atom they’re attached to.

  • requireDummyMatch: if the molecule has side chains that attach at points not flagged with a dummy, it will be rejected (None is returned)

RETURNS: a new molecule with the core removed

NOTES:

  • The original molecule is not modified.

EXAMPLES:

>>> from rdkit.Chem import MolToSmiles, MolFromSmiles, ReplaceCore
>>> mol = MolFromSmiles('C1ONNCC1')
>>> core = MolFromSmiles('NN')
>>> MolToSmiles(ReplaceCore(mol, core, mol.GetSubstructMatch(core)))
'[1*]OCCC[2*]'

Since NN is symmetric, we should actually get two matches here if we don’t uniquify the matches.

>>> [MolToSmiles(ReplaceCore(mol, core, match))
...     for match in mol.GetSubstructMatches(core, uniquify=False)]
['[1*]OCCC[2*]', '[1*]CCCO[2*]']
C++ signature :

RDKit::ROMol* ReplaceCore(RDKit::ROMol,RDKit::ROMol,boost::python::api::object [,bool=True [,bool=False [,bool=False]]])

ReplaceCore( (Mol)mol, (Mol)coreQuery [, (bool)replaceDummies=True [, (bool)labelByIndex=False [, (bool)requireDummyMatch=False [, (bool)useChirality=False]]]]) -> 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

  • labelByIndex: toggles labeling the attachment point dummy atoms with the index of the core atom they’re attached to.

  • requireDummyMatch: if the molecule has side chains that attach at points not flagged with a dummy, it will be rejected (None is returned)

  • useChirality: use chirality matching in the coreQuery

RETURNS: a new molecule with the core removed

NOTES:

  • The original molecule is not modified.

EXAMPLES:

>>> from rdkit.Chem import MolToSmiles, MolFromSmiles, MolFromSmarts, ReplaceCore

Basic usage: remove a core as specified by SMILES (or another molecule). To get the atom labels which are stored as an isotope of the matched atom, the output must be written as isomeric smiles. A small confusion is that atom isotopes of 0 aren’t shown in smiles strings.

Here we remove a ring and leave the decoration (r-group) behind.

>>> MolToSmiles(ReplaceCore(MolFromSmiles('CCCC1CCC1'),MolFromSmiles('C1CCC1')),
...             isomericSmiles=True)
'[1*]CCC'

The isotope label by default is matched by the first connection found. In order to indicate which atom the decoration is attached in the core query, use labelByIndex=True. Here the attachment is from the third atom in the smiles string, which is indexed by 3 in the core, like all good computer scientists expect, atoms indices start at 0.

>>> MolToSmiles(ReplaceCore(MolFromSmiles('CCN1CCC1'),MolFromSmiles('C1CCN1'),
...                         labelByIndex=True),
...   isomericSmiles=True)
'[3*]CC'

Non-core matches just return None

>>> ReplaceCore(MolFromSmiles('CCC1CC1'),MolFromSmiles('C1CCC1'))

The bond between atoms are considered part of the core and are removed as well

>>> MolToSmiles(ReplaceCore(MolFromSmiles('C1CC2C1CCC2'),MolFromSmiles('C1CCC1')),
...             isomericSmiles=True)
'[1*]CCC[2*]'
>>> MolToSmiles(ReplaceCore(MolFromSmiles('C1CNCC1'),MolFromSmiles('N')),
...             isomericSmiles=True)
'[1*]CCCC[2*]'

When using dummy atoms, cores should be read in as SMARTS. When read as SMILES dummy atoms only match other dummy atoms. The replaceDummies flag indicates whether matches to the dummy atoms should be considered as part of the core or as part of the decoration (r-group)

>>> MolToSmiles(ReplaceCore(MolFromSmiles('C1CNCC1'),MolFromSmarts('[*]N[*]'),
...                         replaceDummies=True),
...             isomericSmiles=True)
'[1*]CC[2*]'
>>> MolToSmiles(ReplaceCore(MolFromSmiles('C1CNCC1'),MolFromSmarts('[*]N[*]'),
...                         replaceDummies=False),
...             isomericSmiles=True)
'[1*]CCCC[2*]'
>>> MolToSmiles(ReplaceCore(MolFromSmiles('C1CCC1CN'),MolFromSmarts('C1CCC1[*]'),
...                         replaceDummies=False),
...             isomericSmiles=True)
'[1*]CN'
C++ signature :

RDKit::ROMol* ReplaceCore(RDKit::ROMol,RDKit::ROMol [,bool=True [,bool=False [,bool=False [,bool=False]]]])

rdkit.Chem.rdmolops.ReplaceSidechains((Mol)mol, (Mol)coreQuery[, (bool)useChirality=False]) 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

  • useChirality: (optional) match the substructure query using chirality

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 [,bool=False])

rdkit.Chem.rdmolops.ReplaceSubstructs((Mol)mol, (Mol)query, (Mol)replacement[, (bool)replaceAll=False[, (int)replacementConnectionPoint=0[, (bool)useChirality=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)

  • replacementConnectionPoint: (optional) index of the atom in the replacement that the bond should be made to.

  • useChirality: (optional) match the substructure query using chirality

RETURNS: a tuple of new molecules with the substructures replaced removed

NOTES:

  • The original molecule is not modified.

  • A bond is only formed to the remaining atoms, if any, that were bonded to the first atom in the substructure query. (For finer control over substructure replacement, consider using ChemicalReaction.)

EXAMPLES:

The following examples substitute SMILES/SMARTS strings for molecules, you’d have to actually use molecules:

  • ReplaceSubstructs(‘CCOC’,’O[CH3]’,’NC’) -> (‘CCNC’,)

  • ReplaceSubstructs(‘COCCOC’,’O[CH3]’,’NC’) -> (‘COCCNC’,’CNCCOC’)

  • ReplaceSubstructs(‘COCCOC’,’O[CH3]’,’NC’,True) -> (‘CNCCNC’,)

  • ReplaceSubstructs(‘COCCOC’,’O[CH3]’,’CN’,True,1) -> (‘CNCCNC’,)

  • ReplaceSubstructs(‘CCOC’,’[CH3]O’,’NC’) -> (‘CC.CN’,)

C++ signature :

_object* ReplaceSubstructs(RDKit::ROMol,RDKit::ROMol,RDKit::ROMol [,bool=False [,unsigned int=0 [,bool=False]]])

class rdkit.Chem.rdmolops.SanitizeFlags

Bases: enum

SANITIZE_ADJUSTHS = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_ADJUSTHS
SANITIZE_ALL = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_ALL
SANITIZE_CLEANUP = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUP
SANITIZE_CLEANUPATROPISOMERS = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUPATROPISOMERS
SANITIZE_CLEANUPCHIRALITY = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUPCHIRALITY
SANITIZE_CLEANUP_ORGANOMETALLICS = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUP_ORGANOMETALLICS
SANITIZE_FINDRADICALS = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_FINDRADICALS
SANITIZE_KEKULIZE = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_KEKULIZE
SANITIZE_NONE = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_NONE
SANITIZE_PROPERTIES = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_PROPERTIES
SANITIZE_SETAROMATICITY = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SETAROMATICITY
SANITIZE_SETCONJUGATION = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SETCONJUGATION
SANITIZE_SETHYBRIDIZATION = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SETHYBRIDIZATION
SANITIZE_SYMMRINGS = rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SYMMRINGS
names = {'SANITIZE_ADJUSTHS': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_ADJUSTHS, 'SANITIZE_ALL': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_ALL, 'SANITIZE_CLEANUP': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUP, 'SANITIZE_CLEANUPATROPISOMERS': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUPATROPISOMERS, 'SANITIZE_CLEANUPCHIRALITY': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUPCHIRALITY, 'SANITIZE_CLEANUP_ORGANOMETALLICS': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUP_ORGANOMETALLICS, 'SANITIZE_FINDRADICALS': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_FINDRADICALS, 'SANITIZE_KEKULIZE': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_KEKULIZE, 'SANITIZE_NONE': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_NONE, 'SANITIZE_PROPERTIES': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_PROPERTIES, 'SANITIZE_SETAROMATICITY': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SETAROMATICITY, 'SANITIZE_SETCONJUGATION': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SETCONJUGATION, 'SANITIZE_SETHYBRIDIZATION': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SETHYBRIDIZATION, 'SANITIZE_SYMMRINGS': rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SYMMRINGS}
values = {0: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_NONE, 1: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUP, 2: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_PROPERTIES, 4: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SYMMRINGS, 8: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_KEKULIZE, 16: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_FINDRADICALS, 32: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SETAROMATICITY, 64: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SETCONJUGATION, 128: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_SETHYBRIDIZATION, 256: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUPCHIRALITY, 512: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_ADJUSTHS, 1024: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUP_ORGANOMETALLICS, 2048: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_CLEANUPATROPISOMERS, 268435455: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_ALL}
rdkit.Chem.rdmolops.SanitizeMol((Mol)mol[, (int)sanitizeOps=rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_ALL[, (bool)catchErrors=False]]) SanitizeFlags :

Kekulize, check valencies, set aromaticity, conjugation and hybridization

  • The molecule is modified in place.

  • If sanitization fails, an exception will be thrown unless catchErrors is set

ARGUMENTS:

  • mol: the molecule to be modified

  • sanitizeOps: (optional) sanitization operations to be carried out these should be constructed by or’ing together the operations in rdkit.Chem.SanitizeFlags

  • catchErrors: (optional) if provided, instead of raising an exception when sanitization fails (the default behavior), the first operation that failed (as defined in rdkit.Chem.SanitizeFlags) is returned. Zero is returned on success.

C++ signature :

RDKit::MolOps::SanitizeFlags SanitizeMol(RDKit::ROMol {lvalue} [,unsigned long=rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_ALL [,bool=False]])

rdkit.Chem.rdmolops.SetAllowNontetrahedralChirality((bool)val) None :

toggles recognition of non-tetrahedral chirality from 3D structures

C++ signature :

void SetAllowNontetrahedralChirality(bool)

rdkit.Chem.rdmolops.SetAromaticity((Mol)mol[, (AromaticityModel)model=rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_DEFAULT]) None :

does aromaticity perception

ARGUMENTS:

  • mol: the molecule to use

  • model: the model to use

NOTES:

  • The molecule is modified in place.

C++ signature :

void SetAromaticity(RDKit::ROMol {lvalue} [,RDKit::MolOps::AromaticityModel=rdkit.Chem.rdmolops.AromaticityModel.AROMATICITY_DEFAULT])

rdkit.Chem.rdmolops.SetBondStereoFromDirections((Mol)mol) None :

Uses the directions of neighboring bonds to set cis/trans stereo on double bonds.

ARGUMENTS:

  • mol: the molecule to be modified

C++ signature :

void SetBondStereoFromDirections(RDKit::ROMol {lvalue})

rdkit.Chem.rdmolops.SetConjugation((Mol)mol) None :

finds conjugated bonds

ARGUMENTS:

  • mol: the molecule to use

NOTES:

  • The molecule is modified in place.

C++ signature :

void SetConjugation(RDKit::ROMol {lvalue})

rdkit.Chem.rdmolops.SetDoubleBondNeighborDirections((Mol)mol[, (AtomPairsParameters)conf=None]) None :

Uses the stereo info on double bonds to set the directions of neighboring single bonds

ARGUMENTS:

  • mol: the molecule to be modified

C++ signature :

void SetDoubleBondNeighborDirections(RDKit::ROMol {lvalue} [,boost::python::api::object=None])

rdkit.Chem.rdmolops.SetGenericQueriesFromProperties((Mol)mol[, (bool)useAtomLabels=True[, (bool)useSGroups=True]]) None :

documentation

C++ signature :

void SetGenericQueriesFromProperties(RDKit::ROMol {lvalue} [,bool=True [,bool=True]])

rdkit.Chem.rdmolops.SetHybridization((Mol)mol) None :

Assigns hybridization states to atoms

ARGUMENTS:

  • mol: the molecule to use

NOTES:

  • The molecule is modified in place.

C++ signature :

void SetHybridization(RDKit::ROMol {lvalue})

rdkit.Chem.rdmolops.SetTerminalAtomCoords((Mol)mol, (int)idx, (int)otherIdx) None :

Sets Cartesian coordinates for a terminal atom.

Useful for growing an atom off a molecule with sensible coordinates based on the geometry of the neighbor.

NOTE: this sets the appropriate coordinates in all of the molecule’s conformers ARGUMENTS:

  • mol: the molecule the atoms belong to.

  • idx: index of the terminal atom whose coordinates are set.

  • mol: index of the bonded neighbor atom.

RETURNS: Nothing

C++ signature :

void SetTerminalAtomCoords(RDKit::ROMol {lvalue},unsigned int,unsigned int)

rdkit.Chem.rdmolops.SetUseLegacyStereoPerception((bool)val) None :

sets usage of the legacy stereo perception code

C++ signature :

void SetUseLegacyStereoPerception(bool)

rdkit.Chem.rdmolops.SortMatchesByDegreeOfCoreSubstitution((Mol)mol, (Mol)core, (AtomPairsParameters)matches) object :

Postprocesses the results of a mol.GetSubstructMatches(core) call where mol has explicit Hs and core bears terminal dummy atoms (i.e., R groups). It returns a copy of matches sorted by decreasing number of non-hydrogen matches to the terminal dummy atoms.

ARGUMENTS:

  • mol: the molecule GetSubstructMatches was run on

  • core: the molecule used as a substructure query

  • matches: the result returned by GetSubstructMatches

RETURNS: a copy of matches sorted by decreasing number of non-hydrogen matches

to the terminal dummy atoms

C++ signature :

_object* SortMatchesByDegreeOfCoreSubstitution(RDKit::ROMol,RDKit::ROMol,boost::python::api::object)

rdkit.Chem.rdmolops.SplitMolByPDBChainId((Mol)mol[, (AtomPairsParameters)whiteList=None[, (bool)negateList=False]]) dict :

Splits a molecule into pieces based on PDB chain information.

ARGUMENTS:

  • mol: the molecule to use

  • whiteList: only residues in this list will be returned

  • negateList: if set, negates the white list inclusion logic

RETURNS: a dictionary keyed by chain id with molecules as the values

C++ signature :

boost::python::dict SplitMolByPDBChainId(RDKit::ROMol [,boost::python::api::object=None [,bool=False]])

rdkit.Chem.rdmolops.SplitMolByPDBResidues((Mol)mol[, (AtomPairsParameters)whiteList=None[, (bool)negateList=False]]) dict :

Splits a molecule into pieces based on PDB residue information.

ARGUMENTS:

  • mol: the molecule to use

  • whiteList: only residues in this list will be returned

  • negateList: if set, negates the white list inclusion logic

RETURNS: a dictionary keyed by residue name with molecules as the values

C++ signature :

boost::python::dict SplitMolByPDBResidues(RDKit::ROMol [,boost::python::api::object=None [,bool=False]])

class rdkit.Chem.rdmolops.StereoBondThresholds

Bases: instance

Constants used to set the thresholds for which single bonds can be made wavy.

Raises an exception This class cannot be instantiated from Python

CHIRAL_ATOM = 100000
DBL_BOND_NO_STEREO = 1000
DBL_BOND_SPECIFIED_STEREO = 10000
DIRECTION_SET = 1000000
rdkit.Chem.rdmolops.TranslateChiralFlagToStereoGroups((Mol)mol[, (StereoGroupType)zeroFlagGroupType=rdkit.Chem.rdchem.StereoGroupType.STEREO_AND]) None :

Generate enhanced stereo groups based on the status of the chiral flag property.

Arguments:
  • mol: molecule to be modified

  • zeroFlagGroupType: how to handle non-grouped stereo centers when the

    chiral flag is set to zero

If the chiral flag is set to a value of 1 then all specified tetrahedral chiral centers which are not already in StereoGroups will be added to an ABS StereoGroup.

If the chiral flag is set to a value of 0 then all specified tetrahedral chiral centers will be added to a StereoGroup of the type zeroFlagGroupType

If there is no chiral flag set (i.e. the property is not present), the molecule will not be modified.

C++ signature :

void TranslateChiralFlagToStereoGroups(RDKit::ROMol {lvalue} [,RDKit::StereoGroupType=rdkit.Chem.rdchem.StereoGroupType.STEREO_AND])

rdkit.Chem.rdmolops.UnfoldedRDKFingerprintCountBased((Mol)mol[, (int)minPath=1[, (int)maxPath=7[, (bool)useHs=True[, (bool)branchedPaths=True[, (bool)useBondOrder=True[, (AtomPairsParameters)atomInvariants=0[, (AtomPairsParameters)fromAtoms=0[, (AtomPairsParameters)atomBits=None[, (AtomPairsParameters)bitInfo=None]]]]]]]]]) ULongSparseIntVect :

Returns an unfolded count-based version of the RDKit fingerprint for a molecule

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.

  • useHs: (optional) include paths involving Hs in the fingerprint if the molecule has explicit Hs. Defaults to True.

  • branchedPaths: (optional) if set both branched and unbranched paths will be used in the fingerprint. Defaults to True.

  • useBondOrder: (optional) if set both bond orders will be used in the path hashes Defaults to True.

  • atomInvariants: (optional) a sequence of atom invariants to use in the path hashes Defaults to empty.

  • fromAtoms: (optional) a sequence of atom indices. If provided, only paths/subgraphs starting from these atoms will be used. Defaults to empty.

  • atomBits: (optional) an empty list. If provided, the result will contain a list containing the bits each atom sets. Defaults to empty.

  • bitInfo: (optional) an empty dict. If provided, the result will contain a dict with bits as keys and corresponding bond paths as values. Defaults to empty.

C++ signature :

RDKit::SparseIntVect<unsigned long>* UnfoldedRDKFingerprintCountBased(RDKit::ROMol [,unsigned int=1 [,unsigned int=7 [,bool=True [,bool=True [,bool=True [,boost::python::api::object=0 [,boost::python::api::object=0 [,boost::python::api::object=None [,boost::python::api::object=None]]]]]]]]])

rdkit.Chem.rdmolops.WedgeBond((Bond)bond, (int)fromAtomIdx, (Conformer)conf) None :
Set the wedging on an individual bond from a molecule.

The wedging scheme used is that from Mol files.

ARGUMENTS:
  • bond: the bond to update

  • atom ID: the atom from which to do the wedging

  • conformer: the conformer to use to determine wedge direction

C++ signature :

void WedgeBond(RDKit::Bond*,unsigned int,RDKit::Conformer const*)

rdkit.Chem.rdmolops.WedgeMolBonds((Mol)mol, (Conformer)conformer[, (BondWedgingParameters)params=None]) 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

  • conformer: the conformer to use to determine wedge direction

C++ signature :

void WedgeMolBonds(RDKit::ROMol {lvalue},RDKit::Conformer const* [,RDKit::Chirality::BondWedgingParameters const*=None])

rdkit.Chem.rdmolops.molzip((Mol)a, (Mol)b[, (MolzipParams)params=<rdkit.Chem.rdmolops.MolzipParams object at 0x7461926b5170>]) Mol :

zip together two molecules using the given matching parameters

C++ signature :

RDKit::ROMol* molzip(RDKit::ROMol,RDKit::ROMol [,RDKit::MolzipParams=<rdkit.Chem.rdmolops.MolzipParams object at 0x7461926b5170>])

molzip( (Mol)a [, (MolzipParams)params=<rdkit.Chem.rdmolops.MolzipParams object at 0x7461926b5220>]) -> Mol :

zip together two molecules using the given matching parameters

C++ signature :

RDKit::ROMol* molzip(RDKit::ROMol [,RDKit::MolzipParams=<rdkit.Chem.rdmolops.MolzipParams object at 0x7461926b5220>])

molzip( (dict)row [, (MolzipParams)params=<rdkit.Chem.rdmolops.MolzipParams object at 0x7461926b5380>]) -> Mol :

zip an RGroupRow together to recreate the original molecule. This correctly handles broken cycles that can occur in decompositions.

example:

>>> from rdkit import Chem
>>> from rdkit.Chem import rdRGroupDecomposition as rgd
>>> core = Chem.MolFromSmiles('CO')
>>> mols = [Chem.MolFromSmiles('C1NNO1')]
>>> rgroups, unmatched = rgd.RGroupDecompose(core, mols)
>>> for rgroup in rgroups:
...     mol = rgd.molzip(rgroup)
C++ signature :

RDKit::ROMol* molzip(boost::python::dict [,RDKit::MolzipParams=<rdkit.Chem.rdmolops.MolzipParams object at 0x7461926b5380>])

rdkit.Chem.rdmolops.molzipFragments((AtomPairsParameters)mols[, (MolzipParams)params=<rdkit.Chem.rdmolops.MolzipParams object at 0x7461926b52d0>]) Mol :

zip together multiple molecules from an R group decomposition using the given matching parameters. The first molecule in the list must be the core

C++ signature :

RDKit::ROMol* molzipFragments(boost::python::api::object {lvalue} [,RDKit::MolzipParams=<rdkit.Chem.rdmolops.MolzipParams object at 0x7461926b52d0>])