/
polygon.py
205 lines (162 loc) · 6.49 KB
/
polygon.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
# 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 astropy.wcs.utils import skycoord_to_pixel, pixel_to_skycoord
from ..core import PixelRegion, SkyRegion, RegionMask, BoundingBox, PixCoord
from .._geometry import polygonal_overlap_grid
from .._geometry.pnpoly import points_in_polygon
from ..core.attributes import OneDPix, OneDSky, RegionMeta, RegionVisual
__all__ = ['PolygonPixelRegion', 'PolygonSkyRegion']
class PolygonPixelRegion(PixelRegion):
"""
A polygon in pixel coordinates.
Parameters
----------
vertices : `~regions.PixCoord`
The vertices of the polygon
meta : `~regions.RegionMeta` object, optional
A dictionary which stores the meta attributes of this region.
visual : `~regions.RegionVisual` object, optional
A dictionary which stores the visual meta attributes of this region.
Examples
--------
.. plot::
:include-source:
import numpy as np
from astropy.coordinates import Angle
from regions import PixCoord, PolygonPixelRegion
import matplotlib.pyplot as plt
x, y = [45, 45, 55, 60], [75, 70, 65, 75]
fig, ax = plt.subplots(1, 1)
vertices = PixCoord(x=x, y=y)
reg = PolygonPixelRegion(vertices=vertices)
patch = reg.as_artist(facecolor='none', edgecolor='red', lw=2)
ax.add_patch(patch)
plt.xlim(30, 80)
plt.ylim(50, 80)
ax.set_aspect('equal')
"""
_params = ('vertices',)
vertices = OneDPix('vertices')
def __init__(self, vertices, meta=None, visual=None):
self.vertices = vertices
self.meta = meta or RegionMeta()
self.visual = visual or RegionVisual()
@property
def area(self):
"""Region area (float)."""
# See https://stackoverflow.com/questions/24467972
# Use offsets to improve numerical precision
x_ = self.vertices.x - self.vertices.x.mean()
y_ = self.vertices.y - self.vertices.y.mean()
# Shoelace formula, for our case where the start vertex
# isn't duplicated at the end, written to avoid an array copy
area_main = np.dot(x_[:-1], y_[1:]) - np.dot(y_[:-1], x_[1:])
area_last = x_[-1] * y_[0] - y_[-1] * x_[0]
return 0.5 * np.abs(area_main + area_last)
def contains(self, pixcoord):
pixcoord = PixCoord._validate(pixcoord, 'pixcoord')
x = np.atleast_1d(np.asarray(pixcoord.x, dtype=float))
y = np.atleast_1d(np.asarray(pixcoord.y, dtype=float))
vx = np.asarray(self.vertices.x, dtype=float)
vy = np.asarray(self.vertices.y, dtype=float)
shape = x.shape
mask = points_in_polygon(x.flatten(), y.flatten(), vx, vy).astype(bool)
in_poly = mask.reshape(shape)
if self.meta.get('include', True):
return in_poly
else:
return np.logical_not(in_poly)
def to_sky(self, wcs):
vertices_sky = pixel_to_skycoord(self.vertices.x, self.vertices.y, wcs)
return PolygonSkyRegion(vertices=vertices_sky)
@property
def bounding_box(self):
xmin = self.vertices.x.min()
xmax = self.vertices.x.max()
ymin = self.vertices.y.min()
ymax = self.vertices.y.max()
return BoundingBox.from_float(xmin, xmax, ymin, ymax)
def to_mask(self, mode='center', subpixels=5):
self._validate_mode(mode, subpixels)
if mode == 'center':
mode = 'subpixels'
subpixels = 1
if mode == 'subpixels':
use_exact = 0
else:
use_exact = 1
# Find bounding box and mask size
bbox = self.bounding_box
ny, nx = bbox.shape
# Find position of pixel edges and recenter so that circle is at origin
xmin = float(bbox.ixmin) - 0.5
xmax = float(bbox.ixmax) - 0.5
ymin = float(bbox.iymin) - 0.5
ymax = float(bbox.iymax) - 0.5
vx = np.asarray(self.vertices.x, dtype=float)
vy = np.asarray(self.vertices.y, dtype=float)
fraction = polygonal_overlap_grid(
xmin, xmax, ymin, ymax,
nx, ny, vx, vy, use_exact, subpixels,
)
return RegionMask(fraction, bbox=bbox)
def as_artist(self, origin=(0, 0), **kwargs):
"""
Matplotlib patch object for this region (`matplotlib.patches.Polygon`).
Parameters
----------
origin : array_like, optional
The ``(x, y)`` pixel position of the origin of the displayed image.
Default is (0, 0).
kwargs : `dict`
All keywords that a `~matplotlib.patches.Polygon` object accepts
Returns
-------
patch : `~matplotlib.patches.Polygon`
Matplotlib polygon patch
"""
from matplotlib.patches import Polygon
xy = np.vstack([self.vertices.x - origin[0],
self.vertices.y - origin[1]]).transpose()
mpl_params = self.mpl_properties_default('patch')
mpl_params.update(kwargs)
return Polygon(xy=xy, **mpl_params)
def rotate(self, center, angle):
"""Make a rotated region.
Rotates counter-clockwise for positive ``angle``.
Parameters
----------
center : `PixCoord`
Rotation center point
angle : `~astropy.coordinates.Angle`
Rotation angle
Returns
-------
region : `PolygonPixelRegion`
Rotated region (an independent copy)
"""
vertices = self.vertices.rotate(center, angle)
return self.copy(vertices=vertices)
class PolygonSkyRegion(SkyRegion):
"""
A polygon defined using vertices in sky coordinates.
Parameters
----------
vertices : `~astropy.coordinates.SkyCoord`
The vertices of the polygon
meta : `~regions.RegionMeta` object, optional
A dictionary which stores the meta attributes of this region.
visual : `~regions.RegionVisual` object, optional
A dictionary which stores the visual meta attributes of this region.
"""
_params = ('vertices',)
vertices = OneDSky('vertices')
def __init__(self, vertices, meta=None, visual=None):
self.vertices = vertices
self.meta = meta or RegionMeta()
self.visual = visual or RegionVisual()
def to_pixel(self, wcs):
x, y = skycoord_to_pixel(self.vertices, wcs)
vertices_pix = PixCoord(x, y)
return PolygonPixelRegion(vertices_pix)