Skip to content

Commit

Permalink
Feat: Add n_initial_contacts attribute and document Contacts
Browse files Browse the repository at this point in the history
…class attributes (#4415)
  • Loading branch information
HeetVekariya authored Jan 30, 2024
1 parent e189a90 commit 5a992d7
Show file tree
Hide file tree
Showing 3 changed files with 22 additions and 0 deletions.
2 changes: 2 additions & 0 deletions package/CHANGELOG
Original file line number Diff line number Diff line change
Expand Up @@ -27,6 +27,8 @@ Fixes
* Fix deploy action to use the correct version of the pypi upload action.

Enhancements
* Documented the r0 attribute in the `Contacts` class and added the
`n_initial_contacts` attribute, with documentation. (Issue #2604, PR #4415)

Changes
* As per NEP29, the minimum version of numpy has been raised to 1.23.
Expand Down
10 changes: 10 additions & 0 deletions package/MDAnalysis/analysis/contacts.py
Original file line number Diff line number Diff line change
Expand Up @@ -405,6 +405,13 @@ def __init__(self, u, select, refgroup, method="hard_cut", radius=4.5,
Show detailed progress of the calculation if set to ``True``; the
default is ``False``.
Attributes
----------
n_initial_contacts : int
Total number of initial contacts.
r0 : list[numpy.ndarray]
List of distance arrays between reference groups.
Notes
-----
Expand Down Expand Up @@ -455,6 +462,9 @@ def __init__(self, u, select, refgroup, method="hard_cut", radius=4.5,
box=self._get_box(refA.universe)))
self.initial_contacts.append(contact_matrix(self.r0[-1], radius))

self.n_initial_contacts = self.initial_contacts[0].sum()


@staticmethod
def _get_atomgroup(u, sel):
select_error_message = ("selection must be either string or a "
Expand Down
10 changes: 10 additions & 0 deletions testsuite/MDAnalysisTests/analysis/test_contacts.py
Original file line number Diff line number Diff line change
Expand Up @@ -375,6 +375,16 @@ def test_warn_deprecated_attr(self, universe):
with pytest.warns(DeprecationWarning, match=wmsg):
assert_equal(CA1.timeseries, CA1.results.timeseries)

@pytest.mark.parametrize("datafiles, expected", [((PSF, DCD), 0),
([TPR, XTC], 41814)])
def test_n_initial_contacts(self, datafiles, expected):
"""Test for n_initial_contacts attribute"""
u = mda.Universe(*datafiles)
select = ('protein', 'not protein')
refgroup = (u.select_atoms('protein'), u.select_atoms('not protein'))

r = contacts.Contacts(u, select=select, refgroup=refgroup)
assert_equal(r.n_initial_contacts, expected)

def test_q1q2():
u = mda.Universe(PSF, DCD)
Expand Down

0 comments on commit 5a992d7

Please sign in to comment.