PointPixelRegion

class regions.PointPixelRegion(center, meta=None, visual=None)[source]

Bases: regions.core.core.PixelRegion

A point position in pixel coordinates.

Parameters
centerPixCoord

The position of the point.

metaRegionMeta, optional

A dictionary that stores the meta attributes of this region.

visualRegionVisual, optional

A dictionary that stores the visual meta attributes of this region.

Examples

from regions import PixCoord, PointPixelRegion, RegionVisual
import matplotlib.pyplot as plt

fig, ax = plt.subplots(1, 1)
regs = []
regs.append(PointPixelRegion(PixCoord(2, 2),
            visual=RegionVisual(symbol='D')))
regs.append(PointPixelRegion(PixCoord(2, 3),
            visual=RegionVisual(symbol='*')))
regs.append(PointPixelRegion(PixCoord(3, 3),
            visual=RegionVisual(symbol='^')))
regs.append(PointPixelRegion(PixCoord(3, 2),
            visual=RegionVisual(symbol='*')))
regs.append(PointPixelRegion(PixCoord(2, 4),
            visual=RegionVisual(symbol='x')))
regs.append(PointPixelRegion(PixCoord(4, 2)))

for reg in regs:
    reg.plot(ax=ax)

plt.xlim(0, 6)
plt.ylim(0, 6)
ax.set_aspect('equal')

(Source code, png, hires.png, pdf, svg)

../_images/regions-PointPixelRegion-1.png

Attributes Summary

area

The area of the region.

bounding_box

The minimal bounding box (in integer pixel coordinates) that contains the region.

center

Descriptor class for PixelRegion, which takes a scalar PixCoord object.

mpl_artist

Methods Summary

as_artist([origin])

Return a matplotlib Line2D object for this region (matplotlib.lines.Line2D).

contains(pixcoord)

Check whether a position or positions fall inside the region.

rotate(center, angle)

Rotate the region.

to_mask([mode, subpixels])

Return a mask for the region.

to_sky(wcs)

Return a region defined in sky coordinates.

Attributes Documentation

area
bounding_box
center

Descriptor class for PixelRegion, which takes a scalar PixCoord object.

mpl_artist = 'Line2D'

Methods Documentation

as_artist(origin=(0, 0), **kwargs)[source]

Return a matplotlib Line2D object for this region (matplotlib.lines.Line2D).

Parameters
originarray_like, optional

The (x, y) pixel position of the origin of the displayed image.

**kwargsdict

Any keyword arguments accepted by Line2D. These keywords will override any visual meta attributes of this region.

Returns
artistLine2D

A matplotlib Line2D object.

contains(pixcoord)[source]

Check whether a position or positions fall inside the region.

Parameters
pixcoordPixCoord

The position or positions to check.

rotate(center, angle)[source]

Rotate the region.

Positive angle corresponds to counter-clockwise rotation.

Parameters
centerPixCoord

The rotation center point.

angleAngle

The rotation angle.

Returns
regionPointPixelRegion

The rotated region (which is an independent copy).

to_mask(mode='center', subpixels=5)[source]

Return a mask for the region.

Parameters
mode{‘center’, ‘exact’, ‘subpixels’}, optional

The method used to determine the overlap of the region on the pixel grid. Not all options are available for all region types. Note that the more precise methods are generally slower. The following methods are available:

  • 'center': A pixel is considered to be entirely in or out of the region depending on whether its center is in or out of the region. The returned mask will contain values only of 0 (out) and 1 (in).

  • 'exact' (default): The exact fractional overlap of the region and each pixel is calculated. The returned mask will contain values between 0 and 1.

  • 'subpixel': A pixel is divided into subpixels (see the subpixels keyword), each of which are considered to be entirely in or out of the region depending on whether its center is in or out of the region. If subpixels=1, this method is equivalent to 'center'. The returned mask will contain values between 0 and 1.

subpixelsint, optional

For the 'subpixel' mode, resample pixels by this factor in each dimension. That is, each pixel is divided into subpixels ** 2 subpixels.

Returns
maskRegionMask

A mask for the region.

to_sky(wcs)[source]

Return a region defined in sky coordinates.

Parameters
wcsWCS

The world coordinate system transformation to use to convert from pixels to sky coordinates.

Returns
sky_regionSkyRegion

The sky region.