Skip to content

Commit

Permalink
Merge pull request #620 from astrofrog/fix-mask-channels
Browse files Browse the repository at this point in the history
 Fix bug in VaryingResolutionSpectralCube.mask_channels
  • Loading branch information
astrofrog committed Apr 10, 2020
2 parents 3aa15a1 + fae39fa commit 4c23f5b
Show file tree
Hide file tree
Showing 3 changed files with 27 additions and 7 deletions.
2 changes: 2 additions & 0 deletions CHANGES.rst
@@ -1,6 +1,8 @@
0.5 (unreleased)
----------------
- Bugfix: subcubes from compound regions previously did not work. #601
- Bugfix: VaryingResolutionSpectralCube.mask_channels now preserves
previous mask. #620
- Refactor tests to use fixtures for accessing data instead of needing to
run a script to generate test files. #598
- Refactor package infrastructure to no longer use astropy-helpers. #599
Expand Down
10 changes: 3 additions & 7 deletions spectral_cube/spectral_cube.py
Expand Up @@ -3971,14 +3971,10 @@ def mask_channels(self, goodchannels):
raise ValueError("goodchannels must have a length equal to the "
"cube's spectral dimension.")

mask = BooleanArrayMask(goodchannels[:,None,None], self._wcs,
shape=self._data.shape)

return self._new_cube_with(mask=mask,
goodbeams_mask=np.logical_and(goodchannels,
self.goodbeams_mask))

cube = self.with_mask(goodchannels[:,None,None])
cube.goodbeams_mask = np.logical_and(goodchannels, self.goodbeams_mask)

return cube

def spectral_interpolate(self, *args, **kwargs):
raise AttributeError("VaryingResolutionSpectralCubes can't be "
Expand Down
22 changes: 22 additions & 0 deletions spectral_cube/tests/test_spectral_cube.py
Expand Up @@ -2224,3 +2224,25 @@ def test_mask_none():
[[0, 1, 2, 3], [4, 5, 6, 7], [8, 9, 10, 11]] * u.Jy / u.beam)
assert_quantity_allclose(cube[:, 0, 0],
[0, 12] * u.Jy / u.beam)


@pytest.mark.parametrize('filename', ['data_vda', 'data_vda_beams'],
indirect=['filename'])
def test_mask_channels_preserve_mask(filename):

# Regression test for a bug that caused the mask to not be preserved.

cube, data = cube_and_raw(filename)

# Add a mask to the cube
mask = np.ones(cube.shape, dtype=bool)
mask[:, ::2, ::2] = False
cube = cube.with_mask(mask)

# Mask by channels
cube = cube.mask_channels([False, True, False, True])

# Check final mask is a combination of both
expected_mask = mask.copy()
expected_mask[::2] = False
np.testing.assert_equal(cube.mask.include(), expected_mask)

0 comments on commit 4c23f5b

Please sign in to comment.