Package rdkit :: Package Chem :: Module Lipinski
[hide private]
[frames] | no frames]

Source Code for Module rdkit.Chem.Lipinski

 1  # $Id: Lipinski.py 997 2009-02-25 06:12:43Z glandrum $ 
 2  # 
 3  # Copyright (C) 2001-2006 greg Landrum and Rational Discovery LLC 
 4  # 
 5  #   @@ All Rights Reserved  @@ 
 6  # 
 7  """ Calculation of Lipinski parameters for molecules 
 8   
 9  """ 
10  #from Chem import rdchem 
11  from rdkit import Chem 
12   
13  #----------------------------------- 
14  # on import build the SMARTS patterns so we only have to do it once 
15  #----------------------------------- 
16   
17  # The Daylight SMARTS expressions for 
18  # recognizing H-bond donors and acceptors in the Lipinski scheme. 
19  # HDonor     '[!#6;!H0;-0]' 
20  # HAcceptor  '[$([!#6;+0]);!$([F,Cl,Br,I]); 
21  #             !$([o,s,nX3]);!$([Nv5,Pv5,Sv4,Sv6])]' 
22  # Heteroatom '[B,N,O,P,S,F,Cl,Br,I]' 
23   
24   
25  # 2 definitions adapted from those in the Gobbi Paper 
26  #  NOTE: if you want traditional Lipinski numbers, you 
27  #  should use NOCounts (below) instead of HAcceptor 
28  # 
29  HDonorSmarts = Chem.MolFromSmarts('[$([N;!H0;v3]),$([N;!H0;+1;v4]),$([O,S;H1;+0]),$([n;H1;+0])]') 
30  # changes log for HAcceptorSmarts: 
31  #  v2, 1-Nov-2008, GL : fix amide-N exclusion; remove Fs from definition 
32  HAcceptorSmarts = Chem.MolFromSmarts('[$([O,S;H1;v2]-[!$(*=[O,N,P,S])]),\ 
33  $([O,S;H0;v2]),$([O,S;-]),\ 
34  $([N;v3;!$(N-*=!@[O,N,P,S])]),\ 
35  $([nH0,o,s;+0])\ 
36  ]') 
37  HeteroatomSmarts = Chem.MolFromSmarts('[!#6;!#1]') 
38  #  NOTE: the Rotatable bond smarts here doesn't treat deuteriums (which are left in the graph 
39  #  and therefore contribute to the degree of a carbon) the same as hydrogens (which are removed 
40  #  from the graph). So the bond in [2H]C([2H])([2H])C([2H])([2H])[2H] *is* considered 
41  #  rotatable. 
42  RotatableBondSmarts = Chem.MolFromSmarts('[!$(*#*)&!D1]-&!@[!$(*#*)&!D1]') 
43  NHOHSmarts = Chem.MolFromSmarts('[#8H1,#7H1,#7H2,#7H3]') 
44  NOCountSmarts = Chem.MolFromSmarts('[#7,#8]') 
45   
46  # this little trick saves duplicated code 
47 -def _NumMatches(mol,smarts):
48 return len(mol.GetSubstructMatches(smarts,uniquify=1))
49 50 NumHDonors = lambda x,y=HDonorSmarts:_NumMatches(x,y) 51 NumHDonors.__doc__="Number of Hydrogen Bond Donors" 52 NumHDonors.version="1.0.0" 53 _HDonors = lambda x,y=HDonorSmarts:x.GetSubstructMatches(y,uniquify=1) 54 NumHAcceptors = lambda x,y=HAcceptorSmarts:_NumMatches(x,y) 55 NumHAcceptors.__doc__="Number of Hydrogen Bond Acceptors" 56 NumHAcceptors.version="2.0.0" 57 _HAcceptors = lambda x,y=HAcceptorSmarts:x.GetSubstructMatches(y,uniquify=1) 58 NumHeteroatoms = lambda x,y=HeteroatomSmarts:_NumMatches(x,y) 59 NumHeteroatoms.__doc__="Number of Heteroatoms" 60 NumHeteroatoms.version="1.0.0" 61 _Heteroatoms = lambda x,y=HeteroatomSmarts:x.GetSubstructMatches(y,uniquify=1) 62 NumRotatableBonds = lambda x,y=RotatableBondSmarts:_NumMatches(x,y) 63 NumRotatableBonds.__doc__="Number of Rotatable Bonds" 64 NumRotatableBonds.version="1.0.0" 65 _RotatableBonds = lambda x,y=RotatableBondSmarts:x.GetSubstructMatches(y,uniquify=1) 66 NOCount = lambda x,y=NOCountSmarts:_NumMatches(x,y) 67 NOCount.__doc__="Number of Nitrogens and Oxygens" 68 NOCount.version="1.0.0" 69 NHOHCount = lambda x,y=NHOHSmarts:_NumMatches(x,y) 70 NHOHCount.__doc__="Number of NHs or OHs" 71 NHOHCount.version="1.0.0" 72 73
74 -def RingCount(mol):
75 " Number of rings a molecule has." 76 return mol.GetRingInfo().NumRings()
77 RingCount.version = "1.0.0" 78
79 -def HeavyAtomCount(mol):
80 " Number of heavy atoms a molecule." 81 return mol.GetNumAtoms(onlyHeavy=True)
82 HeavyAtomCount.version = "1.0.0" 83