Atom-wise distances between matching AtomGroups¶
Here we compare the distances between alpha-carbons of the enzyme adenylate kinase in its open and closed conformations. distances.dist
can be used to calculate distances between atom groups with the same number of atoms within them.
Last executed: May 18, 2021 with MDAnalysis 1.1.1
Last updated: January 2020
Minimum version of MDAnalysis: 0.19.0
Packages required:
Optional packages for visualisation:
[1]:
import MDAnalysis as mda
from MDAnalysis.tests.datafiles import PDB_small, PDB_closed
from MDAnalysis.analysis import distances
import matplotlib.pyplot as plt
%matplotlib inline
Loading files¶
The test files we will be working with here feature adenylate kinase (AdK), a phosophotransferase enzyme. ([BDPW09]) AdK has three domains:
CORE
LID: an ATP-binding domain (residues 122-159)
NMP: an AMP-binding domain (residues 30-59)
The LID and NMP domains move around the stable CORE as the enzyme transitions between the opened and closed conformations.
[2]:
u1 = mda.Universe(PDB_small) # open AdK
u2 = mda.Universe(PDB_closed) # closed AdK
Calculating the distance between CA atoms¶
We select the atoms named ‘CA’ of each Universe.
[3]:
ca1 = u1.select_atoms('name CA')
ca2 = u2.select_atoms('name CA')
distances.dist
(API docs) returns the residue numbers of both selections given. The offset
keyword adds an offset to these residue numbers to help with comparison to each other and other file formats. Here we are happy with our residue numbers, so we use the default offset of 0. (See the documentation of distances.dist
for more information.)
[4]:
resids1, resids2, dist = distances.dist(ca1, ca2,
offset=0) # for residue numbers
Plotting¶
Below, we plot the distance over the residue numbers and highlight the LID and NMP domains of the protein. The LID domain in particular moves a significant distance between its opened and closed conformations.
[5]:
plt.plot(resids1, dist)
plt.ylabel('Ca distance (Angstrom)')
plt.axvspan(122, 159, zorder=0, alpha=0.2, color='orange', label='LID')
plt.axvspan(30, 59, zorder=0, alpha=0.2, color='green', label='NMP')
plt.legend()
[5]:
<matplotlib.legend.Legend at 0x7f8f74682910>
Calculating the distance with periodic boundary conditions¶
It is common to want to calculate distances with the minimum image convention. To do this, you must pass the unitcell dimensions of the system to the box
keyword, even if your Universe has dimensions defined.
This should have the format: [lx, ly, lz, alpha, beta, gamma]
, where the first three numbers are the box lengths along each axis and the last three are the angles of the box.
[6]:
resids1_box, resids2_box, dist_box = distances.dist(ca1, ca2,
box=[10, 10, 10, 90, 90, 90])
Plotting¶
[7]:
plt.plot(resids1_box, dist_box)
plt.ylabel('Ca distance (Angstrom)')
plt.axvspan(122, 159, zorder=0, alpha=0.2, color='orange', label='LID')
plt.axvspan(30, 59, zorder=0, alpha=0.2, color='green', label='NMP')
plt.legend()
[7]:
<matplotlib.legend.Legend at 0x7f8f74699b90>
References¶
[1] Oliver Beckstein, Elizabeth J. Denning, Juan R. Perilla, and Thomas B. Woolf. Zipping and Unzipping of Adenylate Kinase: Atomistic Insights into the Ensemble of Open↔Closed Transitions. Journal of Molecular Biology, 394(1):160–176, November 2009. 00107. URL: https://linkinghub.elsevier.com/retrieve/pii/S0022283609011164, doi:10.1016/j.jmb.2009.09.009.
[2] Richard J. Gowers, Max Linke, Jonathan Barnoud, Tyler J. E. Reddy, Manuel N. Melo, Sean L. Seyler, Jan Domański, David L. Dotson, Sébastien Buchoux, Ian M. Kenney, and Oliver Beckstein. MDAnalysis: A Python Package for the Rapid Analysis of Molecular Dynamics Simulations. Proceedings of the 15th Python in Science Conference, pages 98–105, 2016. 00152. URL: https://conference.scipy.org/proceedings/scipy2016/oliver_beckstein.html, doi:10.25080/Majora-629e541a-00e.
[3] Naveen Michaud-Agrawal, Elizabeth J. Denning, Thomas B. Woolf, and Oliver Beckstein. MDAnalysis: A toolkit for the analysis of molecular dynamics simulations. Journal of Computational Chemistry, 32(10):2319–2327, July 2011. 00778. URL: http://doi.wiley.com/10.1002/jcc.21787, doi:10.1002/jcc.21787.