1
2
3
4
5
6 """Half-sphere exposure and coordination number calculation."""
7
8 import warnings
9 from math import pi
10
11 from Bio.PDB.AbstractPropertyMap import AbstractPropertyMap
12 from Bio.PDB.PDBParser import PDBParser
13 from Bio.PDB.Polypeptide import CaPPBuilder, is_aa
14 from Bio.PDB.Vector import rotaxis
15
16
18 """
19 Abstract class to calculate Half-Sphere Exposure (HSE).
20
21 The HSE can be calculated based on the CA-CB vector, or the pseudo CB-CA
22 vector based on three consecutive CA atoms. This is done by two separate
23 subclasses.
24 """
25 - def __init__(self, model, radius, offset, hse_up_key, hse_down_key,
26 angle_key=None):
27 """
28 @param model: model
29 @type model: L{Model}
30
31 @param radius: HSE radius
32 @type radius: float
33
34 @param offset: number of flanking residues that are ignored in the calculation
35 of the number of neighbors
36 @type offset: int
37
38 @param hse_up_key: key used to store HSEup in the entity.xtra attribute
39 @type hse_up_key: string
40
41 @param hse_down_key: key used to store HSEdown in the entity.xtra attribute
42 @type hse_down_key: string
43
44 @param angle_key: key used to store the angle between CA-CB and CA-pCB in
45 the entity.xtra attribute
46 @type angle_key: string
47 """
48 assert(offset>=0)
49
50 self.ca_cb_list=[]
51 ppb=CaPPBuilder()
52 ppl=ppb.build_peptides(model)
53 hse_map={}
54 hse_list=[]
55 hse_keys=[]
56 for pp1 in ppl:
57 for i in range(0, len(pp1)):
58 if i==0:
59 r1=None
60 else:
61 r1=pp1[i-1]
62 r2=pp1[i]
63 if i==len(pp1)-1:
64 r3=None
65 else:
66 r3=pp1[i+1]
67
68 result=self._get_cb(r1, r2, r3)
69 if result is None:
70
71 continue
72 pcb, angle=result
73 hse_u=0
74 hse_d=0
75 ca2=r2['CA'].get_vector()
76 for pp2 in ppl:
77 for j in range(0, len(pp2)):
78 if pp1 is pp2 and abs(i-j)<=offset:
79
80 continue
81 ro=pp2[j]
82 if not is_aa(ro) or not ro.has_id('CA'):
83 continue
84 cao=ro['CA'].get_vector()
85 d=(cao-ca2)
86 if d.norm()<radius:
87 if d.angle(pcb)<(pi/2):
88 hse_u+=1
89 else:
90 hse_d+=1
91 res_id=r2.get_id()
92 chain_id=r2.get_parent().get_id()
93
94 hse_map[(chain_id, res_id)]=(hse_u, hse_d, angle)
95 hse_list.append((r2, (hse_u, hse_d, angle)))
96 hse_keys.append((chain_id, res_id))
97
98 r2.xtra[hse_up_key]=hse_u
99 r2.xtra[hse_down_key]=hse_d
100 if angle_key:
101 r2.xtra[angle_key]=angle
102 AbstractPropertyMap.__init__(self, hse_map, hse_keys, hse_list)
103
105 """This method is provided by the subclasses to calculate HSE."""
106 return NotImplemented
107
109 """
110 Return a pseudo CB vector for a Gly residue.
111 The pseudoCB vector is centered at the origin.
112
113 CB coord=N coord rotated over -120 degrees
114 along the CA-C axis.
115 """
116 try:
117 n_v=residue["N"].get_vector()
118 c_v=residue["C"].get_vector()
119 ca_v=residue["CA"].get_vector()
120 except:
121 return None
122
123 n_v=n_v-ca_v
124 c_v=c_v-ca_v
125
126 rot=rotaxis(-pi*120.0/180.0, c_v)
127 cb_at_origin_v=n_v.left_multiply(rot)
128
129 cb_v=cb_at_origin_v+ca_v
130
131 self.ca_cb_list.append((ca_v, cb_v))
132 return cb_at_origin_v
133
134
136 """
137 Class to calculate HSE based on the approximate CA-CB vectors,
138 using three consecutive CA positions.
139 """
140 - def __init__(self, model, radius=12, offset=0):
141 """
142 @param model: the model that contains the residues
143 @type model: L{Model}
144
145 @param radius: radius of the sphere (centred at the CA atom)
146 @type radius: float
147
148 @param offset: number of flanking residues that are ignored in the calculation of the number of neighbors
149 @type offset: int
150 """
151 _AbstractHSExposure.__init__(self, model, radius, offset,
152 'EXP_HSE_A_U', 'EXP_HSE_A_D', 'EXP_CB_PCB_ANGLE')
153
155 """
156 Calculate the approximate CA-CB direction for a central
157 CA atom based on the two flanking CA positions, and the angle
158 with the real CA-CB vector.
159
160 The CA-CB vector is centered at the origin.
161
162 @param r1, r2, r3: three consecutive residues
163 @type r1, r2, r3: L{Residue}
164 """
165 if r1 is None or r3 is None:
166 return None
167 try:
168 ca1=r1['CA'].get_vector()
169 ca2=r2['CA'].get_vector()
170 ca3=r3['CA'].get_vector()
171 except:
172 return None
173
174 d1=ca2-ca1
175 d3=ca2-ca3
176 d1.normalize()
177 d3.normalize()
178
179 b=(d1+d3)
180 b.normalize()
181
182 self.ca_cb_list.append((ca2, b+ca2))
183 if r2.has_id('CB'):
184 cb=r2['CB'].get_vector()
185 cb_ca=cb-ca2
186 cb_ca.normalize()
187 angle=cb_ca.angle(b)
188 elif r2.get_resname()=='GLY':
189 cb_ca=self._get_gly_cb_vector(r2)
190 if cb_ca is None:
191 angle=None
192 else:
193 angle=cb_ca.angle(b)
194 else:
195 angle=None
196
197 return b, angle
198
200 """
201 Write a PyMol script that visualizes the pseudo CB-CA directions
202 at the CA coordinates.
203
204 @param filename: the name of the pymol script file
205 @type filename: string
206 """
207 if len(self.ca_cb_list)==0:
208 warnings.warn("Nothing to draw.", RuntimeWarning)
209 return
210 fp=open(filename, "w")
211 fp.write("from pymol.cgo import *\n")
212 fp.write("from pymol import cmd\n")
213 fp.write("obj=[\n")
214 fp.write("BEGIN, LINES,\n")
215 fp.write("COLOR, %.2f, %.2f, %.2f,\n" % (1.0, 1.0, 1.0))
216 for (ca, cb) in self.ca_cb_list:
217 x,y,z=ca.get_array()
218 fp.write("VERTEX, %.2f, %.2f, %.2f,\n" % (x,y,z))
219 x,y,z=cb.get_array()
220 fp.write("VERTEX, %.2f, %.2f, %.2f,\n" % (x,y,z))
221 fp.write("END]\n")
222 fp.write("cmd.load_cgo(obj, 'HS')\n")
223 fp.close()
224
225
227 """
228 Class to calculate HSE based on the real CA-CB vectors.
229 """
230 - def __init__(self, model, radius=12, offset=0):
231 """
232 @param model: the model that contains the residues
233 @type model: L{Model}
234
235 @param radius: radius of the sphere (centred at the CA atom)
236 @type radius: float
237
238 @param offset: number of flanking residues that are ignored in the calculation of the number of neighbors
239 @type offset: int
240 """
241 _AbstractHSExposure.__init__(self, model, radius, offset,
242 'EXP_HSE_B_U', 'EXP_HSE_B_D')
243
245 """
246 Method to calculate CB-CA vector.
247
248 @param r1, r2, r3: three consecutive residues (only r2 is used)
249 @type r1, r2, r3: L{Residue}
250 """
251 if r2.get_resname()=='GLY':
252 return self._get_gly_cb_vector(r2), 0.0
253 else:
254 if r2.has_id('CB') and r2.has_id('CA'):
255 vcb=r2['CB'].get_vector()
256 vca=r2['CA'].get_vector()
257 return (vcb-vca), 0.0
258 return None
259
260
262 - def __init__(self, model, radius=12.0, offset=0):
263 """
264 A residue's exposure is defined as the number of CA atoms around
265 that residues CA atom. A dictionary is returned that uses a L{Residue}
266 object as key, and the residue exposure as corresponding value.
267
268 @param model: the model that contains the residues
269 @type model: L{Model}
270
271 @param radius: radius of the sphere (centred at the CA atom)
272 @type radius: float
273
274 @param offset: number of flanking residues that are ignored in the calculation of the number of neighbors
275 @type offset: int
276
277 """
278 assert(offset>=0)
279 ppb=CaPPBuilder()
280 ppl=ppb.build_peptides(model)
281 fs_map={}
282 fs_list=[]
283 fs_keys=[]
284 for pp1 in ppl:
285 for i in range(0, len(pp1)):
286 fs=0
287 r1=pp1[i]
288 if not is_aa(r1) or not r1.has_id('CA'):
289 continue
290 ca1=r1['CA']
291 for pp2 in ppl:
292 for j in range(0, len(pp2)):
293 if pp1 is pp2 and abs(i-j)<=offset:
294 continue
295 r2=pp2[j]
296 if not is_aa(r2) or not r2.has_id('CA'):
297 continue
298 ca2=r2['CA']
299 d=(ca2-ca1)
300 if d<radius:
301 fs+=1
302 res_id=r1.get_id()
303 chain_id=r1.get_parent().get_id()
304
305 fs_map[(chain_id, res_id)]=fs
306 fs_list.append((r1, fs))
307 fs_keys.append((chain_id, res_id))
308
309 r1.xtra['EXP_CN']=fs
310 AbstractPropertyMap.__init__(self, fs_map, fs_keys, fs_list)
311
312
313 if __name__=="__main__":
314
315 import sys
316
317 p=PDBParser()
318 s=p.get_structure('X', sys.argv[1])
319 model=s[0]
320
321
322 RADIUS=13.0
323 OFFSET=0
324
325 hse=HSExposureCA(model, radius=RADIUS, offset=OFFSET)
326 for l in hse:
327 print l
328 print
329
330 hse=HSExposureCB(model, radius=RADIUS, offset=OFFSET)
331 for l in hse:
332 print l
333 print
334
335 hse=ExposureCN(model, radius=RADIUS, offset=OFFSET)
336 for l in hse:
337 print l
338 print
339
340 for c in model:
341 for r in c:
342 try:
343 print r.xtra['PCB_CB_ANGLE']
344 except:
345 pass
346