Skip to content

Commit

Permalink
Move angles to results.angles for dihedral classes (#3266)
Browse files Browse the repository at this point in the history
Towards #3261

## Work done in this PR
* Move angles to results.angles for dihedral classes

Co-authored-by: Oliver Beckstein <[email protected]>
  • Loading branch information
IAlibay and orbeckst authored May 6, 2021
1 parent b6ae72f commit a43c0a0
Show file tree
Hide file tree
Showing 3 changed files with 122 additions and 30 deletions.
5 changes: 5 additions & 0 deletions package/CHANGELOG
Original file line number Diff line number Diff line change
Expand Up @@ -174,6 +174,8 @@ Enhancements
checking if it can be used in parallel analysis. (Issue #2996, PR #2950)

Changes
* `analysis.dihedrals` classes now store angle data using the
`results.angles` attribute (Issue #3261)
* `msd.EinsteinMSD` now uses the `analysis.base.Results` class to store
analysis results (Issue #3261)
* `bat.BAT` now uses the `analysis.base.Results` class to store the `bat`
Expand Down Expand Up @@ -231,6 +233,9 @@ Changes
* Added OpenMM coordinate and topology converters (Issue #2863, PR #2917)

Deprecations
* The `angles` attribute for the `dihedrals` classes (Dihedral,
Ramachandran, Janin) is now deprecated in favour of `results.angles`. It
will be removed in 3.0.0 (Issue #3261)
* The `analysis.Contacts.timeseries` attribute is now deprecated in favour of
`analysis.Contacts.results.timeseries`. It will be removed in 3.0.0
(Issue #3261)
Expand Down
104 changes: 84 additions & 20 deletions package/MDAnalysis/analysis/dihedrals.py
Original file line number Diff line number Diff line change
Expand Up @@ -66,7 +66,7 @@
from MDAnalysis.analysis.dihedrals import Dihedral
R = Dihedral(ags).run()
The angles can then be accessed with :attr:`Dihedral.angles`.
The angles can then be accessed with :attr:`Dihedral.results.angles`.
Ramachandran analysis
~~~~~~~~~~~~~~~~~~~~~
Expand Down Expand Up @@ -103,7 +103,7 @@
and the "marginally allowed" regions.
To plot the data yourself, the angles can be accessed using
:attr:`Ramachandran.angles`.
:attr:`Ramachandran.results.angles`.
.. Note::
Expand All @@ -124,7 +124,7 @@
<figure-janin>` as an example.
The data for the angles can be accessed in the attribute
:attr:`Janin.angles`.
:attr:`Janin.results.angles`.
.. _figure-janin:
Expand Down Expand Up @@ -153,7 +153,8 @@
~~~~~~~~~~~~~~~
Reference plots can be added to the axes for both the Ramachandran and Janin
classes using the kwarg ``ref=True``. The Ramachandran reference data
classes using the kwarg ``ref=True`` for the :meth:`Ramachandran.plot`
and :meth:`Janin.plot` methods. The Ramachandran reference data
(:data:`~MDAnalysis.analysis.data.filenames.Rama_ref`) and Janin reference data
(:data:`~MDAnalysis.analysis.data.filenames.Janin_ref`) were made using data
obtained from a large selection of 500 PDB files, and were analyzed using these
Expand All @@ -173,32 +174,65 @@
:members:
:inherited-members:
.. attribute:: angles
.. attribute:: results.angles
Contains the time steps of the angles for each atomgroup in the list as
an ``n_frames×len(atomgroups)`` :class:`numpy.ndarray` with content
``[[angle 1, angle 2, ...], [time step 2], ...]``.
.. versionadded:: 2.0.0
.. attribute:: angles
Alias to the :attr:`results.angles` attribute.
.. deprecated:: 2.0.0
Will be removed in MDAnalysis 3.0.0. Please use
:attr:`results.angles` instead.
.. autoclass:: Ramachandran
:members:
:inherited-members:
.. attribute:: angles
.. attribute:: results.angles
Contains the time steps of the :math:`\phi` and :math:`\psi` angles for
each residue as an ``n_frames×n_residues×2`` :class:`numpy.ndarray` with
content ``[[[phi, psi], [residue 2], ...], [time step 2], ...]``.
.. versionadded:: 2.0.0
.. attribute:: angles
Alias to the :attr:`results.angles` attribute.
.. deprecated:: 2.0.0
Will be removed in MDAnalysis 3.0.0. Please use
:attr:`results.angles` instead.
.. autoclass:: Janin
:members:
:inherited-members:
.. attribute:: angles
.. attribute:: results.angles
Contains the time steps of the :math:`\chi_1` and :math:`\chi_2` angles
for each residue as an ``n_frames×n_residues×2`` :class:`numpy.ndarray`
with content ``[[[chi1, chi2], [residue 2], ...], [time step 2], ...]``.
.. versionadded:: 2.0.0
.. attribute:: angles
Alias to the :attr:`results.angles` attribute.
.. deprecated:: 2.0.0
Will be removed in MDAnalysis 3.0.0. Please use
:attr:`results.angles` instead.
References
----------
Expand Down Expand Up @@ -251,6 +285,11 @@ class Dihedral(AnalysisBase):
selection of atomgroups. If there is only one atomgroup of interest, then
it must be given as a list of one atomgroup.
.. versionchanged:: 2.0.0
:attr:`angles` results are now stored in a
:class:`MDAnalysis.analysis.base.Results` instance.
"""

def __init__(self, atomgroups, **kwargs):
Expand Down Expand Up @@ -278,16 +317,25 @@ def __init__(self, atomgroups, **kwargs):
self.ag4 = mda.AtomGroup([ag[3] for ag in atomgroups])

def _prepare(self):
self.angles = []
self.results.angles = []

def _single_frame(self):
angle = calc_dihedrals(self.ag1.positions, self.ag2.positions,
self.ag3.positions, self.ag4.positions,
box=self.ag1.dimensions)
self.angles.append(angle)
self.results.angles.append(angle)

def _conclude(self):
self.angles = np.rad2deg(np.array(self.angles))
self.results.angles = np.rad2deg(np.array(self.results.angles))

@property
def angles(self):
wmsg = ("The `angle` attribute was deprecated in MDAnalysis 2.0.0 "
"and will be removed in MDAnalysis 3.0.0. Please use "
"`results.angles` instead")
warnings.warn(wmsg, DeprecationWarning)
return self.results.angles


class Ramachandran(AnalysisBase):
r"""Calculate :math:`\phi` and :math:`\psi` dihedral angles of selected
Expand Down Expand Up @@ -349,6 +397,9 @@ class Ramachandran(AnalysisBase):
.. versionchanged:: 1.0.0
added c_name, n_name, ca_name, and check_protein keyword arguments
.. versionchanged:: 2.0.0
:attr:`angles` results are now stored in a
:class:`MDAnalysis.analysis.base.Results` instance.
"""

Expand Down Expand Up @@ -405,7 +456,7 @@ def __init__(self, atomgroup, c_name='C', n_name='N', ca_name='CA',


def _prepare(self):
self.angles = []
self.results.angles = []

def _single_frame(self):
phi_angles = calc_dihedrals(self.ag1.positions, self.ag2.positions,
Expand All @@ -415,16 +466,16 @@ def _single_frame(self):
self.ag4.positions, self.ag5.positions,
box=self.ag1.dimensions)
phi_psi = [(phi, psi) for phi, psi in zip(phi_angles, psi_angles)]
self.angles.append(phi_psi)
self.results.angles.append(phi_psi)

def _conclude(self):
self.angles = np.rad2deg(np.array(self.angles))
self.results.angles = np.rad2deg(np.array(self.results.angles))

def plot(self, ax=None, ref=False, **kwargs):
"""Plots data into standard Ramachandran plot.
Each time step in :attr:`Ramachandran.angles` is plotted onto the same
graph.
Each time step in :attr:`Ramachandran.results.angles` is plotted onto
the same graph.
Parameters
----------
Expand Down Expand Up @@ -463,10 +514,19 @@ def plot(self, ax=None, ref=False, **kwargs):
levels = [1, 17, 15000]
colors = ['#A1D4FF', '#35A1FF']
ax.contourf(X, Y, np.load(Rama_ref), levels=levels, colors=colors)
a = self.angles.reshape(np.prod(self.angles.shape[:2]), 2)
a = self.results.angles.reshape(
np.prod(self.results.angles.shape[:2]), 2)
ax.scatter(a[:, 0], a[:, 1], **kwargs)
return ax

@property
def angles(self):
wmsg = ("The `angle` attribute was deprecated in MDAnalysis 2.0.0 "
"and will be removed in MDAnalysis 3.0.0. Please use "
"`results.angles` instead")
warnings.warn(wmsg, DeprecationWarning)
return self.results.angles


class Janin(Ramachandran):
r"""Calculate :math:`\chi_1` and :math:`\chi_2` dihedral angles of selected
Expand Down Expand Up @@ -523,7 +583,9 @@ def __init__(self, atomgroup,
.. versionchanged:: 2.0.0
`select_remove` and `select_protein` keywords were added
`select_remove` and `select_protein` keywords were added.
:attr:`angles` results are now stored in a
:class:`MDAnalysis.analysis.base.Results` instance.
"""
super(Ramachandran, self).__init__(
atomgroup.universe.trajectory, **kwargs)
Expand Down Expand Up @@ -557,12 +619,13 @@ def __init__(self, atomgroup,
"missing or duplicate atoms in topology.")

def _conclude(self):
self.angles = (np.rad2deg(np.array(self.angles)) + 360) % 360
self.results.angles = (np.rad2deg(np.array(
self.results.angles)) + 360) % 360

def plot(self, ax=None, ref=False, **kwargs):
"""Plots data into standard Janin plot.
Each time step in :attr:`Janin.angles` is plotted onto the
Each time step in :attr:`Janin.results.angles` is plotted onto the
same graph.
Parameters
Expand Down Expand Up @@ -601,6 +664,7 @@ def plot(self, ax=None, ref=False, **kwargs):
levels = [1, 6, 600]
colors = ['#A1D4FF', '#35A1FF']
ax.contourf(X, Y, np.load(Janin_ref), levels=levels, colors=colors)
a = self.angles.reshape(np.prod(self.angles.shape[:2]), 2)
a = self.results.angles.reshape(np.prod(
self.results.angles.shape[:2]), 2)
ax.scatter(a[:, 0], a[:, 1], **kwargs)
return ax
43 changes: 33 additions & 10 deletions testsuite/MDAnalysisTests/analysis/test_dihedrals.py
Original file line number Diff line number Diff line change
Expand Up @@ -21,7 +21,7 @@
# J. Comput. Chem. 32 (2011), 2319--2327, doi:10.1002/jcc.21787
#
import numpy as np
from numpy.testing import assert_almost_equal
from numpy.testing import assert_almost_equal, assert_equal
import matplotlib
import pytest

Expand All @@ -46,30 +46,38 @@ def test_dihedral(self, atomgroup):
dihedral = Dihedral([atomgroup]).run()
test_dihedral = np.load(DihedralArray)

assert_almost_equal(dihedral.angles, test_dihedral, 5,
assert_almost_equal(dihedral.results.angles, test_dihedral, 5,
err_msg="error: dihedral angles should "
"match test values")

def test_dihedral_single_frame(self, atomgroup):
dihedral = Dihedral([atomgroup]).run(start=5, stop=6)
test_dihedral = [np.load(DihedralArray)[5]]

assert_almost_equal(dihedral.angles, test_dihedral, 5,
assert_almost_equal(dihedral.results.angles, test_dihedral, 5,
err_msg="error: dihedral angles should "
"match test vales")

def test_atomgroup_list(self, atomgroup):
dihedral = Dihedral([atomgroup, atomgroup]).run()
test_dihedral = np.load(DihedralsArray)

assert_almost_equal(dihedral.angles, test_dihedral, 5,
assert_almost_equal(dihedral.results.angles, test_dihedral, 5,
err_msg="error: dihedral angles should "
"match test values")

def test_enough_atoms(self, atomgroup):
with pytest.raises(ValueError):
dihedral = Dihedral([atomgroup[:2]]).run()

def test_dihedral_attr_warning(self, atomgroup):
dihedral = Dihedral([atomgroup]).run(stop=2)

wmsg = "The `angle` attribute was deprecated in MDAnalysis 2.0.0"
with pytest.warns(DeprecationWarning, match=wmsg):
assert_equal(dihedral.angles, dihedral.results.angles)


class TestRamachandran(object):

@pytest.fixture()
Expand All @@ -83,23 +91,23 @@ def rama_ref_array(self):
def test_ramachandran(self, universe, rama_ref_array):
rama = Ramachandran(universe.select_atoms("protein")).run()

assert_almost_equal(rama.angles, rama_ref_array, 5,
assert_almost_equal(rama.results.angles, rama_ref_array, 5,
err_msg="error: dihedral angles should "
"match test values")

def test_ramachandran_single_frame(self, universe, rama_ref_array):
rama = Ramachandran(universe.select_atoms("protein")).run(
start=5, stop=6)

assert_almost_equal(rama.angles[0], rama_ref_array[5], 5,
assert_almost_equal(rama.results.angles[0], rama_ref_array[5], 5,
err_msg="error: dihedral angles should "
"match test values")

def test_ramachandran_residue_selections(self, universe):
rama = Ramachandran(universe.select_atoms("resname GLY")).run()
test_rama = np.load(GLYRamaArray)

assert_almost_equal(rama.angles, test_rama, 5,
assert_almost_equal(rama.results.angles, test_rama, 5,
err_msg="error: dihedral angles should "
"match test values")

Expand Down Expand Up @@ -128,6 +136,14 @@ def test_plot(self, universe):
assert isinstance(ax, matplotlib.axes.Axes), \
"Ramachandran.plot() did not return and Axes instance"

def test_ramachandran_attr_warning(self, universe):
rama = Ramachandran(universe.select_atoms("protein")).run(stop=2)

wmsg = "The `angle` attribute was deprecated in MDAnalysis 2.0.0"
with pytest.warns(DeprecationWarning, match=wmsg):
assert_equal(rama.angles, rama.results.angles)


class TestJanin(object):

@pytest.fixture()
Expand All @@ -153,22 +169,22 @@ def _test_janin(self, u, ref_array):
janin = Janin(u.select_atoms("protein")).run()

# Test precision lowered to account for platform differences with osx
assert_almost_equal(janin.angles, ref_array, 3,
assert_almost_equal(janin.results.angles, ref_array, 3,
err_msg="error: dihedral angles should "
"match test values")

def test_janin_single_frame(self, universe, janin_ref_array):
janin = Janin(universe.select_atoms("protein")).run(start=5, stop=6)

assert_almost_equal(janin.angles[0], janin_ref_array[5], 3,
assert_almost_equal(janin.results.angles[0], janin_ref_array[5], 3,
err_msg="error: dihedral angles should "
"match test values")

def test_janin_residue_selections(self, universe):
janin = Janin(universe.select_atoms("resname LYS")).run()
test_janin = np.load(LYSJaninArray)

assert_almost_equal(janin.angles, test_janin, 3,
assert_almost_equal(janin.results.angles, test_janin, 3,
err_msg="error: dihedral angles should "
"match test values")

Expand All @@ -190,3 +206,10 @@ def test_plot(self, universe):
ax = Janin(universe.select_atoms("resid 5-10")).run().plot(ref=True)
assert isinstance(ax, matplotlib.axes.Axes), \
"Ramachandran.plot() did not return and Axes instance"

def test_janin_attr_warning(self, universe):
janin = Janin(universe.select_atoms("protein")).run(stop=2)

wmsg = "The `angle` attribute was deprecated in MDAnalysis 2.0.0"
with pytest.warns(DeprecationWarning, match=wmsg):
assert_equal(janin.angles, janin.results.angles)

0 comments on commit a43c0a0

Please sign in to comment.