Hi Brian, The warning is actually because you have double bonds with unspecified stereochemistry. You are skipping sanitization of the molecules. When you do this no stereochemistry perception is done, so the InChI code is called without any stereochemistry information and you get the warning. If you construct the molecule "normally" (i.e. with sanitization) you get the correct InChI and no warning:
In [57]: m = Chem.MolFromSmiles(r'O=C(/C=C/c1ccccc1)c1ccc(OC/C=C(/CC/C=C(\C)/C)\C)cc1') In [58]: Chem.MolToInchi(m) Out[58]: 'InChI=1S/C25H28O2/c1-20(2)8-7-9-21(3)18-19-27-24-15-13-23(14-16-24)25(26)17-12-22-10-5-4-6-11-22/h4-6,8,10-18H,7,9,19H2,1-3H3/b17-12+,21-18+' If you really want to call the InChI code without sanitizing the molecules and want the stereochemistry to be correct, you have to do a bit more work: In [63]: m = Chem.MolFromSmiles(r'O=C(/C=C/c1ccccc1)c1ccc(OC/C=C(/CC/C=C(\C)/C)\C)cc1',sanitize=False) In [64]: m.UpdatePropertyCache(strict=False) In [65]: Chem.AssignStereochemistry(m) In [66]: Chem.MolToInchi(m) Out[66]: 'InChI=1S/C25H28O2/c1-20(2)8-7-9-21(3)18-19-27-24-15-13-23(14-16-24)25(26)17-12-22-10-5-4-6-11-22/h4-6,8,10-18H,7,9,19H2,1-3H3/b17-12+,21-18+' Best, -greg On Thu, Jun 11, 2020 at 3:46 AM Bennion, Brian via Rdkit-discuss < rdkit-discuss@lists.sourceforge.net> wrote: > Hello, > Below I show a smiles string from MOE and the smiles string calculated > from RDKit and the InChI string calculated by RDkit(2020_1). > > The error on conversion to inchi string is confusing me after entering > both smiles strings into a viewer I don't see any undefined stereo center. > > O=C(/C=C/c1ccccc1)c1ccc(OC/C=C(/CC/C=C(\C)/C)\C)cc1 > CC(C)=CCC/C(C)=C/COc1ccc(C(=O)/C=C/c2ccccc2)cc1 > [18:10:42] WARNING: Omitted undefined stereo > > InChI=1S/C25H28O2/c1-20(2)8-7-9-21(3)18-19-27-24-15-13-23(14-16-24)25(26)17-12-22-10-5-4-6-11-22/h4-6,8,10-18H,7,9,19H2,1-3H3 > > > while len(line) != 0: > fields = line.replace('","',' ').split() > mol1 = fields[0].replace('"','') > mol_name = fields[1] > > try: > mol = Chem.MolFromSmiles(mol1,sanitize=False) #, > removeHs=False) > except: > mol = None > if mol is None: > print("mol1 failed:",mol1) > output.write("mol1 failes:",mol1) > else: > rkditsmiout.write('\"'+Chem.MolToSmiles(mol, > isomericSmiles=True)+'\"\n') > print(Chem.MolToSmiles(mol, isomericSmiles=True)) > rkditsmiout.write('\"'+Chem.inchi.MolToInchi(mol)+'\"\n') > print(Chem.inchi.MolToInchi(mol)) > count += 1 > print(count) > > _______________________________________________ > 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