Skip to content

Commit

Permalink
Merge pull request #195 from sushobhana/annulus
Browse files Browse the repository at this point in the history
Improvements in masking for compound regions.
  • Loading branch information
keflavich committed Jul 7, 2018
2 parents dba1b4a + 780d8fd commit d9cd14c
Show file tree
Hide file tree
Showing 2 changed files with 38 additions and 15 deletions.
15 changes: 5 additions & 10 deletions regions/core/compound.py
Original file line number Diff line number Diff line change
Expand Up @@ -43,6 +43,7 @@ def contains(self, pixcoord):
return in_reg

def to_mask(self, mode='center', subpixels=1):

if mode != 'center':
raise NotImplementedError

Expand All @@ -57,19 +58,13 @@ def to_mask(self, mode='center', subpixels=1):
iymax=max(mask1.bbox.iymax, mask2.bbox.iymax)
)

# Bounding boxes must not extend over array, see #168
bbox_borders = np.array([bbox.ixmin, bbox.ixmax, bbox.iymin, bbox.iymax])
if (bbox_borders < 0).any():
raise NotImplementedError("Bounding box must be within array for "
"compound regions, see #168")

# Pad mask1.data and mask2.data to get the same shape
padded_data = list()
for mask in (mask1, mask2):
pleft = mask.bbox.ixmin - bbox.ixmin
pright = bbox.ixmax - mask.bbox.ixmax
ptop = bbox.iymax - mask.bbox.iymax
pbottom = mask.bbox.iymin - bbox.iymin
pleft = abs(mask.bbox.ixmin - bbox.ixmin)
pright = abs(bbox.ixmax - mask.bbox.ixmax)
ptop = abs(bbox.iymax - mask.bbox.iymax)
pbottom = abs(mask.bbox.iymin - bbox.iymin)
padded_data.append(np.pad(mask.data,
((ptop, pbottom), (pleft, pright)),
'constant'))
Expand Down
38 changes: 33 additions & 5 deletions regions/core/tests/test_compound.py
Original file line number Diff line number Diff line change
@@ -1,12 +1,16 @@
# Licensed under a 3-clause BSD style license - see LICENSE.rst
from __future__ import absolute_import, division, print_function, unicode_literals

import numpy as np
from numpy.testing import assert_allclose

import astropy.units as u
from astropy.coordinates import SkyCoord

from ...shapes import CircleSkyRegion, CirclePixelRegion
from ...core import PixCoord, CompoundPixelRegion
from ...tests.helpers import make_simple_wcs
import pytest
from numpy.testing import assert_allclose


def test_compound_pixel():
# Two circles that overlap in one column
Expand All @@ -31,12 +35,36 @@ def test_compound_pixel():
assert_allclose(mask.data[:,7], [0, 0, 0, 0, 0, 0, 0, 0, 0])
assert_allclose(mask.data[:,6], [0, 1, 1, 1, 1, 1, 1, 1, 0])

# Circle bigger than map, see #168
# fixes #168
pixcoord3 = PixCoord(1, 1)
c3 = CirclePixelRegion(pixcoord3, 4)
union2 = c1 | c3
with pytest.raises(NotImplementedError):
mask = union2.to_mask()
mask1 = union2.to_mask()

pixcoord4 = PixCoord(9, 9)
c4 = CirclePixelRegion(pixcoord4, 4)
union3 = c1 | c4
mask2 = union3.to_mask()

# mask1 and mask2 should be equal
assert_allclose(mask1.data, mask2.data)

ref_data = np.array([[0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0],
[0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0],
[0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 1, 1, 0],
[0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 1, 1, 0],
[0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 1, 1, 0],
[0, 0, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 0],
[0, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 0],
[0, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 0, 0],
[0, 1, 1, 1, 1, 1, 1, 1, 0, 0, 0, 0, 0],
[0, 1, 1, 1, 1, 1, 1, 1, 0, 0, 0, 0, 0],
[0, 1, 1, 1, 1, 1, 1, 1, 0, 0, 0, 0, 0],
[0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 0, 0],
[0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0]]
)
assert_allclose(mask1.data, ref_data)


def test_compound_sky():

Expand Down

0 comments on commit d9cd14c

Please sign in to comment.