diff --git a/package/CHANGELOG b/package/CHANGELOG index 4fa11c2b58a..422b10feea5 100644 --- a/package/CHANGELOG +++ b/package/CHANGELOG @@ -14,6 +14,7 @@ The rules for this file: ------------------------------------------------------------------------------ + mm/dd/17 utkbansal, kain88-de, xiki-tempula, kaplajon, wouterboomsma @@ -29,6 +30,8 @@ Fixes * Various documentation sphinx errors (PR #1312) * Bugfix in confdistmatrix.get_distance_matrix; now works on all trajectory types. (issue #1324) + * Fix hbond_analysis cannot deal with Universe where no two atoms are with 3A. + (PR #1325) Changes * Enable various pylint warnings to increase python 3 compatibility diff --git a/package/MDAnalysis/analysis/hbonds/hbond_analysis.py b/package/MDAnalysis/analysis/hbonds/hbond_analysis.py index fc118762632..3065dbfbad8 100644 --- a/package/MDAnalysis/analysis/hbonds/hbond_analysis.py +++ b/package/MDAnalysis/analysis/hbonds/hbond_analysis.py @@ -823,6 +823,8 @@ def _update_selection_2(self): self._s2_donors = {} self._s2_donors_h = {} self._s2_acceptors = {} + if not self._s2: + return None if self.selection1_type in ('donor', 'both'): self._s2_acceptors = self._s2.select_atoms( 'name {0}'.format(' '.join(self.acceptors))) diff --git a/testsuite/MDAnalysisTests/analysis/test_hbonds.py b/testsuite/MDAnalysisTests/analysis/test_hbonds.py index b6cabd17c11..aca854145f2 100644 --- a/testsuite/MDAnalysisTests/analysis/test_hbonds.py +++ b/testsuite/MDAnalysisTests/analysis/test_hbonds.py @@ -31,6 +31,7 @@ import itertools import warnings +from six import StringIO from MDAnalysisTests.datafiles import PDB_helix, GRO, XTC # For type guessing: @@ -87,6 +88,17 @@ def test_generate_table(self): assert_array_equal(h.table.donor_resid, self.values['donor_resid']) assert_array_equal(h.table.acceptor_resnm, self.values['acceptor_resnm']) + @staticmethod + def test_atoms_too_far(): + pdb = ''' +ATOM 1 N LEU 1 32.310 13.778 14.372 1.00 0.00 SYST N 0 +ATOM 2 OW SOL 2 3.024 4.456 4.147 1.00 0.00 SYST H 0''' + + u = MDAnalysis.Universe(StringIO(pdb), format="pdb") + h = MDAnalysis.analysis.hbonds.HydrogenBondAnalysis(u, 'resname SOL', 'protein') + h.run(verbose=False) + assert_equal(h.timeseries, [[]]) + @staticmethod def test_true_traj(): u = MDAnalysis.Universe(GRO, XTC)