Hi,

On Sun, Apr 14, 2019 at 6:10 PM KENNETH PAUL RIVADENEIRA GUADAMUD <
kac...@usal.es> wrote:

> Working with some molecules and reactions. I just encounter that chiral
> centers in smiles might not being found after applying reactions.
>
> What I get after applying some reactions on a molecule is this smile:
> C[C](C)[C]1[CH+]/C=C(\\C)CC/C=C(\\C)CC1
>
> which actually seems to a have a chiral center in carbon 3 [C]. If I use
> Chem.FindMolChiralCenters(n,force=True,includeUnassigned=True) I get an
> empty list which means that there is no chiral center.
>
> The thing is that if I add H to that Carbon 3 so it becomes [CH] it is
> recognized as chiral center but with unassigned type (R or S). I tried
> adding Hs using Chem.AddHs(mol) and then try again
> Chem.FindMolChiralCenters() but didn't get any chiral center.
>
>
I think you might have misunderstood the SMILES syntax. At atom in square
brackets, like [C], has no implicit Hs assigned to it. In your case it
means that that Carbon 3 only has three neighbors (and a radical electron).
The RDKit does not consider this to be a potential stereocenter:

In [7]: m = Chem.MolFromSmiles('C[C](C)[C]1[CH+]/C=C(\\C)CC/C=C(\\C)CC1')

In [8]: Chem.FindMolChiralCenters(m,includeUnassigned=True)

Out[8]: []

If you avoid putting the atom in square brackets (and thus telling the
RDKit that it has no implicit Hs), things work as you would expect:

In [9]: m = Chem.MolFromSmiles('CC(C)C1[CH+]/C=C(\\C)CC/C=C(\\C)CC1')

In [10]: Chem.FindMolChiralCenters(m,includeUnassigned=True)

Out[10]: [(3, '?')]

Afer applying two 1,2 hydride shift to my initial mol (
> Chem.MolFromSmiles('C/C1=C\\C[C@H]([C+](C)C)CC/C(C)=C/CC1')) I get the
> smiles mentioned before. So given that I had some initial chiral tag I want
> to know if there is a way to recover lost chirality after reactions.
>
> smarts used for 1,2 hydride shift: [Ch:1]-[C+1:2]>>[C+1:1]-[Ch+0:2]
>
> mol = Chem.MolFromSmiles('C/C1=C\\C[C@H]([C+](C)C)CC/C(C)=C/CC1')
> rxn  = AllChem.ReactionFromSmarts('[Ch:1]-[C+1:2]>>[C+1:1]-[Ch+0:2]')
> products = list()for product in rxn.RunReactant(mol, 0):
>     Chem.SanitizeMol(product[0])
>     products.append(product[0])print(Chem.MolToSmiles(products[0]))
>
> After applying this reaction twice to the product created I eventually get
> this smile.
>
> Output:'C[C](C)[C]1[CH+]/C=C(\\C)CC/C=C(\\C)CC1'
>
> which actually is where it is supposed to be a chiral center in carbon 3
>
> Any idea or should I report it as a bug?
>
You haven't quite defined the reaction that you think you have.
This turns out to be somewhat tricky, mainly due to the somewhat
unsatisfying way Hs are handled in the RDKit's Reaction SMARTS.
The first thing is that I would almost never use the 'h' primitive in a
reaction SMARTS. What you almost always care about the is the number of Hs
on an atom, not the number of implicit Hs, so use the `H` primitive instead.

That still doesn't make your reaction quite right though:

In [26]: mol = Chem.MolFromSmiles('OC(F)[C+](Cl)Br')


In [27]: rxn  =
AllChem.ReactionFromSmarts('[C!H0:1]-[C+1:2]>>[C+1:1]-[C+0:2]')


In [28]: for product in rxn.RunReactant(mol,0):
    ...:     Chem.SanitizeMol(product[0])
    ...:     print(Chem.MolToSmiles(product[0]))
    ...:

O[C+](F)[C](Cl)Br


Here the charge has been transferred, but the number of Hs on C:2 hasn't
been updated.

I can fix it for this example as follows:

In [29]: rxn2  =
AllChem.ReactionFromSmarts('[C!H0:1]-[C+1:2]>>[C+1:1]-[C+0H:2]')


In [30]: for product in rxn2.RunReactant(mol,0):
    ...:     Chem.SanitizeMol(product[0])
    ...:     print(Chem.MolToSmiles(product[0]))
    ...:

O[C+](F)C(Cl)Br



But that only works correctly if C:2 doesn't start out with any Hs on it:

In [31]: mol2 = Chem.MolFromSmiles('OC[CH+]Br')


In [32]: for product in rxn2.RunReactant(mol2,0):
    ...:     Chem.SanitizeMol(product[0])
    ...:     print(Chem.MolToSmiles(product[0]))
    ...:

O[CH+][CH]Br


This is solvable by telling the RDKit that it's allowed to recalculate the
number of implicit Hs on carbon atoms:

In [45]: for product in rxn2.RunReactant(mol,0):
    ...:     p = product[0]
    ...:     for atom in p.GetAtoms():
    ...:         if atom.GetAtomicNum()==6:
    ...:             atom.SetNoImplicit(False)
    ...:     Chem.SanitizeMol(p)
    ...:     print(Chem.MolToSmiles(p))
    ...:
    ...:

O[C+](F)C(Cl)Br

In [46]: mol = Chem.MolFromSmiles('OC(F)[CH+]Br')


In [47]: for product in rxn2.RunReactant(mol,0):
    ...:     p = product[0]
    ...:     for atom in p.GetAtoms():
    ...:         if atom.GetAtomicNum()==6:
    ...:             atom.SetNoImplicit(False)
    ...:     Chem.SanitizeMol(p)
    ...:     print(Chem.MolToSmiles(p))
    ...:

O[C+](F)CBr


This will fail (I think) though if you have a reaction that's intended to
produce radicals. Ah well... you can't win them all.

Another approach is to work with molecules that have Hs added to them:

In [33]: rxn2  =
AllChem.ReactionFromSmarts('[#1:3][C:1]-[C+1:2]>>[C+1:1]-[C+0:2]-[#1:3]')


In [34]: mol3 = Chem.AddHs(Chem.MolFromSmiles('OC(F)[C+](Cl)Br'))


In [35]: for product in rxn2.RunReactant(mol3,0):
    ...:     p = Chem.RemoveHs(product[0])
    ...:     print(Chem.MolToSmiles(p))
    ...:
    ...:

O[C+](F)C(Cl)Br

In [36]: mol4 = Chem.AddHs(Chem.MolFromSmiles('OC(F)[CH+]Br'))


In [37]: for product in rxn2.RunReactant(mol4,0):
    ...:     p = Chem.RemoveHs(product[0])
    ...:     print(Chem.MolToSmiles(p))
    ...:

O[C+](F)CBr


If you wanted to apply the reaction multiple times, you'd call
Chem.SanitizeMol() on the products instead of RemoveHs().

Note that this is not the greatest general strategy for dealing with
reactions since it can very quickly lead to combinatorial explosions in the
number of products that you get.

-greg
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss

Reply via email to