diff --git a/package/AUTHORS b/package/AUTHORS index 6a32bc9ce33..00076883aa4 100644 --- a/package/AUTHORS +++ b/package/AUTHORS @@ -157,6 +157,7 @@ Chronological list of authors - Leonardo Barneschi - Henrik Jäger - Jan Stevens + - Orion Cohen External code ------------- diff --git a/package/CHANGELOG b/package/CHANGELOG index bf10ef51bc4..0de5b4a35a4 100644 --- a/package/CHANGELOG +++ b/package/CHANGELOG @@ -17,11 +17,13 @@ The rules for this file: lilyminium, daveminh, jbarnoud, yuxuanzhuang, VOD555, ianmkenney, calcraven,xiki-tempula, mieczyslaw, manuel.nuno.melo, PicoCentauri, hanatok, rmeli, aditya-kamath, tirkarthi, LeonardoBarneschi, hejamu, - biogen98 + biogen98, orioncohen * 2.0.0 Fixes + * Fixed 'sphzone', 'sphlayer', 'cyzone', and 'cylayer' to return empty if the + zone/layer is empty, consistent with 'around' (Issue #2915) * A Universe created from an ROMol with no atoms returns now a Universe with 0 atoms (Issue #3142) * ValueError raised when empty atomgroup is given to DensityAnalysis diff --git a/package/MDAnalysis/core/selection.py b/package/MDAnalysis/core/selection.py index adff7ab3eeb..d72ec4f6c03 100644 --- a/package/MDAnalysis/core/selection.py +++ b/package/MDAnalysis/core/selection.py @@ -217,6 +217,7 @@ def apply(self, group): return func(self, group) return apply + class _Selectionmeta(type): def __init__(cls, name, bases, classdict): type.__init__(type, name, bases, classdict) @@ -352,6 +353,8 @@ def __init__(self, parser, tokens): def apply(self, group): indices = [] sel = self.sel.apply(group) + if len(sel) == 0: + return group[[]] box = self.validate_dimensions(group.dimensions) periodic = box is not None ref = sel.center_of_geometry().reshape(1, 3).astype(np.float32) @@ -378,6 +381,8 @@ def __init__(self, parser, tokens): def apply(self, group): indices = [] sel = self.sel.apply(group) + if len(sel) == 0: + return group[[]] box = self.validate_dimensions(group.dimensions) periodic = box is not None ref = sel.center_of_geometry().reshape(1, 3).astype(np.float32) @@ -394,7 +399,8 @@ class CylindricalSelection(Selection): @return_empty_on_apply def apply(self, group): sel = self.sel.apply(group) - + if len(sel) == 0: + return group[[]] # Calculate vectors between point of interest and our group vecs = group.positions - sel.center_of_geometry() diff --git a/testsuite/MDAnalysisTests/core/test_atomselections.py b/testsuite/MDAnalysisTests/core/test_atomselections.py index 6d52ccf0a8f..eb02538e7b0 100644 --- a/testsuite/MDAnalysisTests/core/test_atomselections.py +++ b/testsuite/MDAnalysisTests/core/test_atomselections.py @@ -210,6 +210,10 @@ def test_cylayer(self, universe, selstr): sel = universe.select_atoms(selstr) assert_equal(len(sel), 88) + def test_empty_cylayer(self, universe): + empty = universe.select_atoms('cylayer 4.0 6.0 10 -10 name NOT_A_NAME') + assert_equal(len(empty), 0) + @pytest.mark.parametrize('selstr', [ 'cyzone 6.0 10 -10 bynum 1281', 'cyzone 6.0 10 -10 index 1280' @@ -218,6 +222,10 @@ def test_cyzone(self, universe, selstr): sel = universe.select_atoms(selstr) assert_equal(len(sel), 166) + def test_empty_cyzone(self, universe): + empty = universe.select_atoms('cyzone 6.0 10 -10 name NOT_A_NAME') + assert_equal(len(empty), 0) + def test_point(self, universe): ag = universe.select_atoms('point 5.0 5.0 5.0 3.5') @@ -770,6 +778,10 @@ def test_sphlayer(self, u): assert idx == set(ag.indices) + def test_empty_sphlayer(self, u): + empty = u.select_atoms('sphlayer 2.4 6.0 name NOT_A_NAME') + assert len(empty) == 0 + def test_sphzone(self, u): r1 = u.select_atoms('resid 1') cog = r1.center_of_geometry().reshape(1, 3) @@ -781,6 +793,10 @@ def test_sphzone(self, u): assert idx == set(ag.indices) + def test_empty_sphzone(self, u): + empty = u.select_atoms('sphzone 5.0 name NOT_A_NAME') + assert len(empty) == 0 + def test_point_1(self, u): # The example selection ag = u.select_atoms('point 5.0 5.0 5.0 3.5')