Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
3 changes: 3 additions & 0 deletions package/CHANGELOG
Original file line number Diff line number Diff line change
Expand Up @@ -20,6 +20,9 @@ The rules for this file:
* 2.11.0

Fixes
* DSSP now explicitly checks for a minimum of 6 residues and raises a clear
error message, unlike the previous behavior where it would fail with an
incomprehensible broadcasting error at execution time (Issue #5046, PR #5163)
* Fixes the verbose=False in EinsteinMSD, and only shows progress bar when
verbose=True (Issue #5144, PR #5153)
* Fix incorrect assignment of topology_format to format (and vice versa)
Expand Down
10 changes: 9 additions & 1 deletion package/MDAnalysis/analysis/dssp/dssp.py
Original file line number Diff line number Diff line change
Expand Up @@ -210,7 +210,7 @@ class DSSP(AnalysisBase):
Parameters
----------
atoms : Union[Universe, AtomGroup]
input Universe or AtomGroup. In both cases, only protein residues will
input Universe or AtomGroup with at least 6 protein residues. In both cases, only protein residues will
be chosen prior to the analysis via `select_atoms('protein')`.
Heavy atoms of the protein are then selected by name
`heavyatom_names`, and hydrogens are selected by name
Expand Down Expand Up @@ -239,6 +239,8 @@ class DSSP(AnalysisBase):
Raises
------
ValueError
If fewer than 6 residues are provided in the selection.
ValueError
if ``guess_hydrogens`` is True but some non-PRO hydrogens are missing.
Expand Down Expand Up @@ -357,6 +359,12 @@ def __init__(
)
)

if len(ag.residues) < 6:
raise ValueError(
"DSSP requires at least 6 residues for secondary structure analysis, "
f"but only {len(ag.residues)} residue(s) were provided in the selection."
)

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

Expand Down
20 changes: 20 additions & 0 deletions testsuite/MDAnalysisTests/analysis/test_dssp.py
Original file line number Diff line number Diff line change
Expand Up @@ -89,3 +89,23 @@ def test_exception_raises_with_atom_index(pdb_filename, client_DSSP):
match="Residue <Residue SER, 298> contains*",
):
DSSP(u, guess_hydrogens=False).run(**client_DSSP)


def test_insufficient_residues_raises_error(client_DSSP):
"""Test that DSSP raises clear error for insufficient residues."""
u = mda.Universe(TPR, XTC)

protein = u.select_atoms("protein")

with pytest.raises(ValueError, match="DSSP requires at least 6 residues"):
res2 = protein.residues[:2].atoms
DSSP(res2)

with pytest.raises(ValueError, match="DSSP requires at least 6 residues"):
res4 = protein.residues[:4].atoms
DSSP(res4)

res6 = protein.residues[:6].atoms
dssp = DSSP(res6)
result = dssp.run(**client_DSSP, stop=1)
assert result.results.dssp.shape[1] == 6
Loading