moleculekit.rdkittools module#

moleculekit.rdkittools.extend_residue_from_smiles(mol: Molecule, sel: str, extension_smiles: str | None = None, target_atom_sel: str | None = None, new_smiles: str | None = None, sanitizeSmiles: bool = True, _logger: bool = True)#

Extends a residue by attaching new atoms, preserving original 3D coordinates and relaxing new atoms via constrained embedding and MMFF minimization.

Two modes are supported (mutually exclusive):

  1. Extension SMILES mode: provide extension_smiles (containing a dummy atom * at the attachment point) together with target_atom_sel (the atom on the residue to attach to).

  2. New SMILES mode: provide new_smiles with the complete SMILES of the modified molecule. The function uses Maximum Common Substructure (MCS) matching to identify unchanged atoms and generates 3D coordinates for the new ones.

Parameters:
  • mol (Molecule) – The molecule containing the residue to extend

  • sel (str) – Atom selection for the residue to extend

  • extension_smiles (str, optional) – SMILES of the extension moiety with a dummy atom *

  • target_atom_sel (str, optional) – Atom selection of the attachment point (required with extension_smiles)

  • new_smiles (str, optional) – Complete SMILES of the modified molecule (alternative to extension_smiles)

  • sanitizeSmiles (bool) – If True the SMILES string will be sanitized

  • _logger (bool) – If True the logger will be used to print information

Examples

>>> mol = Molecule('3ptb')
>>> mol.templateResidueFromSmiles("resname BEN", "[NH2+]=C(N)c1ccccc1", addHs=True)
>>> mol.extendResidueFromSmiles("resname BEN", extension_smiles="*C(C)(C)C", target_atom_sel="resname BEN and name H6")

Or equivalently using the full modified SMILES:

>>> mol = Molecule('3ptb')
>>> mol.templateResidueFromSmiles("resname BEN", "[NH2+]=C(N)c1ccccc1", addHs=True)
>>> mol.extendResidueFromSmiles("resname BEN", new_smiles="[NH2+]=C(N)c1cc(C(C)(C)C)ccc1")
moleculekit.rdkittools.molecule_to_rdkitmol(mol: Molecule, sanitize=False, kekulize=False, assignStereo=True, guessBonds=False, _logger=True)#

Converts the Molecule to an RDKit molecule

Parameters:
  • mol (Molecule) – The molecule to convert to an RDKit molecule

  • sanitize (bool) – If True the molecule will be sanitized

  • kekulize (bool) – If True the molecule will be kekulized

  • assignStereo (bool) – If True the molecule will have stereochemistry assigned from its 3D coordinates

  • guessBonds (bool) – If True the molecule will have bonds guessed

  • _logger (bool) – If True the logger will be used to print information

moleculekit.rdkittools.rdkitmol_to_molecule(rmol)#

Converts an RDKit molecule to a Molecule object

Parameters:

rmol (rdkit.Chem.rdchem.Mol) – The RDKit molecule to convert

Examples

>>> from rdkit import Chem
>>> from rdkit.Chem import AllChem
>>> rmol = Chem.MolFromSmiles("C1CCCCC1")
>>> rmol = Chem.AddHs(rmol)
>>> AllChem.EmbedMolecule(rmol)
>>> AllChem.MMFFOptimizeMolecule(rmol)
>>> mol = Molecule.fromRDKitMol(rmol)
moleculekit.rdkittools.template_residue_from_smiles(mol: Molecule, sel: str, smiles: str, sanitizeSmiles: bool = True, addHs: bool = False, onlyOnAtoms: str | None = None, guessBonds: bool = False, _logger: bool = True)#

Template a residue from a SMILES string

This function will assign bonds, bond orders and formal charges to a residue according to a corresponding SMILES string. In addition it can also protonate the residue.

Parameters:
  • mol (Molecule) – The molecule to template the residue in

  • sel (str) – The atom selection of the residue which we want to template

  • smiles (str) – The SMILES string of the template residue

  • sanitizeSmiles (bool) – If True the SMILES string will be sanitized

  • addHs (bool) – If True the residue will be protonated

  • onlyOnAtoms (str) – If not None, only the atoms in this atom selection will be protonated

  • guessBonds (bool) – Set to True to guess bonds for the residue we are templating

  • _logger (bool) – If True the logger will be used to print information

Examples

>>> mol = Molecule("3ptb")
>>> mol.templateResidueFromSmiles("resname BEN", "[NH2+]=C(N)c1ccccc1", addHs=True)
>>> mol.templateResidueFromSmiles("resname GLY and resid 18", "C(C(=O))N", addHs=True)