Skip to content
Draft
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
15 changes: 15 additions & 0 deletions docs/shapes.rst
Original file line number Diff line number Diff line change
Expand Up @@ -59,6 +59,21 @@ Circle
... outer_radius=5.2)


`~regions.CircleSectorPixelRegion`

.. code-block:: python

>>> from astropy.coordinates import SkyCoord
>>> from astropy import units as u
>>> from regions import PixCoord
>>> from regions import CircleSectorPixelRegion

>>> region_pix = CircleSectorPixelRegion(center=PixCoord(x=42, y=43),
... radius=4.2,
... angle_start=0 * u.deg,
... angle_stop=120 * u.deg)


Ellipse
*******

Expand Down
179 changes: 178 additions & 1 deletion regions/shapes/circle.py
Original file line number Diff line number Diff line change
Expand Up @@ -20,7 +20,8 @@
from regions.core.metadata import RegionMeta, RegionVisual
from regions.core.pixcoord import PixCoord

__all__ = ['CirclePixelRegion', 'CircleSkyRegion']

__all__ = ['CirclePixelRegion', 'CircleSkyRegion', 'CircleSectorPixelRegion']


class CirclePixelRegion(PixelRegion):
Expand Down Expand Up @@ -215,3 +216,179 @@ def to_pixel(self, wcs):
radius = (self.radius / pixscale).to(u.pix).value
return CirclePixelRegion(center, radius, meta=self.meta.copy(),
visual=self.visual.copy())


class CircleSectorPixelRegion(PixelRegion):
"""
A circle sector defined using pixel coordinates.

Parameters
----------
center : `~regions.PixCoord`
The center position.
radius : float
The radius in pixels.
angle_start: `~astropy.units.Quantity`, optional
The start angle of the sector, measured anti-clockwise.
angle_stop : `~astropy.units.Quantity`, optional
The stop angle of the sector, measured anti-clockwise.
meta : `~regions.RegionMeta` or `dict`, optional
A dictionary that stores the meta attributes of the region.
visual : `~regions.RegionVisual` or `dict`, optional
A dictionary that stores the visual meta attributes of the
region.

Examples
--------
.. plot::
:include-source:

from astropy import units as u
from regions import PixCoord, CircleSectorPixelRegion
import matplotlib.pyplot as plt

fig, ax = plt.subplots(1, 1)

reg = CircleSectorPixelRegion(PixCoord(x=8, y=7), radius=3.5, angle_start=0 * u.deg,
angle_stop=120 * u.deg)
patch = reg.plot(ax=ax, facecolor='none', edgecolor='red', lw=2,
label='Circle')

ax.legend(handles=(patch,), loc='upper center')
ax.set_xlim(0, 15)
ax.set_ylim(0, 15)
ax.set_aspect('equal')
"""

_params = ('center', 'radius', 'angle_start', 'angle_stop')
_mpl_artist = 'Patch'
center = ScalarPixCoord('The center pixel position as a |PixCoord|.')
radius = PositiveScalar('The radius in pixels as a float.')
angle_start = ScalarAngle('The start angle measured anti-clockwise as a '
'|Quantity| angle.')
angle_stop = ScalarAngle('The stop angle measured anti-clockwise as a '
'|Quantity| angle.')
meta = RegionMetaDescr('The meta attributes as a |RegionMeta|')
visual = RegionVisualDescr('The visual attributes as a |RegionVisual|.')

def __init__(self, center, radius, angle_start=0 * u.deg, angle_stop=360 * u.deg,
meta=None, visual=None):
self.center = center
self.radius = radius

if angle_start >= angle_stop:
raise ValueError('angle_stop must be greater than angle_start')

self.angle_start = angle_start
self.angle_stop = angle_stop
self.meta = meta or RegionMeta()
self.visual = visual or RegionVisual()

@property
def theta(self):
"""Opening angle of the sector (`~astropy.coordinates.Angle`)"""
return self.angle_stop - self.angle_start

@property
def area(self):
return self.radius ** 2 * self.theta.to_value("rad") / 2.

def contains(self, pixcoord):
pixcoord = PixCoord._validate(pixcoord, name='pixcoord')
in_circle = self.center.separation(pixcoord) < self.radius

dx = pixcoord.x - self.center.x
dy = pixcoord.y - self.center.y
angle = (Angle(np.arctan2(dy, dx), "rad") - self.angle_start).wrap_at("360d")

in_angle = (angle > 0 * u.deg) & (angle < self.theta)
in_sector = in_circle & in_angle

if self.meta.get('include', True):
return in_sector
else:
return np.logical_not(in_sector)

def to_sky(self, wcs):
raise NotImplementedError

def to_mask(self, **kwargs):
raise NotImplementedError

@property
def bounding_box(self):
"""Bounding box (`~regions.RegionBoundingBox`)."""
x_start = self.radius * np.cos(self.angle_start)
y_start = self.radius * np.sin(self.angle_start)

x_stop = self.radius * np.cos(self.angle_stop)
y_stop = self.radius * np.sin(self.angle_stop)

def wrap(angle):
return Angle(angle).wrap_at("360d")

cross_0 = wrap(self.angle_start) > wrap(self.angle_stop)
cross_90 = wrap(self.angle_start - 90 * u.deg) > wrap(self.angle_stop - 90 * u.deg)
cross_180 = wrap(self.angle_start - 180 * u.deg) > wrap(self.angle_stop - 180 * u.deg)
cross_270 = wrap(self.angle_start - 270 * u.deg) > wrap(self.angle_stop - 270 * u.deg)

xmin = self.center.x + min(np.where(cross_180, -self.radius, min(x_start, x_stop)), 0)
xmax = self.center.x + max(np.where(cross_0, self.radius, max(x_start, x_stop)), 0)
ymin = self.center.y + min(np.where(cross_270, -self.radius, min(y_start, y_stop)), 0)
ymax = self.center.y + max(np.where(cross_90, self.radius, max(y_start, y_stop)), 0)

return RegionBoundingBox.from_float(xmin, xmax, ymin, ymax)

def as_artist(self, origin=(0, 0), **kwargs):
"""
Return a matplotlib patch object for this region
(`matplotlib.patches.Wedge).

Parameters
----------
origin : array_like, optional
The ``(x, y)`` pixel position of the origin of the displayed
image.

**kwargs : dict
Any keyword arguments accepted by
`~matplotlib.patches.Circle`. These keywords will override
any visual meta attributes of this region.

Returns
-------
artist : `~matplotlib.patches.Wedge`
A matplotlib circle patch.
"""
from matplotlib.patches import Wedge

center = self.center.x - origin[0], self.center.y - origin[1]
radius = self.radius
mpl_kwargs = self.visual.define_mpl_kwargs(self._mpl_artist)
mpl_kwargs.update(kwargs)

return Wedge(center=center, r=radius, theta1=self.angle_start.to_value("deg"),
theta2=self.angle_stop.to_value("deg"), **mpl_kwargs)

def rotate(self, center, angle):
"""
Rotate the region.

Positive ``angle`` corresponds to counter-clockwise rotation.

Parameters
----------
center : `~regions.PixCoord`
The rotation center point.
angle : `~astropy.coordinates.Angle`
The rotation angle.

Returns
-------
region : `CirclePixelRegion`
The rotated region (which is an independent copy).
"""
center = self.center.rotate(center, angle)
angle_start = self.angle_start + angle
angle_stop = self.angle_stop + angle
return self.copy(center=center, angle_start=angle_start, angle_stop=angle_stop)
27 changes: 27 additions & 0 deletions regions/shapes/tests/test_circle.py
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@
from regions.tests.helpers import make_simple_wcs



@pytest.fixture(scope='session', name='wcs')
def wcs_fixture():
filename = get_pkg_data_filename('data/example_header.fits')
Expand Down Expand Up @@ -140,3 +141,29 @@ def test_eq(self):
def test_zero_size(self):
with pytest.raises(ValueError):
CircleSkyRegion(SkyCoord(3 * u.deg, 4 * u.deg), 0. * u.arcsec)


class TestCircleSectorPixelRegion(BaseTestPixelRegion):
meta = RegionMeta({'text': 'test'})
visual = RegionVisual({'color': 'blue'})
reg = CircleSectorPixelRegion(center=PixCoord(3, 4), radius=2, angle_start=30 * u.deg,
angle_stop=120 * u.deg, meta=meta, visual=visual)
sample_box = [0, 6, 1, 7]
inside = [(3, 5)]
outside = [(2, 4)]
expected_area = np.pi
expected_repr = ('<CircleSectorPixelRegion(center=PixCoord(x=3, y=4), radius=2, '
'angle_start=30.0 deg, angle_stop=120.0 deg)>')
expected_str = ('Region: CircleSectorPixelRegion\n'
'center: PixCoord(x=3, y=4)\n'
'radius: 2\n'
'angle_start: 30.0 deg\n'
'angle_stop: 120.0 deg')

@pytest.mark.skipif('not HAS_MATPLOTLIB')
def test_as_artist(self):
patch = self.reg.as_artist()
assert_allclose(patch.center, (3, 4))
assert_allclose(patch.r, 2)
assert_allclose(patch.theta1, 30)
assert_allclose(patch.theta2, 120)
Loading