World Coordinate System (astropy.wcs
)¶
Introduction¶
astropy.wcs
contains utilities for managing World Coordinate System
(WCS) transformations in FITS files. These transformations map the
pixel locations in an image to their real-world units, such as their
position on the sky sphere. These transformations can work both
forward (from pixel to sky) and backward (from sky to pixel).
It performs three separate classes of WCS transformations:
- Core WCS, as defined in the FITS WCS standard, based on Mark
Calabretta’s wcslib. (Also includes
TPV
andTPD
distortion, but notSIP
). - Simple Imaging Polynomial (SIP) convention.
- table lookup distortions as defined in the FITS WCS distortion paper.
Each of these transformations can be used independently or together in a standard pipeline.
Getting Started¶
The basic workflow is as follows:
from astropy import wcs
- Call the
WCS
constructor with anastropy.io.fits
Header
and/orHDUList
object.- Optionally, if the FITS file uses any deprecated or non-standard features, you may need to call one of the
fix
methods on the object.- Use one of the following transformation methods:
- From pixels to world coordinates:
all_pix2world
: Perform all three transformations in series (core WCS, SIP and table lookup distortions) from pixel to world coordinates. Use this one if you’re not sure which to use.wcs_pix2world
: Perform just the core WCS transformation from pixel to world coordinates.- From world to pixel coordinates:
all_world2pix
: Perform all three transformations (core WCS, SIP and table lookup distortions) from world to pixel coordinates, using an iterative method if necessary.wcs_world2pix
: Perform just the core WCS transformation from world to pixel coordinates.- Performing SIP transformations only:
sip_pix2foc
: Convert from pixel to focal plane coordinates using the SIP polynomial coefficients.sip_foc2pix
: Convert from focal plane to pixel coordinates using the SIP polynomial coefficients.- Performing distortion paper transformations only:
p4_pix2foc
: Convert from pixel to focal plane coordinates using the table lookup distortion method described in the FITS WCS distortion paper.det2im
: Convert from detector coordinates to image coordinates. Commonly used for narrow column correction.
For example, to convert pixel coordinates to world coordinates:
>>> from astropy.wcs import WCS
>>> w = WCS('image.fits')
>>> lon, lat = w.all_pix2world(30, 40, 0)
>>> print(lon, lat)
Using astropy.wcs
¶
Loading WCS information from a FITS file¶
This example loads a FITS file (supplied on the commandline) and uses the WCS cards in its primary header to transform.
# Load the WCS information from a fits header, and use it
# to convert pixel coordinates to world coordinates.
from __future__ import division, print_function
import numpy
from astropy import wcs
from astropy.io import fits
import sys
def load_wcs_from_file(filename):
# Load the FITS hdulist using astropy.io.fits
hdulist = fits.open(filename)
# Parse the WCS keywords in the primary HDU
w = wcs.WCS(hdulist[0].header)
# Print out the "name" of the WCS, as defined in the FITS header
print(w.wcs.name)
# Print out all of the settings that were parsed from the header
w.wcs.print_contents()
# Some pixel coordinates of interest.
pixcrd = numpy.array([[0, 0], [24, 38], [45, 98]], numpy.float_)
# Convert pixel coordinates to world coordinates
# The second argument is "origin" -- in this case we're declaring we
# have 1-based (Fortran-like) coordinates.
world = w.wcs_pix2world(pixcrd, 1)
print(world)
# Convert the same coordinates back to pixel coordinates.
pixcrd2 = w.wcs_world2pix(world, 1)
print(pixcrd2)
# These should be the same as the original pixel coordinates, modulo
# some floating-point error.
assert numpy.max(numpy.abs(pixcrd - pixcrd2)) < 1e-6
if __name__ == '__main__':
load_wcs_from_file(sys.argv[-1])
Building a WCS structure programmatically¶
This example, rather than starting from a FITS header, sets WCS values programmatically, uses those settings to transform some points, and then saves those settings to a new FITS header.
# Set the WCS information manually by setting properties of the WCS
# object.
from __future__ import division, print_function
import numpy
from astropy import wcs
from astropy.io import fits
# Create a new WCS object. The number of axes must be set
# from the start
w = wcs.WCS(naxis=2)
# Set up an "Airy's zenithal" projection
# Vector properties may be set with Python lists, or Numpy arrays
w.wcs.crpix = [-234.75, 8.3393]
w.wcs.cdelt = numpy.array([-0.066667, 0.066667])
w.wcs.crval = [0, -90]
w.wcs.ctype = ["RA---AIR", "DEC--AIR"]
w.wcs.set_pv([(2, 1, 45.0)])
# Some pixel coordinates of interest.
pixcrd = numpy.array([[0, 0], [24, 38], [45, 98]], numpy.float_)
# Convert pixel coordinates to world coordinates
world = w.wcs_pix2world(pixcrd, 1)
print(world)
# Convert the same coordinates back to pixel coordinates.
pixcrd2 = w.wcs_world2pix(world, 1)
print(pixcrd2)
# These should be the same as the original pixel coordinates, modulo
# some floating-point error.
assert numpy.max(numpy.abs(pixcrd - pixcrd2)) < 1e-6
# Now, write out the WCS object as a FITS header
header = w.to_header()
# header is an astropy.io.fits.Header object. We can use it to create a new
# PrimaryHDU and write it to a file.
hdu = fits.PrimaryHDU(header=header)
# Save to FITS file
# hdu.writeto('test.fits')
Note
The members of the WCS object correspond roughly to the key/value
pairs in the FITS header. However, they are adjusted and
normalized in a number of ways that make performing the WCS
transformation easier. Therefore, they can not be relied upon to
get the original values in the header. To build up a FITS header
directly and specifically, use astropy.io.fits.Header
directly.
Validating the WCS keywords in a FITS file¶
Astropy includes a commandline tool, wcslint
to check the WCS
keywords in a FITS file:
> wcslint invalid.fits
HDU 1:
WCS key ' ':
- RADECSYS= 'ICRS ' / Astrometric system
RADECSYS is non-standard, use RADESYSa.
- The WCS transformation has more axes (2) than the image it is
associated with (0)
- 'celfix' made the change 'PV1_5 : Unrecognized coordinate
transformation parameter'.
HDU 2:
WCS key ' ':
- The WCS transformation has more axes (3) than the image it is
associated with (0)
- 'celfix' made the change 'In CUNIT2 : Mismatched units type
'length': have 'Hz', want 'm''.
- 'unitfix' made the change 'Changed units: 'HZ ' -> 'Hz''.
Bounds checking¶
Bounds checking is enabled by default, and any computed world
coordinates outside of [-180°, 180°] for longitude and [-90°, 90°] in
latitude are marked as invalid. To disable this behavior, use
astropy.wcs.Wcsprm.bounds_check
.
Supported projections¶
As astropy.wcs
is based on wcslib, it supports the standard
projections defined in the FITS WCS standard. These projection
codes are specified in the second part of the CTYPEn
keywords
(accessible through Wcsprm.ctype
), for
example, RA---TAN-SIP
. The supported projection codes are:
AZP
: zenithal/azimuthal perspectiveSZP
: slant zenithal perspectiveTAN
: gnomonicSTG
: stereographicSIN
: orthographic/synthesisARC
: zenithal/azimuthal equidistantZPN
: zenithal/azimuthal polynomialZEA
: zenithal/azimuthal equal areaAIR
: Airy’s projectionCYP
: cylindrical perspectiveCEA
: cylindrical equal areaCAR
: plate carréeMER
: Mercator’s projectionCOP
: conic perspectiveCOE
: conic equal areaCOD
: conic equidistantCOO
: conic orthomorphicSFL
: Sanson-Flamsteed (“global sinusoid”)PAR
: parabolicMOL
: Mollweide’s projectionAIT
: Hammer-AitoffBON
: Bonne’s projectionPCO
: polyconicTSC
: tangential spherical cubeCSC
: COBE quadrilateralized spherical cubeQSC
: quadrilateralized spherical cubeHPX
: HEALPixXPH
: HEALPix polar, aka “butterfly”
And, if built with wcslib 5.0 or later, the following polynomial distortions are supported:
TPV
: Polynomial distortionTUV
: Polynomial distortion
Note
Though wcslib 5.4 and later handles SIP
polynomial distortion,
for backward compatibility, SIP
is handled by astropy itself
and methods exist to handle it specially.
Subsetting and Pixel Scales¶
WCS objects can be broken apart into their constituent axes using the
sub
function. There is also a celestial
convenience function that will return a WCS object with only the celestial axes
included.
The pixel scales of a celestial image or the pixel dimensions of a non-celestial
image can be extracted with the utility functions
proj_plane_pixel_scales
and
non_celestial_pixel_scales
. Likewise, celestial pixel
area can be extracted with the utility function
proj_plane_pixel_area
.
Matplotlib plots with correct WCS projection¶
The WCSAxes affiliated package adds the
ability to use the WCS
to define projections in
Matplotlib. More information on installing and using WCSAxes can be found here.
from matplotlib import pyplot as plt
from astropy.io import fits
from astropy.wcs import WCS
from astropy.utils.data import download_file
fits_file = 'http://data.astropy.org/tutorials/FITS-images/HorseHead.fits'
image_file = download_file(fits_file, cache=True )
hdu = fits.open(image_file)[0]
wcs = WCS(hdu.header)
fig = plt.figure()
fig.add_subplot(111, projection=wcs)
plt.imshow(hdu.data, origin='lower', cmap='cubehelix')
plt.xlabel('RA')
plt.ylabel('Dec')
plt.show()
(Source code, png, hires.png, pdf)
Other information¶
Reference/API¶
astropy.wcs Package¶
astropy.wcs
contains utilities for managing World Coordinate System
(WCS) transformations in FITS files. These transformations map the
pixel locations in an image to their real-world units, such as their
position on the sky sphere.
It performs three separate classes of WCS transformations:
- Core WCS, as defined in the FITS WCS standard, based on Mark
Calabretta’s wcslib. See
Wcsprm
. - Simple Imaging Polynomial (SIP) convention. See
Sip
. - table lookup distortions as defined in WCS distortion paper. See
DistortionLookupTable
.
Each of these transformations can be used independently or together in a standard pipeline.
Functions¶
find_all_wcs (header[, relax, keysel, fix, ...]) |
Find all the WCS transformations in the given header. |
get_include () |
Get the path to astropy.wcs’s C header files. |
validate (source) |
Prints a WCS validation report for the given FITS file. |
Classes¶
DistortionLookupTable |
Represents a single lookup table for a distortion paper transformation. |
FITSFixedWarning |
The warning raised when the contents of the FITS header have been modified to be standards compliant. |
InconsistentAxisTypesError |
The WCS header inconsistent or unrecognized coordinate axis type(s). |
InvalidCoordinateError |
One or more of the world coordinates is invalid. |
InvalidSubimageSpecificationError |
The subimage specification is invalid. |
InvalidTabularParametersError |
The given tabular parameters are invalid. |
InvalidTransformError |
The WCS transformation is invalid, or the transformation parameters are invalid. |
NoSolutionError |
No solution can be found in the given interval. |
NoWcsKeywordsFoundError |
No WCS keywords were found in the given header. |
NonseparableSubimageCoordinateSystemError |
Non-separable subimage coordinate system. |
SingularMatrixError |
The linear transformation matrix is singular. |
Sip |
The Sip class performs polynomial distortion correction using the SIP convention in both directions. |
Tabprm |
A class to store the information related to tabular coordinates, i.e., coordinates that are defined via a lookup table. |
WCS ([header, fobj, key, minerr, relax, ...]) |
WCS objects perform standard WCS transformations, and correct for SIP and distortion paper table-lookup transformations, based on the WCS keywords and supplementary data read from a FITS file. |
WCSBase |
Wcs objects amalgamate basic WCS (as provided by wcslib), with SIP and distortion paper operations. |
WcsError |
Base class of all invalid WCS errors. |
Wcsprm |
Wcsprm performs the core WCS transformations. |
Class Inheritance Diagram¶
astropy.wcs.utils Module¶
Functions¶
add_stokes_axis_to_wcs (wcs, add_before_ind) |
Add a new Stokes axis that is uncorrelated with any other axes. |
wcs_to_celestial_frame (wcs) |
For a given WCS, return the coordinate frame that matches the celestial component of the WCS. |
proj_plane_pixel_scales (wcs) |
For a WCS returns pixel scales along each axis of the image pixel at the CRPIX location once it is projected onto the “plane of intermediate world coordinates” as defined in Greisen & Calabretta 2002, A&A, 395, 1061. |
proj_plane_pixel_area (wcs) |
For a celestial WCS (see astropy.wcs.WCS.celestial ) returns pixel area of the image pixel at the CRPIX location once it is projected onto the “plane of intermediate world coordinates” as defined in Greisen & Calabretta 2002, A&A, 395, 1061. |
is_proj_plane_distorted (wcs[, maxerr]) |
For a WCS returns False if square image (detector) pixels stay square when projected onto the “plane of intermediate world coordinates” as defined in Greisen & Calabretta 2002, A&A, 395, 1061. |
non_celestial_pixel_scales (inwcs) |
For a non-celestial WCS, e.g. |
skycoord_to_pixel (coords, wcs[, origin, mode]) |
Convert a set of SkyCoord coordinates into pixels. |
pixel_to_skycoord (xp, yp, wcs[, origin, ...]) |
Convert a set of pixel coordinates into a SkyCoord coordinate. |
Classes¶
custom_frame_mappings ([mappings]) |
Class Inheritance Diagram¶
Acknowledgments and Licenses¶
wcslib is licenced under the GNU Lesser General Public License.