Skip to content

Commit d14978f

Browse files
committed
Merge branch 'inchiKeyAug'.
Thanks @nickvandewiele. This closes #75.
2 parents 9d6b030 + fdd117e commit d14978f

File tree

3 files changed

+82
-1
lines changed

3 files changed

+82
-1
lines changed

rmgpy/molecule.pxd

Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -146,6 +146,12 @@ cdef class Molecule(Graph):
146146

147147
cpdef str toInChI(self)
148148

149+
cpdef str toAugmentedInChI(self)
150+
151+
cpdef str toInChIKey(self)
152+
153+
cpdef str toAugmentedInChIKey(self)
154+
149155
cpdef str toSMILES(self)
150156

151157
cpdef toOBMol(self)

rmgpy/molecule.py

Lines changed: 46 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1035,7 +1035,52 @@ def toInChI(self):
10351035
obConversion.SetOutFormat('inchi')
10361036
obConversion.SetOptions('w', openbabel.OBConversion.OUTOPTIONS)
10371037
return obConversion.WriteString(obmol).strip()
1038-
1038+
1039+
def toAugmentedInChI(self):
1040+
"""
1041+
Adds an extra layer to the InChI denoting the number of unpaired electrons in case
1042+
more than 1 ( >= 2) unpaired electrons are present in the molecule.
1043+
"""
1044+
inchi = self.toInChI()
1045+
1046+
radicalNumber = sum([i.radicalElectrons for i in self.atoms])
1047+
1048+
if radicalNumber >= 2:
1049+
return inchi+'/mult'+str(radicalNumber+1)
1050+
else:
1051+
return inchi
1052+
1053+
def toInChIKey(self):
1054+
"""
1055+
Convert a molecular structure to an InChI Key string. Uses
1056+
`OpenBabel <http://openbabel.org/>`_ to perform the conversion.
1057+
1058+
Removes check-sum dash (-) and character so that only
1059+
the 14 + 9 characters remain.
1060+
"""
1061+
import openbabel
1062+
# This version does not write a warning to stderr if stereochemistry is undefined
1063+
obmol = self.toOBMol()
1064+
obConversion = openbabel.OBConversion()
1065+
obConversion.SetOutFormat('inchi')
1066+
obConversion.SetOptions('w', openbabel.OBConversion.OUTOPTIONS)
1067+
obConversion.SetOptions('K', openbabel.OBConversion.OUTOPTIONS)
1068+
return obConversion.WriteString(obmol).strip()[:-2]
1069+
1070+
def toAugmentedInChIKey(self):
1071+
"""
1072+
Adds an extra layer to the InChIKey denoting the number of unpaired electrons in case
1073+
more than 1 ( >= 2) unpaired electrons are present in the molecule.
1074+
"""
1075+
key = self.toInChIKey()
1076+
1077+
radicalNumber = sum([i.radicalElectrons for i in self.atoms])
1078+
1079+
if radicalNumber >= 2:
1080+
return key+'mult'+str(radicalNumber+1)
1081+
else:
1082+
return key
1083+
10391084
def toSMILES(self):
10401085
"""
10411086
Convert a molecular structure to an SMILES string. Uses

unittest/moleculeTest.py

Lines changed: 30 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -872,6 +872,36 @@ def testRadicalCH2(self):
872872
molecule = Molecule().fromSMILES('[CH2]')
873873
self.assertEqual(molecule.atoms[0].radicalElectrons, 2)
874874
self.assertEqual(molecule.atoms[0].spinMultiplicity, 3)
875+
876+
def testInChIKey(self):
877+
"""
878+
Test that InChI Key generation is working properly.
879+
"""
880+
molecule = Molecule().fromInChI('InChI=1S/C7H12/c1-2-7-4-3-6(1)5-7/h6-7H,1-5H2')
881+
key = molecule.toInChIKey()
882+
self.assertEqual(key, 'UMRZSTCPUPJPOJ-UHFFFAOYSA')
883+
884+
def testAugmentedInChI(self):
885+
"""
886+
Test that Augmented InChI generation is printing the /mult layer
887+
"""
888+
mol = Molecule().fromAdjacencyList("""
889+
1 C 1 {2,S}
890+
2 C 1 {1,S}
891+
""")
892+
893+
self.assertEqual(mol.toAugmentedInChI(), 'InChI=1S/C2H4/c1-2/h1-2H2/mult3')
894+
895+
def testAugmentedInChIKey(self):
896+
"""
897+
Test that Augmented InChI Key generation is printing the mult layer
898+
"""
899+
mol = Molecule().fromAdjacencyList("""
900+
1 C 1 {2,S}
901+
2 C 1 {1,S}
902+
""")
903+
904+
self.assertEqual(mol.toAugmentedInChIKey(), 'VGGSQFUCUMXWEO-UHFFFAOYSAmult3')
875905
################################################################################
876906

877907
class TestMoleculeSymmetry(unittest.TestCase):

0 commit comments

Comments
 (0)