On Wed, Apr 3, 2019 at 8:26 PM Jean-Marc Nuzillard < jm.nuzill...@univ-reims.fr> wrote:
> Dear Greg, > > In [10] and Out[10] of your answer below (thanks!) incited me to believe > that the geometry of double bonds > outside of the reactive center in the reaction product is not lost for > ever, contrarily to what my example showed. > So, I went back to my example and printed the SMILES representations of > the reagent and of the product. > The double bond are there and in good shape: > C/C=C(\C=C(C)\C=C(C)\C=C/C=C(/C)C(=O)[C@]12O[C@H]1[C@@](O)(CCO)N=C2O)C(=O)OC > for the reagent > C/C=C(\C=C(C)\C=C(C)\C=C/C=C(/C)C(=O)[C@]12O[C@H]1[C@@](O)(CCO)NC2=O)C(=O)OC > for the product. > So, what I conclude is that MolToInchi does not consider the double bond > geometry in the reaction product > but that MolToSmiles does it. > Apologies... I didn't read your original email closely enough and answered the wrong question. What's going on in your case is that the stereochemistry information in the products of the reaction hasn't been updated. This doesn't happen automatically (just as the products of reactions are not automatically sanitized). If you call Chem.AssignStereochemistry() on the products you get the expected result: In [38]: inchi1 = "InChI=1S/C23H29NO7/c1-6-17(19(27)30-5)13-15(3)12-14(2)8-7-9-16(4)18(26)23-20(31-23)22(29,10-11-25)24-21(23)28/h6-9,1 ...: 2-13,20,25,29H,10-11H2,1-5H3,(H,24,28)/b8-7-,14-12+,15-13+,16-9-,17-6+/t20-,22-,23-/m0/s1" ...: m1 = Chem.MolFromInchi(inchi1) ...: rxn = AllChem.ReactionFromSmarts('[C:1]([OH:2])=[N:3]>>[C:1](=[OH0:2])[NH:3]') ...: ps = rxn.RunReactants((m1,)) ...: m2 = ps[0][0] ...: Chem.AssignStereochemistry(m2,force=True,cleanIt=True) ...: inchi2 = Chem.MolToInchi(m2) ...: print(inchi1) ...: print(inchi2) ...: print(inchi1 == inchi2) InChI=1S/C23H29NO7/c1-6-17(19(27)30-5)13-15(3)12-14(2)8-7-9-16(4)18(26)23-20(31-23)22(29,10-11-25)24-21(23)28/h6-9,12-13,20,25,29H,10-11H2,1-5H3,(H,24,28)/b8-7-,14-12+,15-13+,16-9-,17-6+/t20-,22-,23-/m0/s1 InChI=1S/C23H29NO7/c1-6-17(19(27)30-5)13-15(3)12-14(2)8-7-9-16(4)18(26)23-20(31-23)22(29,10-11-25)24-21(23)28/h6-9,12-13,20,25,29H,10-11H2,1-5H3,(H,24,28)/b8-7-,14-12+,15-13+,16-9-,17-6+/t20-,22-,23-/m0/s1 True The creation of new Mol objects from these two SMILES and the subsequent > generation of InChI strings produces identical results > as expected, because the reagent and the product are tautomers. > > By the way, does the last RDKit release correct the problem you mentioned > about MolFromSmiles > in you email dated February 14, 2019? > Part of it, but not everything. I will file a bug for the problem you reported so that it doesn't get lost. Maybe I can find/fix it before the next patch release. -greg > Best regards, > > Jean-Marc > > > > Le 03/04/2019 à 18:16, Greg Landrum a écrit : > > Hi Jean-Marc, > > Unfortunately not. The RDKit currently ignores stereo information from > "mapped bonds" (bonds for which both atoms are mapped) from the reactants. > If you include stereo info in the reaction itself, that will be copied over > in to the products, but that is not what you're looking for. Here's an > example: > > In [4]: r1 = > AllChem.ReactionFromSmarts('[C:1][C:2]=[C:3][F:4]>>[C:1][C:2]=[C:3][Cl:4]') > > > In [5]: > Chem.MolToSmiles(r1.RunReactants((Chem.MolFromSmiles('CC=CF'),))[0][0]) > > Out[5]: 'CC=CCl' > > In [6]: > Chem.MolToSmiles(r1.RunReactants((Chem.MolFromSmiles('C/C=C/F'),))[0][0]) > > Out[6]: 'CC=CCl' > > In [7]: r2 = > AllChem.ReactionFromSmarts('[C:1]/[C:2]=[C:3]/[F:4]>>[C:1]/[C:2]=[C:3]/[Cl:4]') > > > In [8]: > Chem.MolToSmiles(r2.RunReactants((Chem.MolFromSmiles('C/C=C/F'),))[0][0]) > > Out[8]: 'C/C=C/Cl' > > In [9]: > Chem.MolToSmiles(r2.RunReactants((Chem.MolFromSmiles('CC=CF'),))[0][0]) > > Out[9]: 'C/C=C/Cl' > > In [10]: > Chem.MolToSmiles(r1.RunReactants((Chem.MolFromSmiles('C/C=C/CC=CF'),))[0][0]) > > Out[10]: 'C/C=C/CC=CCl' > > > It's more tolerant of atomic stereochemistry: > > In [11]: r3 = > AllChem.ReactionFromSmarts('[C:1][C:2][F:3]>>[C:1][C:2][O:3]') > > > In [12]: > Chem.MolToSmiles(r3.RunReactants((Chem.MolFromSmiles('N[C@H](Cl)CF'),))[0][0]) > > Out[12]: 'N[C@H](Cl)CO' > > In [13]: > Chem.MolToSmiles(r3.RunReactants((Chem.MolFromSmiles('C[C@H](F)Cl'),))[0][0]) > > Out[13]: 'C[C@H](O)Cl' > > > Figuring out a good solution for the double bond handling would be a nice > challenge for the next release. > > Best, > -greg > > > > On Wed, Apr 3, 2019 at 5:54 PM Jean-Marc Nuzillard < > jm.nuzill...@univ-reims.fr> wrote: > >> Dear all, >> >> I have a molecule that comes from an InChI string, that is decoded >> as an iminol I would like to transform in its amide tautomer. >> The same molecule contains also a series of doubles bonds. >> >> The following code: >> >> inchi1 = >> "InChI=1S/C23H29NO7/c1-6-17(19(27)30-5)13-15(3)12-14(2)8-7-9-16(4)18(26)23-20(31-23)22(29,10-11-25)24-21(23)28/h6-9,12-13,20,25,29H,10-11H2,1-5H3,(H,24,28)/b8-7-,14-12+,15-13+,16-9-,17-6+/t20-,22-,23-/m0/s1" >> m1 = Chem.MolFromInchi(inchi1) >> rxn = >> AllChem.ReactionFromSmarts('[C:1]([OH:2])=[N:3]>>[C:1](=[OH0:2])[NH:3]') >> ps = rxn.RunReactants((m1,)) >> m2 = ps[0][0] >> inchi2 = Chem.MolToInchi(m2) >> print(inchi1) >> print(inchi2) >> print(inchi1 == inchi2) >> >> prints: >> [17:26:42] WARNING: Omitted undefined stereo >> >> InChI=1S/C23H29NO7/c1-6-17(19(27)30-5)13-15(3)12-14(2)8-7-9-16(4)18(26)23-20(31-23)22(29,10-11-25)24-21(23)28/h6-9,12-13,20,25,29H,10-11H2,1-5H3,(H,24,28)/b8-7-,14-12+,15-13+,16-9-,17-6+/t20-,22-,23-/m0/s1 >> >> InChI=1S/C23H29NO7/c1-6-17(19(27)30-5)13-15(3)12-14(2)8-7-9-16(4)18(26)23-20(31-23)22(29,10-11-25)24-21(23)28/h6-9,12-13,20,25,29H,10-11H2,1-5H3,(H,24,28)/t20-,22-,23-/m0/s1 >> False >> >> in which the /b layer has disappeared, leaving the double bond geometry >> undefined. >> >> Is there a way to preserve or to recover the lost information? >> >> Jean-Marc >> -- >> >> Dr. Jean-Marc Nuzillard >> Institute of Molecular Chemistry, CNRS UMR 7312 >> Faculté des Sciences Exactes et Naturelles, Bâtiment 18 >> BP 1039 >> 51687 REIMS Cedex 2 >> France >> >> Tel : 33 3 26 91 82 10 >> Fax : 33 3 26 91 31 >> 66http://www.univ-reims.fr/icmrhttp://eos.univ-reims.fr/LSD/CSNteam.html >> http://www.univ-reims.fr/LSD/http://www.univ-reims.fr/LSD/JmnSoft/ >> >> _______________________________________________ >> Rdkit-discuss mailing list >> Rdkit-discuss@lists.sourceforge.net >> https://lists.sourceforge.net/lists/listinfo/rdkit-discuss >> > > > -- > Jean-Marc Nuzillard > Directeur de Recherches au CNRS > > Institut de Chimie Moléculaire de Reims > CNRS UMR 7312 > Moulin de la Housse > CPCBAI, Bâtiment 18 > BP 1039 > 51687 REIMS Cedex 2 > France > > Tel : 03 26 91 82 10 > Fax : 03 26 91 31 > 66http://www.univ-reims.fr/icmrhttp://eos.univ-reims.fr/LSD/CSNteam.html > http://www.univ-reims.fr/LSD/http://www.univ-reims.fr/LSD/JmnSoft/ > >
_______________________________________________ Rdkit-discuss mailing list Rdkit-discuss@lists.sourceforge.net https://lists.sourceforge.net/lists/listinfo/rdkit-discuss