-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathhbond.py
448 lines (370 loc) · 17.3 KB
/
hbond.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
##############################################################################
# MDTraj: A Python Library for Loading, Saving, and Manipulating
# Molecular Dynamics Trajectories.
# Copyright 2012-2013 Stanford University and the Authors
#
# Authors: Robert McGibbon
# Contributors:
#
# MDTraj is free software: you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as
# published by the Free Software Foundation, either version 2.1
# of the License, or (at your option) any later version.
#
# This library is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public
# License along with MDTraj. If not, see <http://www.gnu.org/licenses/>.
##############################################################################
##############################################################################
# Imports
##############################################################################
from __future__ import print_function, division
import numpy as np
from mdtraj.utils import ensure_type
from mdtraj.geometry import compute_distances, compute_angles
from mdtraj.geometry import _geometry
__all__ = ['wernet_nilsson', 'baker_hubbard', 'kabsch_sander']
##############################################################################
# Functions
##############################################################################
def wernet_nilsson(traj, exclude_water=True, periodic=True, sidechain_only=False):
"""Identify hydrogen bonds based on cutoffs for the Donor-H...Acceptor
distance and angle according to the criterion outlined in [1].
As opposed to Baker-Hubbard, this is a "cone" criterion where the
distance cutoff depends on the angle.
The criterion employed is :math:`r_\\text{DA} < 3.3 A - 0.00044*\\delta_{HDA}*\\delta_{HDA}`,
where :math:`r_\\text{DA}` is the distance between donor and acceptor heavy atoms,
and :math:`\\delta_{HDA}` is the angle made by the hydrogen atom, donor, and acceptor atoms,
measured in degrees (zero in the case of a perfectly straight bond: D-H ... A).
When donor the donor is 'O' and the acceptor is 'O', this corresponds to
the definition established in [1]_. The donors considered by this method
are NH and OH, and the acceptors considered are O and N. In the paper the only
donor considered is OH.
Parameters
----------
traj : md.Trajectory
An mdtraj trajectory. It must contain topology information.
exclude_water : bool, default=True
Exclude solvent molecules from consideration.
periodic : bool, default=True
Set to True to calculate displacements and angles across periodic box boundaries.
sidechain_only : bool, default=False
Set to True to only consider sidechain-sidechain interactions.
Returns
-------
hbonds : list, len=n_frames
A list containing the atom indices involved in each of the identified
hydrogen bonds at each frame. Each element in the list is an array
where each row contains three integer indices, `(d_i, h_i, a_i)`,
such that `d_i` is the index of the donor atom, `h_i` the index
of the hydrogen atom, and `a_i` the index of the acceptor atom involved
in a hydrogen bond which occurs in that frame.
Notes
-----
Each hydrogen bond is distinguished for the purpose of this function by the
indices of the donor, hydrogen, and acceptor atoms. This means that, for
example, when an ARG sidechain makes a hydrogen bond with its NH2 group,
you might see what appear like double counting of the h-bonds, since the
hydrogen bond formed via the H_1 and H_2 are counted separately, despite
their "chemical indistinguishably"
Examples
--------
>>> md.wernet_nilsson(t)
array([[ 0, 10, 8],
[ 0, 11, 7],
[ 69, 73, 54],
[ 76, 82, 65],
[119, 131, 89],
[140, 148, 265],
[166, 177, 122],
[181, 188, 231]])
>>> label = lambda hbond : '%s -- %s' % (t.topology.atom(hbond[0]), t.topology.atom(hbond[2]))
>>> for hbond in hbonds:
>>> print label(hbond)
GLU1-N -- GLU1-OE2
GLU1-N -- GLU1-OE1
GLY6-N -- SER4-O
CYS7-N -- GLY5-O
TYR11-N -- VAL8-O
MET12-N -- LYS20-O
See Also
--------
baker_hubbard, kabsch_sander
References
----------
.. [1] Wernet, Ph., L.G.M. Pettersson, and A. Nilsson, et al.
"The Structure of the First Coordination Shell in Liquid Water." (2004)
Science 304, 995-999.
"""
distance_cutoff = 0.33
angle_const = 0.000044
angle_cutoff = 45
if traj.topology is None:
raise ValueError('wernet_nilsson requires that traj contain topology '
'information')
# Get the possible donor-hydrogen...acceptor triplets
bond_triplets = _get_bond_triplets(traj.topology,
exclude_water=exclude_water, sidechain_only=sidechain_only)
# Compute geometry
mask, distances, angles = _compute_bounded_geometry(traj, bond_triplets,
distance_cutoff, [0, 2], [2, 0, 1], periodic=periodic)
# Update triplets under consideration
bond_triplets = bond_triplets.compress(mask, axis=0)
# Calculate the true cutoffs for distances
cutoffs = distance_cutoff - angle_const * (angles * 180.0 / np.pi) ** 2
# Find triplets that meet the criteria
presence = np.logical_and(distances < cutoffs, angles < angle_cutoff)
return [bond_triplets.compress(present, axis=0) for present in presence]
def baker_hubbard(traj, freq=0.0, exclude_water=True, periodic=True, sidechain_only=False):
"""Identify hydrogen bonds based on cutoffs for the Donor-H...Acceptor
distance and angle.
The criterion employed is :math:`\\theta > 120` and
:math:`r_\\text{H...Acceptor} < 2.5 A`.
When donor the donor is 'N' and the acceptor is 'O', this corresponds to
the definition established in [1]_. The donors considered by this method
are NH and OH, and the acceptors considered are O and N.
Parameters
----------
traj : md.Trajectory
An mdtraj trajectory. It must contain topology information.
freq : float, default=0.1
Return only hydrogen bonds that occur in greater this fraction of the
frames in the trajectory.
exclude_water : bool, default=True
Exclude solvent molecules from consideration
periodic : bool, default=True
Set to True to calculate displacements and angles across periodic box boundaries.
sidechain_only : bool, default=False
Set to True to only consider sidechain-sidechain interactions.
Returns
-------
hbonds : np.array, shape=[n_hbonds, 3], dtype=int
An array containing the indices atoms involved in each of the identified
hydrogen bonds. Each row contains three integer indices, `(d_i, h_i,
a_i)`, such that `d_i` is the index of the donor atom, `h_i` the index
of the hydrogen atom, and `a_i` the index of the acceptor atom involved
in a hydrogen bond which occurs (according to the definition above) in
proportion greater than `freq` of the trajectory.
Notes
-----
Each hydrogen bond is distinguished for the purpose of this function by the
indices of the donor, hydrogen, and acceptor atoms. This means that, for
example, when an ARG sidechain makes a hydrogen bond with its NH2 group,
you might see what appear like double counting of the h-bonds, since the
hydrogen bond formed via the H_1 and H_2 are counted separately, despite
their "chemical indistinguishably"
Examples
--------
>>> md.baker_hubbard(t)
array([[ 0, 10, 8],
[ 0, 11, 7],
[ 69, 73, 54],
[ 76, 82, 65],
[119, 131, 89],
[140, 148, 265],
[166, 177, 122],
[181, 188, 231]])
>>> label = lambda hbond : '%s -- %s' % (t.topology.atom(hbond[0]), t.topology.atom(hbond[2]))
>>> for hbond in hbonds:
>>> print label(hbond)
GLU1-N -- GLU1-OE2
GLU1-N -- GLU1-OE1
GLY6-N -- SER4-O
CYS7-N -- GLY5-O
TYR11-N -- VAL8-O
MET12-N -- LYS20-O
See Also
--------
kabsch_sander
References
----------
.. [1] Baker, E. N., and R. E. Hubbard. "Hydrogen bonding in globular
proteins." Progress in Biophysics and Molecular Biology
44.2 (1984): 97-179.
"""
# Cutoff criteria: these could be exposed as function arguments, or
# modified if there are better definitions than the this one based only
# on distances and angles
distance_cutoff = 2.5 # nanometers
angle_cutoff = 2.0 * np.pi / 3.0 # radians
if traj.topology is None:
raise ValueError('baker_hubbard requires that traj contain topology '
'information')
# Get the possible donor-hydrogen...acceptor triplets
bond_triplets = _get_bond_triplets(traj.topology,
exclude_water=exclude_water, sidechain_only=sidechain_only)
mask, distances, angles = _compute_bounded_geometry(traj, bond_triplets,
distance_cutoff, [1, 2], [0, 1, 2], freq=freq, periodic=periodic)
# Find triplets that meet the criteria
presence = np.logical_and(distances < distance_cutoff, angles > angle_cutoff)
mask[mask] = np.mean(presence, axis=0) > freq
return bond_triplets.compress(mask, axis=0)
def kabsch_sander(traj):
"""Compute the Kabsch-Sander hydrogen bond energy between each pair
of residues in every frame.
Hydrogen bonds are defined using an electrostatic definition, assuming
partial charges of -0.42 e and +0.20 e to the carbonyl oxygen and amide
hydrogen respectively, their opposites assigned to the carbonyl carbon
and amide nitrogen. A hydrogen bond is identified if E in the following
equation is less than -0.5 kcal/mol:
.. math::
E = 0.42 \cdot 0.2 \cdot 33.2 kcal/(mol \cdot nm) * \\
(1/r_{ON} + 1/r_{CH} - 1/r_{OH} - 1/r_{CN})
Parameters
----------
traj : md.Trajectory
An mdtraj trajectory. It must contain topology information.
Returns
-------
matrices : list of scipy.sparse.csr_matrix
The return value is a list of length equal to the number of frames
in the trajectory. Each element is an n_residues x n_residues sparse
matrix, where the existence of an entry at row `i`, column `j` with value
`x` means that there exists a hydrogen bond between a backbone CO
group at residue `i` with a backbone NH group at residue `j` whose
Kabsch-Sander energy is less than -0.5 kcal/mol (the threshold for
existence of the "bond"). The exact value of the energy is given by the
value `x`.
See Also
--------
wernet_nilsson, baker_hubbard
References
----------
.. [1] Kabsch W, Sander C (1983). "Dictionary of protein secondary structure: pattern recognition of hydrogen-bonded and geometrical features". Biopolymers 22 (12): 2577-637. dio:10.1002/bip.360221211
"""
if traj.topology is None:
raise ValueError('kabsch_sander requires topology')
import scipy.sparse
xyz, nco_indices, ca_indices, proline_indices, _ = _prep_kabsch_sander_arrays(traj)
n_residues = len(ca_indices)
hbonds = np.empty((xyz.shape[0], n_residues, 2), np.int32)
henergies = np.empty((xyz.shape[0], n_residues, 2), np.float32)
hbonds.fill(-1)
henergies.fill(np.nan)
_geometry._kabsch_sander(xyz, nco_indices, ca_indices, proline_indices,
hbonds, henergies)
# The C code returns its info in a pretty inconvenient format.
# Let's change it to a list of scipy CSR matrices.
matrices = []
hbonds_mask = (hbonds != -1)
for i in range(xyz.shape[0]):
# appologies for this cryptic code -- we need to deal with the low
# level aspects of the csr matrix format.
hbonds_frame = hbonds[i]
mask = hbonds_mask[i]
henergies_frame = henergies[i]
indptr = np.zeros(n_residues + 1, np.int32)
indptr[1:] = np.cumsum(mask.sum(axis=1))
indices = hbonds_frame[mask].flatten()
data = henergies_frame[mask].flatten()
matrices.append(scipy.sparse.csr_matrix(
(data, indices, indptr), shape=(n_residues, n_residues)).T)
return matrices
#exclude_water=False
def _get_bond_triplets(topology, exclude_water=False, sidechain_only=False):
def can_participate(atom):
# Filter waters
if exclude_water and atom.residue.is_water:
return False
# Filter non-sidechain atoms
if sidechain_only and not atom.is_sidechain:
return False
# Otherwise, accept it
return True
def get_donors(e0, e1):
# Find all matching bonds
elems = set((e0, e1))
atoms = [(one, two) for one, two in topology.bonds
if set((one.element.symbol, two.element.symbol)) == elems]
# Filter non-participating atoms
atoms = [atom for atom in atoms
if can_participate(atom[0]) and can_participate(atom[1])]
# Get indices for the remaining atoms
indices = []
for a0, a1 in atoms:
pair = (a0.index, a1.index)
# make sure to get the pair in the right order, so that the index
# for e0 comes before e1
if a0.element.symbol == e1:
pair = pair[::-1]
indices.append(pair)
return indices
nh_donors = get_donors('N', 'H')
oh_donors = get_donors('O', 'H')
xh_donors = np.array(nh_donors + oh_donors)
if len(xh_donors) == 0:
# if there are no hydrogens or protein in the trajectory, we get
# no possible pairs and return nothing
return np.zeros((0, 3), dtype=int)
acceptor_elements = frozenset(('O', 'N'))
acceptors = [a.index for a in topology.atoms
if a.element.symbol in acceptor_elements and can_participate(a)]
# Make acceptors a 2-D numpy array
acceptors = np.array(acceptors)[:, np.newaxis]
# Generate the cartesian product of the donors and acceptors
xh_donors_repeated = np.repeat(xh_donors, acceptors.shape[0], axis=0)
acceptors_tiled = np.tile(acceptors, (xh_donors.shape[0], 1))
bond_triplets = np.hstack((xh_donors_repeated, acceptors_tiled))
# Filter out self-bonds
self_bond_mask = (bond_triplets[:, 0] == bond_triplets[:, 2])
return bond_triplets[np.logical_not(self_bond_mask), :]
def _compute_bounded_geometry(traj, triplets, distance_cutoff, distance_indices,
angle_indices, freq=0.0, periodic=True):
"""
Returns a tuple include (1) the mask for triplets that fulfill the distance
criteria frequently enough, (2) the actual distances calculated, and (3) the
angles between the triplets specified by angle_indices.
"""
#print(triplets[:, distance_indices])
# First we calculate the requested distances
distances = compute_distances(traj, triplets[:, distance_indices], periodic=periodic)
# Now we discover which triplets meet the distance cutoff often enough
#temp = distances < distance_cutoff
prevalence = np.mean(distances < distance_cutoff, axis=0)
mask = prevalence > freq
# Update data structures to ignore anything that isn't possible anymore
triplets = triplets.compress(mask, axis=0)
distances = distances.compress(mask, axis=1)
# Calculate angles using the law of cosines
abc_pairs = zip(angle_indices, angle_indices[1:] + angle_indices[:1])#[0,1,2],[1,2,0]
abc_distances = []
# Calculate distances (if necessary)
for abc_pair in abc_pairs:
if set(abc_pair) == set(distance_indices):
abc_distances.append(distances)
else:
abc_distances.append(compute_distances(traj, triplets[:, abc_pair],
periodic=periodic))
# Law of cosines calculation
a, b, c = abc_distances
cosines = (a ** 2 + b ** 2 - c ** 2) / (2 * a * b)
np.clip(cosines, -1, 1, out=cosines) # avoid NaN error
angles = np.arccos(cosines)
return mask, distances, angles
def _get_or_minus1(f):
try:
return f()
except IndexError:
return -1
def _prep_kabsch_sander_arrays(traj):
xyz = ensure_type(traj.xyz, dtype=np.float32, ndim=3, name='traj.xyz',
shape=(None, None, 3), warn_on_cast=False)
ca_indices, nco_indices, is_proline, is_protein = [], [], [], []
for residue in traj.topology.residues:
ca = _get_or_minus1(lambda: [a.index for a in residue.atoms if a.name == 'CA'][0])
n = _get_or_minus1(lambda: [a.index for a in residue.atoms if a.name == 'N'][0])
c = _get_or_minus1(lambda: [a.index for a in residue.atoms if a.name == 'C'][0])
o = _get_or_minus1(lambda: [a.index for a in residue.atoms if a.name == 'O'][0])
ca_indices.append(ca)
is_proline.append(residue.name == 'PRO')
nco_indices.append([n, c, o])
is_protein.append(ca != -1 and n != -1 and c != -1 and o != -1)
nco_indices = np.array(nco_indices, np.int32)
ca_indices = np.array(ca_indices, np.int32)
proline_indices = np.array(is_proline, np.int32)
is_protein = np.array(is_protein, np.int32)
return xyz, nco_indices, ca_indices, proline_indices, is_protein