mdtraj.baker_hubbard

mdtraj.baker_hubbard(traj, freq=0.1, exclude_water=True, periodic=True)

Identify hydrogen bonds based on cutoffs for the Donor-H...Acceptor distance and angle.

The criterion employed is \(\theta > 120\) and \(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 [R20]. 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.

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.

See also

kabsch_sander

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”

References

[R20](1, 2) Baker, E. N., and R. E. Hubbard. “Hydrogen bonding in globular proteins.” Progress in Biophysics and Molecular Biology 44.2 (1984): 97-179.

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
Versions