Position-Velocity Slice Extractor¶
The concept of the pvextractor
package is simple - given a path defined
in sky coordinates, and a spectral cube, extract a slice of the cube along
that path, and along the spectral axis, producing a position-velocity or
position-frequency slice.
The path can be defined programmatically in pixel or world coordinates, but it can also be drawn interactively using a simple GUI. We also provide a DS9 analysis plug-in that allows the path to be drawn in DS9, and the resulting slice shown in a new frame in DS9. Finally, the slicing capability is available in Glue.
Extracting slices programmatically¶
Defining a path¶
Pixel coordinates¶
To define a path in pixel coordinates, import the Path
class:
>>> from pvextractor import Path
Then initialize it using a list of (x,y)
tuples. The simplest path that
can be defined is a line connecting two points:
>>> path1 = Path([(0., 0.), (10., 10.)])
Multi-segment paths can also be similarly defined:
>>> path2 = Path([(0., 0.), (4., 6.), (10., 10.)])
By default, slices are extracted using interpolation along the line, but it
is also possible to define a path with a finite width, and to instead measure
the average flux or surface brightness in finite polygons (rather than
strictly along the line). To give a path a non-zero width, simply use the
width=
argument, which is also in pixels by default:
>>> path3 = Path([(0., 0.), (10., 10.)], width=0.5)
World coordinates¶
To define a path in world coordinates, pass a coordinate array to the Path
object. In addition, the width (if passed) should an Astropy
Quantity
object:
>>> from astropy import units as u
>>> from astropy.coordinates import Galactic
>>> g = Galactic([3.4, 3.6] * u.deg, [0.5, 0.56] * u.deg)
>>> path4 = Path(g, width=1 * u.arcsec)
In additon to the Path
class, we provide a convenience
PathFromCenter
class that can be used for cases where the
center and position angle of the path are known (rather than the end points of
the path). This class is used as follows:
>>> from pvextractor import PathFromCenter
>>> from astropy import units as u
>>> from astropy.coordinates import Galactic
>>> g = Galactic(3 * u.deg, 5 * u.deg)
>>> path5 = PathFromCenter(center=g,
... length=1 * u.arcmin,
... angle=30 * u.deg,
... width=1 * u.arcsec)
The position angle is defined counter-clockwise from North, and the direction of the path is such that for a position angle of zero, the path is defined from South to North.
Extracting a slice¶
Once the path has been defined, you can make use of the
extract_pv_slice()
function to extract the PV slice. The
data to slice can be passed to this function as:
- A 3-d Numpy array
- A
SpectralCube
instance - An HDU object containing a spectral cube
- The name of a FITS file containing a spectral cube
For example:
>>> from pvextractor import extract_pv_slice
>>> slice1 = extract_pv_slice(array, path1)
>>> from spectral_cube import SpectralCube
>>> cube = SpectralCube.read('my_cube.fits')
>>> slice2 = extract_pv_slice(cube, path3)
>>> slice3 = extract_pv_slice('my_cube.fits', path3)
Note
If a path is passed in in world coordinates, and the data are passed
as a plain Numpy array, the WCS information should be passed as a
WCS
object to the wcs=
argument.
Saving the slice¶
The returned slice is an Astropy PrimaryHDU
instance,
which you can write to disk using:
>>> slice1.writeto('my_slice.fits')
Using the built-in graphical user interface¶
The extractor GUI provdes the most direct interface available to the pixel-matched version of the position-velocity extractor. It is simple to initialize:
from pvextractor.gui import PVSlicer
pv = PVSlicer('cube.fits')
pv.show()
A window will pop up showing the data on the left hand side. Click to draw the vertices of the path, then press “enter” to expand the width of the slice, then move your mouse away from the path to define how wide it should be. Once you are happy with the width, click, and the slice will be computed and shown on the right.
Once a path has been defined, you can optionally press “y” to see the polygons along which the data has been collapsed to compute the slice.
The following shows an example of a slice derived from a 13CO cube of L1448:

Slicing in Glue¶
Glue has a position-velocity diagram extraction tool that allows arbitrary paths to be specified:

Slicing in DS9¶
Note
This feature is experimental and does not work with all versions of DS9 at this point.
There is a python command-line script that will be installed into your path
along with pvextractor
. You can invoke it from the command line, but the
preferred approach is to load the
tool into DS9. First, determine the path to ds9_pvextract.ans
;
it is in the scripts
subdirectory of the source code. Then start
up DS9 with the analysis tool loaded
ds9 -analysis load /path/to/pvextractor/scripts/ds9_pvextract.ans &
Then load any cube in DS9. You can draw a line, a vector, or a “segment”; only the first one drawn will have any effect. To extract the PV diagram, press the ‘x’ key or click “PV Extractor” in the analysis menu. Be patient - especially for big cubes, it may take a little while, and there is no progress bar.
API Documentation¶
pvextractor Package¶
This is an Astropy affiliated package.
Functions¶
extract_pv_slice (cube, path[, wcs, spacing, ...]) |
Given a position-position-velocity cube with dimensions (nv, ny, nx), and a path, extract a position-velocity slice. |
paths_from_regfile (regfile) |
Given a ds9 region file, extract pv diagrams for each: |
slice_wcs (wcs, spatial_scale) |
Slice a WCS header for a spectral cube to a Position-Velocity WCS, with |
test ([package, test_path, args, plugins, ...]) |
Run the tests using py.test. |
Classes¶
Path (xy_or_coords[, width]) |
A curved path that may have a non-zero width and is used to extract slices from cubes. |
PathFromCenter (center[, length, angle, ...]) |
A simple path defined by a center, length, and position angle. |