Hi Sarah,
The problem is that you are using an iterator while at the same time changing
the underlying data structure. `rings` are precomputed atom indices for `frame`
and also `editMol`. As you are removing atoms from `editMol`, the indices in
`rings` are no longer valid to be used with `editMol`. Continuing using it will
get you in trouble. To illustrate the problem, consider the following
hypothetical molecule with 4 atoms:
>>> from rdkit import Chem
>>> mol = Chem.MolFromSmiles('NCOF')
>>> editable = Chem.EditableMol(mol)
Let's say we'd like to remove the O and the F, at index 2 and 3. After removing
atom #2:
>>> Chem.MolToSmiles(editable.GetMol())
'F.CN'
the molecule has only 3 atoms, and now if we continue removing using index 3,
we hit the Range Error because the editable molecule now has only 3 atoms.
>>> editable.RemoveAtom(3)
[10:27:08]
****
Range Error
idx
Violation occurred on line 143 in file
/Users/caoyi3/Development/RDKit-repo/Code/GraphMol/ROMol.cpp
Failed Expression: 0 <= 3 <= 2
****
Traceback (most recent call last):
File "<stdin>", line 1, in <module>
RuntimeError: Range Error
I suggest the following alternative:
1. first make a copy so you don't change the original molecule
>>> from copy import deepcopy
>>> mol2 = deepcopy(mol)
Now set the atoms you'd like to remove to dummy atoms with atomic number 0
>>> mol2.GetAtomWithIdx(2).SetAtomicNum(0)
>>> mol2.GetAtomWithIdx(3).SetAtomicNum(0)
Now remove dummy atoms using a query
>>> mol3 = Chem.DeleteSubstructs(mol2, Chem.MolFromSmarts('[#0]'))
You get what you are looking for
>>> Chem.MolToSmiles(mol3)
'CN'
Regards,
Eddie
On Sep 28, 2011, at 9:16 AM, Sarah Langdon wrote:
> Hi there,
>
> I've tried to write a function in Python to generate the Murcko
> Framework of a molecule, then remove a ring from the framework. I want
> to remove a ring based on the atom ID of the atoms of the ring, rather
> than as a substructure so that in the case of a molecule containing
> more than one of the same ring, only one ring remains. Therefore I
> have used RemoveAtoms to remove each ring atom one by one. My code is
> as follows.
>
> def removeRings(mol):
>
> frame = MurckoScaffold.GetScaffoldForMol(mol)
>
> getRings = frame.GetRingInfo()
> rings = getRings.AtomRings()
>
> splitFrames = [x for x in range(len(rings))]
>
> for x in range(len(rings)):
> editMol = Chem.EditableMol(frame)
>
> for atom in rings[x]:
> editMol.RemoveAtom(atom)
>
> splitFrames[x] = editMol.GetMol()
>
> for x in splitFrames:
> print Chem.MolToSmiles(x)
>
> However when I use this function I get the following error:
>
> Range Error
> idx
> Violation occurred on line 143 in file /usr/local/bin/src/
> RDKit_2011_03_2/Code/GraphMol/ROMol.cpp
> Failed Expression: 0 <= 22 <= 20
>
> I assume that this means that the atom I am trying to remove is not in
> the range of the atoms in the framework. I have checked the atom IDs
> for all atoms in the framework and all atoms in the rings and they are
> within the same range, so I do not understand why I am getting this
> error. Can anyone help please?
>
> Many thanks,
>
> Sarah Langdon
> PhD student
> Cancer Research UK Cancer Therapeutics Unit
> Institute of Cancer Research
> Haddow Laboratories
> 15 Cotswold Road
> Sutton
> Surrey SM2 5NG
>
> Tel: 0208 722 4139
> Email: sarah.lang...@icr.ac.uk
>
>
> The Institute of Cancer Research: Royal Cancer Hospital, a charitable Company
> Limited by Guarantee, Registered in England under Company No. 534147 with its
> Registered Office at 123 Old Brompton Road, London SW7 3RP.
>
> This e-mail message is confidential and for use by the addressee only. If
> the message is received by anyone other than the addressee, please return the
> message to the sender by replying to it and then delete the message from your
> computer and network.
>
> ------------------------------------------------------------------------------
> All the data continuously generated in your IT infrastructure contains a
> definitive record of customers, application performance, security
> threats, fraudulent activity and more. Splunk takes this data and makes
> sense of it. Business sense. IT sense. Common sense.
> http://p.sf.net/sfu/splunk-d2dcopy1
> _______________________________________________
> Rdkit-discuss mailing list
> Rdkit-discuss@lists.sourceforge.net
> https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
------------------------------------------------------------------------------
All the data continuously generated in your IT infrastructure contains a
definitive record of customers, application performance, security
threats, fraudulent activity and more. Splunk takes this data and makes
sense of it. Business sense. IT sense. Common sense.
http://p.sf.net/sfu/splunk-d2dcopy1
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss