1
2
3
4
5
6
7
8
9
10
11 """ uses pymol to interact with molecules
12
13 """
14 from rdkit import Chem
15 import xmlrpclib,os
16
17
18 _server=None
20 - def __init__(self,host=None,port=9123,force=0,**kwargs):
21 global _server
22 if not force and _server is not None:
23 self.server=_server
24 else:
25 if not host:
26 host=os.environ.get('PYMOL_RPCHOST','localhost')
27 _server=None
28 serv = xmlrpclib.Server('http://%s:%d'%(host,port))
29 serv.ping()
30 _server = serv
31 self.server=serv
32 self.InitializePyMol()
33
35 """ does some initializations to set up PyMol according to our
36 tastes
37
38 """
39 self.server.do('set valence,1')
40 self.server.do('set stick_rad,0.15')
41 self.server.do('set mouse_selection_mode,0')
42 self.server.do('set line_width,2')
43 self.server.do('set selection_width,10')
44 self.server.do('set auto_zoom,0')
45
46
48 " blows out everything in the viewer "
49 self.server.deleteAll()
50
52 " deletes everything except the items in the provided list of arguments "
53 allNames = self.server.getNames('*',False)
54 for nm in allNames:
55 if nm not in excludes:
56 self.server.deleteObject(nm)
57
58 - def LoadFile(self,filename,name,showOnly=False):
59 """ calls pymol's "load" command on the given filename; the loaded object
60 is assigned the name "name"
61 """
62 if showOnly:
63 self.DeleteAll()
64 id = self.server.loadFile(filename,name)
65 return id
66
67 - def ShowMol(self,mol,name='molecule',showOnly=True,highlightFeatures=[],
68 molB="",confId=-1,zoom=True):
69 """ special case for displaying a molecule or mol block """
70
71 if not molB:
72 molB = Chem.MolToMolBlock(mol,confId=confId)
73 server = self.server
74 if not zoom:
75 self.server.do('view rdinterface,store')
76 if showOnly:
77 self.DeleteAll()
78 id = server.loadMolBlock(molB,name)
79 if highlightFeatures:
80 nm = name+'-features'
81 conf = mol.GetConformer(confId)
82 for feat in highlightFeatures:
83 pt = [0.0,0.0,0.0]
84 for idx in feat:
85 loc = conf.GetAtomPosition(idx)
86 pt[0] += loc[0]/len(feat)
87 pt[1] += loc[1]/len(feat)
88 pt[2] += loc[2]/len(feat)
89 server.sphere(pt,0.2,(1,1,1),nm)
90 if zoom:
91 server.zoom('visible')
92 else:
93 self.server.do('view rdinterface,recall')
94 return id
95
97 " returns the selected atoms "
98 if not whichSelection:
99 sels = self.server.getNames('selections')
100 if sels:
101 whichSelection = sels[-1]
102 else:
103 whichSelection=None
104 if whichSelection:
105 items = self.server.index(whichSelection)
106 else:
107 items = []
108 return items
109
110
111 - def SelectAtoms(self,itemId,atomIndices,selName='selection'):
112 " selects a set of atoms "
113 ids = '(id '
114 ids += ','.join(['%d'%(x+1) for x in atomIndices])
115 ids += ')'
116 cmd = 'select %s,%s and %s'%(selName,ids,itemId)
117 self.server.do(cmd)
118
120 " highlights a set of atoms "
121 if extraHighlight:
122 idxText = ','.join(['%s and (id %d)'%(where,x) for x in indices])
123 self.server.do('edit %s'%idxText)
124 else:
125 idxText = ' or '.join(['id %d'%x for x in indices])
126 self.server.do('select selection, %s and (%s)'%(where,idxText))
127
129 " change the display style of the specified object "
130 self.server.do('hide everything,%s'%(obj,))
131 if style:
132 self.server.do('show %s,%s'%(style,obj))
133
136 """ selects the area of a protein around a specified object/selection name;
137 optionally adds a surface to that """
138 self.server.do('select %(name)s,byres (%(aroundObj)s around %(distance)f) and %(inObj)s'%locals())
139
140
141 if showSurface:
142 self.server.do('show surface,%s'%name)
143 self.server.do('disable %s'%name)
144
146 " adds a set of spheres "
147 self.server.do('view rdinterface,store')
148 self.server.resetCGO(label)
149 for i,loc in enumerate(locs):
150 self.server.sphere(loc,sphereRad,colors[i],label,1)
151 self.server.do('enable %s'%label)
152 self.server.do('view rdinterface,recall')
153
154
156 if not val:
157 self.server.do('set defer_update,1')
158 else:
159 self.server.do('set defer_update,0')
160
162 " returns the coordinates of the selected atoms "
163 res = {}
164 for label,idx in sels:
165 coords = self.server.getAtomCoords('(%s and id %d)'%(label,idx))
166 res[(label,idx)] = coords
167 return res
168
170 self.server.do('disable all')
172 self.server.do('disable %s'%objName)
174 self.server.do('enable %s'%objName)
175
177 self.server.do('refresh')
178 - def Zoom(self,objName):
179 self.server.zoom(objName)
180
181 - def DisplayHBonds(self,objName,molName,proteinName,
182 molSelText='(%(molName)s)',
183 proteinSelText='(%(proteinName)s and not het)'):
184 " toggles display of h bonds between the protein and a specified molecule "
185 cmd = "delete %(objName)s;\n"
186 cmd += "dist %(objName)s," + molSelText+","+proteinSelText+",mode=2;\n"
187 cmd += "enable %(objName)s;"
188 cmd = cmd%locals()
189
190 self.server.do(cmd)
191
192 - def DisplayCollisions(self,objName,molName,proteinName,distCutoff=3.0,
193 color='red',
194 molSelText='(%(molName)s)',
195 proteinSelText='(%(proteinName)s and not het)'):
196 " toggles display of collisions between the protein and a specified molecule "
197 cmd = "delete %(objName)s;\n"
198 cmd += "dist %(objName)s," + molSelText+","+proteinSelText+",%(distCutoff)f,mode=0;\n"
199 cmd += """enable %(objName)s
200 color %(color)s, %(objName)s"""
201 cmd = cmd%locals()
202 self.server.do(cmd)
203