moleculekit.rdkittools module#
- moleculekit.rdkittools.extend_residue_from_smiles(mol, sel, extension_smiles=None, target_atom_sel=None, new_smiles=None, sanitizeSmiles=True, minimize=False, _logger=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):
Extension SMILES mode: provide
extension_smiles(containing a dummy atom*at the attachment point) together withtarget_atom_sel(the atom on the residue to attach to).New SMILES mode: provide
new_smileswith 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 extendsel (
str) – Atom selection for the residue to extendextension_smiles (
str|None) – SMILES of the extension moiety with a dummy atom*target_atom_sel (
str|None) – Atom selection of the attachment point (required with extension_smiles)new_smiles (
str|None) – Complete SMILES of the modified molecule (alternative to extension_smiles)sanitizeSmiles (
bool) – If True the SMILES string will be sanitizedminimize (
bool) – If True and OpenMM is available, run a soft-potential energy minimization of the residue against its surroundings after insertion.
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, 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 moleculesanitize (
bool) – If True the molecule will be sanitizedkekulize (
bool) – If True the molecule will be kekulizedassignStereo (
bool) – If True the molecule will have stereochemistry assigned from its 3D coordinatesguessBonds (
bool) – If True the molecule will have bonds guessed
- 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_molecule(mol, sel, refmol, addHs=False, onlyOnAtoms=None, guessBonds=False, _logger=True)#
Assign bonds, bond orders, formal charges and (optionally) hydrogens to a residue from a reference
Moleculetemplate, matched by atom name.Mirrors
template_residue_from_smiles()but the template source is a reference Molecule (e.g.Molecule("LIG.cif")) that already carriesbonds,bondtypeandformalcharge. Heavy atoms of the selected residue are mapped onto the reference by NAME (not MCS); bond orders and formal charges are transferred through that mapping. The reference is used only as a template and is never appended tomol.molis mutated in place.The residue and the reference must have the same set of heavy-atom names.
- Parameters:
mol (
Molecule) – The molecule containing the residue(s) to template. Mutated in place.sel (
strornumpy.ndarray) – VMD-style atom selection or boolean mask. May span multiple residues with the same chemistry; each residue is templated in sequence with the same reference.refmol (
Molecule) – Reference Molecule used as the template. Must carrybonds,bondtypeandformalcharge, and have unique heavy-atom names matching the residue.addHs (
bool) – If True, add hydrogens after bond orders are transferred. Boundary atoms (those involved in cross-residue covalent bonds) have their explicit H count reduced by the order of the external bond so they are not over-protonated.onlyOnAtoms (
strornumpy.ndarray) – VMD-style selection or boolean mask within the residue restricting which heavy atoms get hydrogens added. Only used whenaddHs=True.guessBonds (
bool) – If True, run distance-based bond guessing on the residue before templating. Use whenmol.bondsis empty.
- Raises:
RuntimeError – If the selection is empty/gapped/multi-residue, has no bonds, the reference has duplicate heavy-atom names, or the residue and reference heavy-atom names do not match.
- moleculekit.rdkittools.template_residue_from_smiles(mol, sel, smiles, sanitizeSmiles=True, addHs=False, onlyOnAtoms=None, guessBonds=False, _logger=True)#
Assign bonds, bond orders, formal charges and (optionally) hydrogens to a residue from a SMILES template.
See
moleculekit.molecule.Molecule.templateResidueFromSmiles()for the full description; this is the underlying implementation. Themolargument is mutated in place.- Parameters:
mol (
Molecule) – The molecule containing the residue(s) to template. Mutated in place.sel (
str) – VMD-style atom selection or boolean mask. May span multiple residues with the same chemistry; each residue is templated in sequence with the same SMILES.smiles (
str) – SMILES string of the template residue. RCSB-style (fully protonated, explicit charges, full heavy-atom set) works best.sanitizeSmiles (
bool) – If True, sanitize the SMILES with RDKit before matching.addHs (
bool) – If True, add hydrogens after bond orders are transferred. Boundary atoms (those involved in cross-residue covalent bonds) have their explicit H count reduced by the order of the external bond so they are not over-protonated.onlyOnAtoms (
str|None) – VMD-style selection within the residue restricting which heavy atoms get hydrogens added. Only used whenaddHs=True.guessBonds (
bool) – If True, run distance-based bond guessing on the residue before templating. Use whenmol.bondsis empty.
- Raises:
RuntimeError – If the selection is empty, has gaps in atom indexes, contains multiple residues with conflicting metadata, has no bonds, or the SMILES contains heavy atoms that cannot be matched (even after stripping recognized terminal atoms).
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)