moleculekit.smallmol.smallmol module

class moleculekit.smallmol.smallmol.SmallMol(mol, ignore_errors=False, force_reading=False, fixHs=True, removeHs=False)

Bases: object

Class to manipulate small molecule structures

Parameters:
  • mol (rdkit.Chem.rdchem.Mol or filename or smile or moleculekit.smallmol.smallmol.SmallMol) – (i) Rdkit molecule or (ii) Location of molecule file (“.pdb”/”.mol2”) or (iii) a smile string or iv) another SmallMol object or v) moleculekit.molecule.Molecule object
  • ignore_errors (bool) – If True, errors will not be raised.
  • force_reading (bool) – If True, and the mol provided is not accepted, the molecule will be initially converted into sdf
  • fixHs (bool) – If True, the missing hydrogens are assigned, the others are correctly assinged into the graph of the molecule
  • removeHs (bool) – If True, remove the hydrogens

Examples

>>> import os
>>> from moleculekit.home import home
>>> SmallMol('CCO')  
>>> SmallMol(os.path.join(home(dataDir='test-smallmol'), 'ligand.pdb'), fixHs=False, removeHs=True )  
>>> sm = SmallMol(os.path.join(home(dataDir='test-smallmol'), 'benzamidine.mol2'))
>>> print(sm)                                     
SmallMol with 18 atoms and 1 conformers
Atom field - bondtype
Atom field - charge
...

Methods

Attributes

align(refmol)
copy()

Create a copy of the molecule object

Returns:newsmallmol – A copy of the object
Return type:SmallMol
depict(sketch=True, filename=None, ipython=False, optimize=False, optimizemode='std', removeHs=True, atomlabels=None, highlightAtoms=None, resolution=(400, 200))

Depicts the molecules. It is possible to save it into an svg file and also generates a jupiter-notebook rendering

Parameters:
  • sketch (bool) – Set to True for 2D depiction
  • filename (str) – Set the filename for the svg file
  • ipython (bool) – Set to True to return the jupiter-notebook rendering
  • optimize (bool) – Set to True to optimize the conformation. Works only with 3D.
  • optimizemode (['std', 'mmff']) – Set the optimization mode for 3D conformation
  • removeHs (bool) – Set to True to hide hydrogens in the depiction
  • atomlabels (str) – Accept any combinations of the following pararemters as unique string ‘%a%i%c%*’ a:atom name, i:atom index, c:atom formal charge (+/-), :chiral ( if atom is chiral)
  • highlightAtoms (list) – List of atom to highlight. It can be also a list of atom list, in this case different colors will be used
  • resolution (tuple of integers) – Resolution in pixels: (X, Y)
Returns:

ipython_svg

Return type:

SVG object if ipython is set to True

Example

>>> sm.depict(ipython=True, optimize=True, optimizemode='std')  
>>> sm.depict(ipython=True, sketch=True)  
>>> sm.depict(ipython=True, sketch=True)  
>>> sm.depict(ipython=True, sketch=True, atomlabels="%a%i%c")  
>>> ids = np.intersect1d(sm.get('idx', 'hybridization SP2'), sm.get('idx', 'element C'))  
>>> sm.depict(ipython=True, sketch=True,highlightAtoms=ids.tolist(), removeHs=False)  
dropFrames(frames='all')
filter(sel)
foundBondBetween(sel1, sel2, bondtype=None)

Returns True if at least a bond is found between the two selections. It is possible to check for specific bond type. A tuple is returned in the form (bool, [ [(idx1,idx2), rdkit.Chem.rdchem.BondType]] ])

Parameters:
  • sel1 (str) – The selection for the first set of atoms
  • sel2 (str) – The selection for the second set of atoms
  • bondtype (str or int) – The bondtype as index or string Default: None
Returns:

  • isbond (bool) – True if a bond was found
  • details (list) – A list of lists with the index of atoms in the bond and its type

frame
generateConformers(num_confs=400, optimizemode='mmff', align=True, append=True)

Generates ligand conformers

Parameters:
  • num_confs (int) – Number of conformers to generate.
  • optimizemode (str) – The optimizemode to use. Can be ‘uff’, ‘mmff’ Default: ‘mmff’
  • align (bool) – If True, the conformer are aligned to the first one Default: True
  • append (bool) – If False, the current conformers are deleted Default: True
get(returnField, sel='all', convertType=True, invert=False)

Returns the property for the atom specified with the selection. The selection is another atom property

Parameters:
  • returnField (str) – The field of the atom to return
  • sel (str) – The selection string. atom field name followed by spaced values for that field
  • convertType (bool) – If True, and where possible the returnField is converted in rdkit object Default: True
  • invert (bool) – If True, the selection is inverted Default: False
Returns:

values – The array of values for the property

Return type:

np.array

Example

>>> sm.get('element', 'idx 0 1 7')  
array(['C', 'C', 'H'],
      dtype='<U1')
>>> sm.get('hybridization', 'element N')  
array([rdkit.Chem.rdchem.HybridizationType.SP2,
       rdkit.Chem.rdchem.HybridizationType.SP2], dtype=object)
>>> sm.get('hybridization', 'element N', convertType=False)
array([3, 3])
>>> sm.get('element', 'hybridization sp2')  
array(['C', 'C', 'C', 'C', 'C', 'C', 'C', 'N', 'N'],
      dtype='<U1')
>>> sm.get('element', 'hybridization S')  
array(['H', 'H', 'H', 'H', 'H', 'H', 'H', 'H', 'H'],
      dtype='<U1')
>>> sm.get('element', 'hybridization 1')  
array(['H', 'H', 'H', 'H', 'H', 'H', 'H', 'H', 'H'],
      dtype='<U1')
>>> sm.get('atomobject', 'element N')  
array([<rdkit.Chem.rdchem.Atom object at 0x7faf616dd120>,
       <rdkit.Chem.rdchem.Atom object at 0x7faf616dd170>], dtype=object)
getAtoms()

Retuns an array with the rdkit.Chem.rdchem.Atom present in the molecule

getCenter()

Returns geometrical center of molecule conformation

isChiral(returnDetails=False)

Returns True if the molecule has at least one chiral atom. If returnDetails is set as True, a list of tuples with the atom idx and chiral type is returned.

Parameters:returnDetails (bool) – If True, returns the chiral atoms and their chiral types Default: False
Returns:
  • ischiral (bool) – True if the atom has at least a chiral atom
  • details (list) – A list of tuple with the chiral atoms and their types

Example

>>> chiralmol.isChiral()  
True
>>> chiralmol.isChiral(returnDetails=True)  
(True, [('C2', 'R')])
ligname
numAtoms
numFrames
toMolecule(formalcharges=False, ids=None)

Return the moleculekit.molecule.Molecule

Parameters:
  • formalcharges (bool) – If True,the formal charges are used instead of partial ones
  • ids (list) – The list of conformer ids to store in the moleculekit Molecule object- If None, all are returned Default: None
Returns:

mol – The moleculekit Molecule object

Return type:

moleculekit.molecule.Molecule

toSMARTS(explicitHs=False)

Returns the smarts string of the molecule

Parameters:explicitHs (bool) – Set as True for keep the hydrogens
Returns:smart – The smarts string
Return type:str
toSMILES(explicitHs=False, kekulizeSmile=True)

Returns the smiles string of the molecule

Parameters:
  • explicitHs (bool) – Set as True to keep the hydrogens
  • kekulizeSmile (bool) – Set as True to return the kekule smile format
Returns:

smi – The smiles string

Return type:

str

view(*args, **kwargs)
write(fname, frames=None, merge=True)