Fix EnumerateStereoisomers for EITHERDOUBLE cis/trans bonds (#4272)

* Fix logic in findPotentialStereoBonds

* Strip EITHERDOUBLE bonds in EnumerateStereoisomers

* Add test for EITHERDOUBLE EnumerateStereoisomers
This commit is contained in:
Matt Swain
2021-06-23 07:20:58 +01:00
committed by GitHub
parent 91d2093dd0
commit 1d3c634691
3 changed files with 14 additions and 8 deletions

View File

@@ -2101,18 +2101,12 @@ void findPotentialStereoBonds(ROMol &mol, bool cleanIt) {
!(mol.getRingInfo()->numBondRings((*bondIt)->getIdx()))) {
// we are ignoring ring bonds here - read the FIX above
Bond *dblBond = *bondIt;
// We ignore bonds flagged as EITHERDOUBLE or STEREOANY which have
// stereo atoms set.
if (dblBond->getBondDir() == Bond::EITHERDOUBLE ||
(dblBond->getStereo() == Bond::STEREOANY &&
dblBond->getStereoAtoms().size() == 2)) {
continue;
}
// proceed only if we either want to clean the stereocode on this bond,
// if none is set on it yet, or it is STEREOANY and we need to find
// stereoatoms
if (cleanIt || dblBond->getStereo() == Bond::STEREONONE ||
dblBond->getStereo() == Bond::STEREOANY) {
(dblBond->getStereo() == Bond::STEREOANY &&
dblBond->getStereoAtoms().size() != 2)) {
dblBond->setStereo(Bond::STEREONONE);
const Atom *begAtom = dblBond->getBeginAtom(),
*endAtom = dblBond->getEndAtom();

View File

@@ -294,6 +294,9 @@ def EnumerateStereoisomers(m, options=StereoEnumerationOptions(), verbose=False)
tm = Chem.Mol(m)
for atom in tm.GetAtoms():
atom.ClearProp("_CIPCode")
for bond in tm.GetBonds():
if bond.GetBondDir() == Chem.BondDir.EITHERDOUBLE:
bond.SetBondDir(Chem.BondDir.NONE)
flippers = _getFlippers(tm, options)
nCenters = len(flippers)
if not nCenters:

View File

@@ -433,5 +433,14 @@ class TestCase(unittest.TestCase):
self.assertTrue(at.HasProp("_ChiralityPossible"))
def testEnumerateEitherDoubleStereo(self):
""" EnumerateStereoisomers from MOL with explicit either cis/trans bond """
rdbase = os.environ["RDBASE"]
filename = os.path.join(rdbase, 'Code/GraphMol/FileParsers/test_data/simple_either.mol')
mol = Chem.MolFromMolFile(filename)
smiles = [Chem.MolToSmiles(m) for m in AllChem.EnumerateStereoisomers(mol)]
self.assertEqual(set(smiles), {"C/C=C/C", "C/C=C\\C"})
if __name__ == '__main__':
unittest.main()