The approach used in Nick's solution - using reactions to put the molecules
together - will certainly work. It's very general and flexible.
For something this specific, you could also use Chem.CombineMols and an
RWMol. Here's a simple demonstration of that:
In [16]: m1 = Chem.MolFromSmiles('c1ccccc1[*:1]')
In [17]: m2 = Chem.MolFromSmiles('CCO[*:1]')
In [18]: m3 = Chem.CombineMols(m1,m2)
In [19]: em = Chem.RWMol(m3)
In [22]: em.AddBond(9,5,Chem.BondType.SINGLE)
Out[22]: 11
In [23]: em.RemoveAtom(10)
In [24]: em.RemoveAtom(6)
In [25]: Chem.SanitizeMol(em)
Out[25]: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_NONE
In [26]: Chem.MolToSmiles(em)
OUT: 0
Out[26]: 'CCOc1ccccc1'
This requires, of course, a bit more code to identify where bonds should be
formed and the atoms that should be deleted (be sure you do this in order
of decreasing index).
Best,
-greg
On Fri, Jul 3, 2015 at 9:19 AM, Nicholas Firth <[email protected]>
wrote:
> Hi Axel,
>
> I was sure there was something in the core that did this but I can't seem
> to find it now! I have a function that uses isotope labelled dummy atoms to
> attach fragments if it's of any use? As with most the code I write, it
> doesn't have the greatest error handling, or rational use of deepcopy.
>
> def connectAllFragsIsomeric(molList):
> numFrags = len(molList)
> isotopeList = []
> allIsos = []
> # This loop is to make a list of molecules along with which isotopes
> they contain
> # so I can keep track of where which R-Groups are later
> for i in range(numFrags):
> molList[i] = processToIsomericMol(molList[i])
> for atom in molList[i].GetAtoms():
> isotope = atom.GetIsotope()
> if(isotope):
> isotopeList.append(isotope)
> allIsos.append(isotope)
> molList[i] = [molList[i], deepcopy(isotopeList)]
> del isotopeList[:]
>
> allIsos.sort()
> allIsos = list(set(allIsos))
>
> for i in allIsos:
> reactantIdx = []
> # Find two fragments with the labeled R-group
> for j in range(len(molList)):
> if(i in molList[j][1]):
> reactantIdx.append(j)
> # Sanity check that there are only two frags with the same R-group
> if(len(reactantIdx) != 2):
> print 'Too many fragments with a R-groups numbered with %i' % i
> exit()
> rxn =
> AllChem.ReactionFromSmarts('[a,A:1]-[%i*:2].[a,A:3]-[%i*:4]>>[a,A:1]-[a,A:3].[%i*:2][%i*:4]'
> % (i, i, i, i))
> ps = rxn.RunReactants((molList[reactantIdx[0]][0],
> molList[reactantIdx[1]][0]))
> newMol = ps[0][0]
> newIsos = list(set(molList[reactantIdx[0]][1]) |
> set(molList[reactantIdx[1]][1]))
> newIsos.remove(i)
> molList.append([newMol, newIsos])
> del molList[reactantIdx[1]]
> del molList[reactantIdx[0]]
> return molList
>
>
>
> Best,
> Nick
>
> Nicholas C. Firth | PhD Student | Cancer Therapeutics
> The Institute of Cancer Research | 15 Cotswold Road | Belmont | Sutton |
> Surrey | SM2 5NG
> T 020 8722 4033 | E [email protected] | W www.icr.ac.uk | Twitter
> @ICRnews
>
> ________________________________________
> From: Axel Pahl [[email protected]]
> Sent: 03 July 2015 07:27
> To: RDKit Discuss
> Subject: [Rdkit-discuss] create a mol from core and mol fragments
>
> Dear fellow RDKitters,
>
> I use ReplaceCore and GetMolFrags to create a 2D SAR table:
>
>
> tmp = Chem.ReplaceCore(mol, core_mol, labelByIndex=True)
>
> frag_mols = list(Chem.GetMolFrags(tmp, asMols=True))
>
>
> Is it also possible to do the reverse? Create a new molecule from the core
> and two fragments?
> I would like to fill the gaps in the SAR table and calculate some
> properties for these molecules.
>
> Kind regards,
> Axel
>
> 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.
>
>
> ------------------------------------------------------------------------------
> Don't Limit Your Business. Reach for the Cloud.
> GigeNET's Cloud Solutions provide you with the tools and support that
> you need to offload your IT needs and focus on growing your business.
> Configured For All Businesses. Start Your Cloud Today.
> https://www.gigenetcloud.com/
> _______________________________________________
> Rdkit-discuss mailing list
> [email protected]
> https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
>
------------------------------------------------------------------------------
Don't Limit Your Business. Reach for the Cloud.
GigeNET's Cloud Solutions provide you with the tools and support that
you need to offload your IT needs and focus on growing your business.
Configured For All Businesses. Start Your Cloud Today.
https://www.gigenetcloud.com/
_______________________________________________
Rdkit-discuss mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss