Skip to content

Commit

Permalink
Split out conversion from CircleSkyRegion to PolygonSkyRegion
Browse files Browse the repository at this point in the history
  • Loading branch information
astrofrog committed Dec 1, 2016
1 parent 888990c commit 4b9bfb1
Show file tree
Hide file tree
Showing 6 changed files with 64 additions and 30 deletions.
11 changes: 11 additions & 0 deletions regions/core/core.py
Original file line number Diff line number Diff line change
Expand Up @@ -284,6 +284,17 @@ def area(self):
Returns the area of the region as a `~astropy.units.Quantity`.
"""

@abc.abstractmethod
def to_polygon(self, points=100):
"""
Convert the region to a polygon.
Parameters
----------
points : int, optional
Number of points in the final polygon.
"""

@abc.abstractmethod
def to_pixel(self, wcs, mode='local', tolerance=None):
"""
Expand Down
53 changes: 26 additions & 27 deletions regions/shapes/circle.py
Original file line number Diff line number Diff line change
Expand Up @@ -7,11 +7,11 @@

from astropy import units as u
from astropy.coordinates import Angle, SkyCoord
from astropy.wcs.utils import pixel_to_skycoord, skycoord_to_pixel
from astropy.wcs.utils import pixel_to_skycoord

from .polygon import PolygonPixelRegion
from .polygon import PolygonSkyRegion

from ..core import PixelRegion, SkyRegion, Mask, BoundingBox, PixCoord
from ..core import PixelRegion, SkyRegion, Mask, BoundingBox
from .._utils.wcs_helpers import skycoord_to_pixel_scale_angle
from .._geometry import circular_overlap_grid, rotate_polygon

Expand Down Expand Up @@ -158,6 +158,28 @@ def area(self):
def contains(self, skycoord):
return self.center.separation(skycoord) < self.radius

def to_polygon(self, points=100):

# TODO: avoid converting to unit spherical or spherical if already
# using a spherical representation

# Extract longitude/latitude, either from a tuple of two quantities, or
# a single 2-element Quantity.
rep = self.center.represent_as('unitspherical')
longitude, latitude = rep.lon, rep.lat

# Start off by generating the circle around the North pole
lon = np.linspace(0., 2 * np.pi, points + 1)[:-1] * u.radian
lat = np.repeat(0.5 * np.pi - self.radius.to(u.radian).value, points) * u.radian

# Now rotate it to the correct longitude/latitude
lon, lat = rotate_polygon(lon, lat, longitude, latitude)

# Make a new SkyCoord
vertices_sky = SkyCoord(lon, lat, frame=self.center)

return PolygonSkyRegion(vertices_sky)

def to_pixel(self, wcs, mode='local', tolerance=100):
"""
Given a WCS, convert the circle to a best-approximation circle in pixel
Expand Down Expand Up @@ -187,30 +209,7 @@ def to_pixel(self, wcs, mode='local', tolerance=100):

elif mode == 'full':

# TODO: avoid converting to unit spherical or spherical if already
# using a spherical representation

# Extract longitude/latitude, either from a tuple of two quantities, or
# a single 2-element Quantity.
rep = self.center.represent_as('unitspherical')
longitude, latitude = rep.lon, rep.lat

# Start off by generating the circle around the North pole
lon = np.linspace(0., 2 * np.pi, tolerance + 1)[:-1] * u.radian
lat = np.repeat(0.5 * np.pi - self.radius.to(u.radian).value, tolerance) * u.radian

# Now rotate it to the correct longitude/latitude
lon, lat = rotate_polygon(lon, lat, longitude, latitude)

# Make a new SkyCoord
vertices_sky = SkyCoord(lon, lat, frame=self.center)

# Convert to PixCoord
x, y = skycoord_to_pixel(vertices_sky, wcs)
vertices_pix = PixCoord(x, y)

# Make polygon
return PolygonPixelRegion(vertices_pix)
return self.to_polygon(points=tolerance).to_pixel(wcs)

else:
raise ValueError('mode should be one of local/affine/full')
4 changes: 4 additions & 0 deletions regions/shapes/ellipse.py
Original file line number Diff line number Diff line change
Expand Up @@ -196,6 +196,10 @@ def contains(self, skycoord):
# TODO: needs to be implemented
raise NotImplementedError

def to_polygon(self, points=100):
# TODO: needs to be implemented
raise NotImplementedError

def to_pixel(self, wcs, mode='local', tolerance=None):
# TODO: needs to be implemented
raise NotImplementedError
4 changes: 4 additions & 0 deletions regions/shapes/point.py
Original file line number Diff line number Diff line change
Expand Up @@ -88,6 +88,10 @@ def area(self):
def contains(self, skycoord):
return False

def to_polygon(self, points=100):
# TODO: needs to be implemented
raise NotImplementedError

def to_pixel(self, wcs, mode='local', tolerance=None):
# TODO: needs to be implemented
raise NotImplementedError
18 changes: 15 additions & 3 deletions regions/shapes/polygon.py
Original file line number Diff line number Diff line change
Expand Up @@ -4,7 +4,9 @@

import numpy as np

from ..core import PixelRegion, SkyRegion, Mask, BoundingBox
from astropy.wcs.utils import skycoord_to_pixel

from ..core import PixelRegion, SkyRegion, Mask, BoundingBox, PixCoord
from .._geometry import polygonal_overlap_grid
from .._geometry.pnpoly import points_in_polygon

Expand Down Expand Up @@ -140,6 +142,16 @@ def contains(self, skycoord):
# TODO: needs to be implemented
raise NotImplementedError

def to_pixel(self, wcs, mode='local', tolerance=None):
# TODO: needs to be implemented
def to_polygon(self, points=100):
# TODO: decide how this should behave. Return self?
raise NotImplementedError

def to_pixel(self, wcs, mode='local', tolerance=None):

if mode == 'local':
x, y = skycoord_to_pixel(self.vertices, wcs)
vertices_pix = PixCoord(x, y)
return PolygonPixelRegion(vertices_pix)
else:

raise NotImplementedError()
4 changes: 4 additions & 0 deletions regions/shapes/rectangle.py
Original file line number Diff line number Diff line change
Expand Up @@ -215,6 +215,10 @@ def contains(self, skycoord):
# TODO: needs to be implemented
raise NotImplementedError

def to_polygon(self, points=100):
# TODO: needs to be implemented
raise NotImplementedError

def to_pixel(self, wcs, mode='local', tolerance=None):
# TODO: needs to be implemented
raise NotImplementedError

0 comments on commit 4b9bfb1

Please sign in to comment.