Dear Wim,

Thanks for your reply!

Apologies for the delay, finally got time to pick up this project again.

Your suggestion works great, though I have found some cases where it
breaks. For instance the molecule:

mol = Chem.MolFromSmiles("C1CC2C3C2C1C3")

It seems, in this case, a bridgehead atom is also a fused-ring atom. Maybe
these looped compounds have too complex topology for this type of analysis.

I don't see a straight way forward to identify just the bridgehead atoms.

Best wishes,
Andreas

On Sat, Dec 3, 2022 at 12:53 PM Wim Dehaen <wimdeh...@gmail.com> wrote:

> Hi Andreas,
> I don't have a good SMARTS pattern available for this but here is a
> function that should return bridgehead idx and not include non bridgehead
> fused ring atoms:
>
> ```
> def return_bridgeheads_idx(mol):
>     bh_list=[]
>     intersections=[]
>     sssr_idx = [set(x) for x in list(Chem.GetSymmSSSR(mol))]
>     for i,ring1 in enumerate(sssr_idx):
>         for j,ring2 in enumerate(sssr_idx):
>             if i>j:
>                 intersections+=[ring1.intersection(ring2)]
>     for iidx in intersections:
>         if len(iidx)>2: #condition for bridgehead
>             for idx in iidx:
>                 neighbors = [a.GetIdx() for a in
> mol.GetAtomWithIdx(idx).GetNeighbors()]
>                 bh_list+=[idx for nidx in neighbors if nidx not in iidx]
>     return tuple(set(bh_list))
> ```
>
> Here are 6 test molecules:
>
> ```
> mol1 = Chem.MolFromSmiles("C1CC2CCC1C2")
> mol2 = Chem.MolFromSmiles("C1CC2C1C1CCC2C1")
> mol3 = Chem.MolFromSmiles("N1(CC2)CCC2CC1")
> mol4 = Chem.MolFromSmiles("C1CCC12CCCCC2")
> mol5 = Chem.MolFromSmiles("C1CC2C1CCCCC2")
> mol6 = Chem.MolFromSmiles("C1CCC(C(CCC3)C23)C12")
> for mol in [mol1,mol2,mol3,mol4,mol5,mol6]:
>     print(return_bridgeheads_idx(mol))
> ```
>
> giving the expected answer:
>
> (2, 5)
> (4, 7)
> (0, 5)
> ()
> ()
> ()
>
> hope this is helpful!
>
> best wishes
> wim
>
> On Sat, Dec 3, 2022 at 8:34 AM Andreas Luttens <andreas.lutt...@gmail.com>
> wrote:
>
>> Dear users,
>>
>> I am trying to identify bridgehead atoms in multi-looped ring systems.
>> The issue I have is that it can be sometimes difficult to distinguish these
>> atoms from ring-fusion atoms. The pattern I used (see below) looks for
>> atoms that are part of three rings but cannot be bonded to an atom that
>> also fits this description, in order to avoid ring-fusion atoms. The code
>> works, except for cases where bridgehead atoms are bonded to a ring-fusion
>> atom.
>>
>> *PASS:*
>> pattern = Chem.MolFromSmarts("[$([x3]);!$([x3][x3])]")
>> rdkit_mol = Chem.MolFromSmiles("C1CC2CCC1C2")
>> print(rdkit_mol.GetSubstructMatches(pattern))
>> >>>((2,),(5,))
>>
>> *FAIL:*
>> pattern = Chem.MolFromSmarts("[$([x3]);!$([x3][x3])]")
>> rdkit_mol = Chem.MolFromSmiles("C1CC2C1C1CCC2C1")
>> print(rdkit_mol.GetSubstructMatches(pattern))
>> >>>()
>>
>> Any hint on what alternative pattern I could use to isolate true
>> bridgeheads would be greatly appreciated. Maybe other strategies are more
>> suitable to find these atoms?
>>
>> Thanks in advance!
>>
>> Best regards,
>> Andreas
>> _______________________________________________
>> Rdkit-discuss mailing list
>> Rdkit-discuss@lists.sourceforge.net
>> https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
>>
>
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss

Reply via email to