Permalink
Fetching contributors…
Cannot retrieve contributors at this time
210 lines (148 sloc) 7.08 KB

Getting started

pyregion is a python module to parse ds9 region files. It also supports ciao region files.

Please note that the main emphasis of the package is to read in the regions files generated by ds9 itself. It reads most of the region files created by ds9. However, it may fail to read some of the user-created (or created by other programs) region files, even if they can be successfully read by ds9. Ruler, Compass and Projection types are ignored.

ds9 pyregion + matplotlib
_static/region_ds9.jpg _static/region_mpl.png

Read Region Files

pyregion.open takes the region name as an argument and returns a ~pyregion.ShapeList object, which is basically a list of ~pyregion.Shape objects (~pyregion.ShapeList is a sub-class of the Python built-in list class).

import pyregion
region_name = "ds9.reg"
r = pyregion.open(region_name)

You may use pyregion.parse if you have a string that defines a region

region = 'fk5;circle(290.96388,14.019167,843.31194")'
r = pyregion.parse(region)

The shape object is a python representation of each region definition. For example,:

import pyregion

region_string = """
# Region file format: DS9 version 4.1
# Filename: test01.fits
global color=green dashlist=8 3 width=1 font="helvetica 10 normal" select=1 highlite=1 dash=0 fixed=0 edit=1 move=1 delete=1 include=1 source=1
fk5
circle(11:24:24.230,-59:15:02.20,18.5108") # color=cyan background
box(11:24:39.213,-59:16:53.91,42.804",23.616",19.0384) # width=4
"""

r = pyregion.parse(region_string)

And you have:

>>> print r[0]
Shape : circle ( HMS(11:24:24.230),DMS(-59:15:02.20),Ang(18.5108") )
>>> print r[1]
Shape : box ( HMS(11:24:39.213),DMS(-59:16:53.91),Ang(42.804"),Ang(23.616"),Number(19.0384) )

The shape object has the following attributes,

  • name : name of the shape. e.g., circle, box, etc..

    >>> print r[0].name
    circle
    
  • coord_format : coordinate format. e.g., "fk5", "image", "physical", etc...

    >>> print r[0].coord_format
    fk5
    
  • coord_list : list of coordinates in coord_format. The coordinate value for sky coordinates is degree.

    >>> print r[0].coord_list
    [171.10095833333332, -59.250611111111112, 0.0051418888888888886]
    
  • comment : comment string associated with the shape (can be None)

    >>> print r[0].comment
    color=cyan background
    
  • attr : attributes of the shape. This includes global attributes defined by the global command and local attributes defined in the comment. The first item is a list of key-only attributes without associated values (e.g., background..) and the second item is a dictionary of attributes of key-value pairs.

    >>> print r[0].attr[0]
    ['background']
    >>> print r[0].attr[1]
    {'color': 'cyan',
     'dash': '0 ',
     'dashlist': '8 3 ',
     'delete': '1 ',
     'edit': '1 ',
     'fixed': '0 ',
     'font': '"helvetica 10 normal"',
     'highlite': '1 ',
     'include': '1 ',
     'move': '1 ',
     'select': '1 ',
     'source': '1',
     'width': '1 '}
    

    Some attributes like "tag" allow multiple items, but this is not currently supported (the last definition overrides any previous ones).

The pyregion.ShapeList class have a few methods that could be useful. ShapeList.as_imagecoord <pyregion.ShapeList.as_imagecoord> returns a new ~pyregion.ShapeList instance with the coordinates converted to the image coordinate system. It requires an astropy.io.fits.Header instance.

from astropy.io import fits
f = fits.open("t1.fits")
r2 = pyregion.parse(region_string).as_imagecoord(f[0].header)

The return value is a new ~pyregion.ShapeList instance, but the coordinate is converted to image coordinates.

>>> print r2[0].coord_format
image

>>> print r2[0].coord_list
[482.27721401429852, 472.76641383805912, 18.811792596807045]

ShapeList.as_imagecoord <pyregion.ShapeList.as_imagecoord> will use the subset of the header defining a celestial coordinate system, ignoring any velocity or channel components.

Draw Regions with Matplotlib

pyregion can help you draw ds9 regions with matplotlib. ShapeList.get_mpl_patches_texts <pyregion.ShapeList.get_mpl_patches_texts> returns a list of matplotlib.artist.Artist objects

r2 = pyregion.parse(region_string).as_imagecoord(f[0].header)
patch_list, artist_list = r2.get_mpl_patches_texts()

The first item is a list of matplotlib.patches.Patch, and the second one is other kinds of artists (usually Text). It is your responsibility to add these to the axes.

# ax is a mpl Axes object
for p in patch_list:
    ax.add_patch(p)
for t in artist_list:
    ax.add_artist(t)
.. plot:: figures/region_drawing.py

The (optional) argument of the get_mpl_patches_texts method is a callable object that takes the shape object as an argument and returns a dictionary object that will be used as a keyword arguments (e.g., colors and line width) for creating the mpl artists. By default, it uses pyregion.mpl_helper.properties_func_default, which tries to respect the ds9 attributes. However, the colors (and other attributes) of some complex shapes are not correctly handled as shown in above example, and you need to manually adjust the associated attributes of patches.

.. plot:: figures/region_drawing2.py
   :include-source:



Use Regions for Spatial Filtering

pyregion includes some basic spatial filter support.

The ShapeList.get_filter <pyregion.ShapeList.get_filter> method returns the filter from the parsed region.

The filter is meant to be used in the image coordinate, thus you need to convert the region to the image coordinate before calling get_filter.

r2 = pyregion.parse(region_string).as_imagecoord(f[0].header)
myfilter = r2.get_filter()
myfilter.inside1(50, 30)

The returned filter has a mask method that creates a 2d mask. You can create the mask directly from the ShapeList object.

r2 = pyregion.parse(region_string)
mymask = r2.get_mask(hdu=f[0])

It will creates an mask in the shape of the given hdu image (the mask will be created after transforming the region to the image coordinate if necessary).

.. plot:: figures/demo_filter_mask.py
   :include-source:

Note that this will fail if your template image is not a simple 2D image. To work around this you may use the shape optional argument of ShapeList.get_mask <pyregion.ShapeList.get_mask>:

mymask = r2.get_mask(hdu=f[0],shape=(1024,1024))