WCS#

class astropy.wcs.WCS(header=None, fobj=None, key=' ', minerr=0.0, relax=True, naxis=None, keysel=None, colsel=None, fix=True, translate_units='', _do_set=True)[source]#

Bases: FITSWCSAPIMixin, WCSBase

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.

See also: https://docs.astropy.org/en/stable/wcs/

Parameters:
headerHeader, PrimaryHDU, ImageHDU, str, dict-like, or None, optional

If header is not provided or None, the object will be initialized to default values.

fobjHDUList, optional

It is needed when header keywords point to a distortion paper lookup table stored in a different extension.

keystr, optional

The name of a particular WCS transform to use. This may be either ' ' or 'A'-'Z' and corresponds to the "a" part of the CTYPEia cards. key may only be provided if header is also provided.

minerrfloat, optional

The minimum value a distortion correction must have in order to be applied. If the value of CQERRja is smaller than minerr, the corresponding distortion is not applied.

relaxbool or int, optional

Degree of permissiveness:

  • True (default): Admit all recognized informal extensions of the WCS standard.

  • False: Recognize only FITS keywords defined by the published WCS standard.

  • int: a bit field selecting specific extensions to accept. See Header-reading relaxation constants for details.

naxisint or sequence, optional

Extracts specific coordinate axes using sub(). If a header is provided, and naxis is not None, naxis will be passed to sub() in order to select specific axes from the header. See sub() for more details about this parameter.

keyselsequence of str, optional

A sequence of flags used to select the keyword types considered by wcslib. When None, only the standard image header keywords are considered (and the underlying wcspih() C function is called). To use binary table image array or pixel list keywords, keysel must be set.

Each element in the list should be one of the following strings:

  • ‘image’: Image header keywords

  • ‘binary’: Binary table image array keywords

  • ‘pixel’: Pixel list keywords

Keywords such as EQUIna or RFRQna that are common to binary table image arrays and pixel lists (including WCSNna and TWCSna) are selected by both ‘binary’ and ‘pixel’.

colselsequence of int, optional

A sequence of table column numbers used to restrict the WCS transformations considered to only those pertaining to the specified columns. If None, there is no restriction.

fixbool, optional

When True (default), call fix on the resulting object to fix any non-standard uses in the header. FITSFixedWarning Warnings will be emitted if any changes were made.

translate_unitsstr, optional

Specify which potentially unsafe translations of non-standard unit strings to perform. By default, performs none. See WCS.fix for more information about this parameter. Only effective when fix is True.

Raises:
MemoryError

Memory allocation failed.

ValueError

Invalid key.

KeyError

Key not found in FITS header.

ValueError

Lookup table distortion present in the header but fobj was not provided.

Notes

  1. astropy.wcs supports arbitrary n dimensions for the core WCS (the transformations handled by WCSLIB). However, the distortion paper lookup table and SIP distortions must be two dimensional. Therefore, if you try to create a WCS object where the core WCS has a different number of dimensions than 2 and that object also contains a distortion paper lookup table or SIP distortion, a ValueError exception will be raised. To avoid this, consider using the naxis kwarg to select two dimensions from the core WCS.

  2. The number of coordinate axes in the transformation is not determined directly from the NAXIS keyword but instead from the highest of:

    • NAXIS keyword

    • WCSAXESa keyword

    • The highest axis number in any parameterized WCS keyword. The keyvalue, as well as the keyword, must be syntactically valid otherwise it will not be considered.

    If none of these keyword types is present, i.e. if the header only contains auxiliary WCS keywords for a particular coordinate representation, then no coordinate description is constructed for it.

    The number of axes, which is set as the naxis member, may differ for different coordinate representations of the same image.

  3. When the header includes duplicate keywords, in most cases the last encountered is used.

  4. set is called immediately after construction, so any invalid keywords or transformations will be raised by the constructor, not when subsequently calling a transformation method.

Attributes Summary

array_shape

The shape of the data that the WCS applies to as a tuple of length pixel_n_dim in (row, column) order (the convention for arrays in Python).

axis_correlation_matrix

Returns an (world_n_dim, pixel_n_dim) matrix that indicates using booleans whether a given world coordinate depends on a given pixel coordinate.

axis_type_names

World names for each coordinate axis.

celestial

A copy of the current WCS with only the celestial axes included.

cpdis1

DistortionLookupTable

cpdis2

DistortionLookupTable

det2im1

A DistortionLookupTable object for detector to image plane correction in the x-axis.

det2im2

A DistortionLookupTable object for detector to image plane correction in the y-axis.

has_celestial

has_distortion

Returns True if any distortion terms are present.

has_spectral

has_temporal

is_celestial

is_spectral

is_temporal

low_level_wcs

Returns a reference to the underlying low-level WCS object.

pixel_axis_names

An iterable of strings describing the name for each pixel axis.

pixel_bounds

The bounds (in pixel coordinates) inside which the WCS is defined, as a list with pixel_n_dim (min, max) tuples.

pixel_n_dim

The number of axes in the pixel coordinate system.

pixel_scale_matrix

pixel_shape

The shape of the data that the WCS applies to as a tuple of length pixel_n_dim in (x, y) order (where for an image, x is the horizontal coordinate and y is the vertical coordinate).

serialized_classes

Indicates whether Python objects are given in serialized form or as actual Python objects.

sip

Get/set the Sip object for performing SIP distortion correction.

spectral

A copy of the current WCS with only the spectral axes included.

temporal

A copy of the current WCS with only the time axes included.

wcs

A Wcsprm object to perform the basic wcslib WCS transformation.

world_axis_names

An iterable of strings describing the name for each world axis.

world_axis_object_classes

A dictionary giving information on constructing high-level objects for the world coordinates.

world_axis_object_components

A list with world_n_dim elements giving information on constructing high-level objects for the world coordinates.

world_axis_physical_types

An iterable of strings describing the physical type for each world axis.

world_axis_units

An iterable of strings given the units of the world coordinates for each axis.

world_n_dim

The number of axes in the world coordinate system.

Methods Summary

all_pix2world(*args, **kwargs)

Transforms pixel coordinates to world coordinates.

all_world2pix(*arg[, tolerance, maxiter, ...])

Transforms world coordinates to pixel coordinates, using numerical iteration to invert the full forward transformation all_pix2world with complete distortion model.

array_index_to_world(*index_arrays)

Convert array indices to world coordinates (represented by Astropy objects).

array_index_to_world_values(*index_arrays)

Convert array indices to world coordinates.

calc_footprint([header, undistort, axes, center])

Calculates the footprint of the image on the sky.

copy()

Return a shallow copy of the object.

deepcopy()

Return a deep copy of the object.

det2im(*args)

Convert detector coordinates to image plane coordinates using distortion paper table-lookup correction.

dropaxis(dropax)

Remove an axis from the WCS.

fix([translate_units, naxis])

Perform the fix operations from wcslib, and warn about any changes it has made.

footprint_contains(coord, **kwargs)

Determines if a given SkyCoord is contained in the wcs footprint.

footprint_to_file([filename, color, width, ...])

Writes out a ds9 style regions file.

get_axis_types()

Similar to self.wcsprm.axis_types but provides the information in a more Python-friendly format.

p4_pix2foc(*args)

Convert pixel coordinates to focal plane coordinates using distortion paper table-lookup correction.

pix2foc(*args)

Convert pixel coordinates to focal plane coordinates using the SIP polynomial distortion convention and distortion paper table-lookup correction.

pixel_to_world(*pixel_arrays)

Convert pixel coordinates to world coordinates (represented by high-level objects).

pixel_to_world_values(*pixel_arrays)

Convert pixel coordinates to world coordinates.

printwcs()

proj_plane_pixel_area()

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.

proj_plane_pixel_scales()

Calculate 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.

reorient_celestial_first()

Reorient the WCS such that the celestial axes are first, followed by the spectral axis, followed by any others.

sip_foc2pix(*args)

Convert focal plane coordinates to pixel coordinates using the SIP polynomial distortion convention.

sip_pix2foc(*args)

Convert pixel coordinates to focal plane coordinates using the SIP polynomial distortion convention.

slice(view[, numpy_order])

Slice a WCS instance using a Numpy slice.

sub(axes)

Extracts the coordinate description for a subimage from a WCS object.

swapaxes(ax0, ax1)

Swap axes in a WCS.

to_fits([relax, key])

Generate an HDUList object with all of the information stored in this object.

to_header([relax, key])

Generate an astropy.io.fits.Header object with the basic WCS and SIP information stored in this object.

to_header_string([relax])

Identical to to_header, but returns a string containing the header cards.

wcs_pix2world(*args, **kwargs)

Transforms pixel coordinates to world coordinates by doing only the basic wcslib transformation.

wcs_world2pix(*args, **kwargs)

Transforms world coordinates to pixel coordinates, using only the basic wcslib WCS transformation.

world_to_array_index(*world_objects)

Convert world coordinates (represented by Astropy objects) to array indices.

world_to_array_index_values(*world_arrays)

Convert world coordinates to array indices.

world_to_pixel(*world_objects)

Convert world coordinates (represented by Astropy objects) to pixel coordinates.

world_to_pixel_values(*world_arrays)

Convert world coordinates to pixel coordinates.

Attributes Documentation

array_shape#
axis_correlation_matrix#
axis_type_names#

World names for each coordinate axis.

Returns:
list of str

A list of names along each axis.

celestial#

A copy of the current WCS with only the celestial axes included.

cpdis1#

DistortionLookupTable

The pre-linear transformation distortion lookup table, CPDIS1.

cpdis2#

DistortionLookupTable

The pre-linear transformation distortion lookup table, CPDIS2.

det2im1#

A DistortionLookupTable object for detector to image plane correction in the x-axis.

det2im2#

A DistortionLookupTable object for detector to image plane correction in the y-axis.

has_celestial#
has_distortion#

Returns True if any distortion terms are present.

has_spectral#
has_temporal#
is_celestial#
is_spectral#
is_temporal#
low_level_wcs#
pixel_axis_names#

An iterable of strings describing the name for each pixel axis.

If an axis does not have a name, an empty string should be returned (this is the default behavior for all axes if a subclass does not override this property). Note that these names are just for display purposes and are not standardized.

pixel_bounds#
pixel_n_dim#
pixel_scale_matrix#
pixel_shape#
serialized_classes#
sip#

Get/set the Sip object for performing SIP distortion correction.

spectral#

A copy of the current WCS with only the spectral axes included.

temporal#

A copy of the current WCS with only the time axes included.

wcs#

A Wcsprm object to perform the basic wcslib WCS transformation.

world_axis_names#
world_axis_object_classes#
world_axis_object_components#
world_axis_physical_types#
world_axis_units#
world_n_dim#

Methods Documentation

all_pix2world(*args, **kwargs)[source]#

Transforms pixel coordinates to world coordinates.

Performs all of the following in series:

  • Detector to image plane correction (if present in the FITS file)

  • SIP distortion correction (if present in the FITS file)

  • distortion paper table-lookup correction (if present in the FITS file)

  • wcslib “core” WCS transformation

Parameters:
*args

There are two accepted forms for the positional arguments:

  • 2 arguments: An N x naxis array of coordinates, and an origin.

  • more than 2 arguments: An array for each axis, followed by an origin. These arrays must be broadcastable to one another.

Here, origin is the coordinate in the upper left corner of the image. In FITS and Fortran standards, this is 1. In Numpy and C standards this is 0.

For a transformation that is not two-dimensional, the two-argument form must be used.

ra_dec_orderbool, optional

When True will ensure that world coordinates are always given and returned in as (ra, dec) pairs, regardless of the order of the axes specified by the in the CTYPE keywords. Default is False.

Returns:
resultarray

Returns the sky coordinates, in degrees. If the input was a single array and origin, a single array is returned, otherwise a tuple of arrays is returned.

Raises:
MemoryError

Memory allocation failed.

SingularMatrixError

Linear transformation matrix is singular.

InconsistentAxisTypesError

Inconsistent or unrecognized coordinate axis types.

ValueError

Invalid parameter value.

ValueError

Invalid coordinate transformation parameters.

ValueError

x- and y-coordinate arrays are not the same size.

InvalidTransformError

Invalid coordinate transformation parameters.

InvalidTransformError

Ill-conditioned coordinate transformation parameters.

Notes

The order of the axes for the result is determined by the CTYPEia keywords in the FITS header, therefore it may not always be of the form (ra, dec). The lat, lng, lattyp and lngtyp members can be used to determine the order of the axes.

all_world2pix(*arg, tolerance=1.0e-4, maxiter=20, adaptive=False, detect_divergence=True, quiet=False)[source]#

Transforms world coordinates to pixel coordinates, using numerical iteration to invert the full forward transformation all_pix2world with complete distortion model.

Parameters:
*args

There are two accepted forms for the positional arguments:

  • 2 arguments: An N x naxis array of coordinates, and an origin.

  • more than 2 arguments: An array for each axis, followed by an origin. These arrays must be broadcastable to one another.

Here, origin is the coordinate in the upper left corner of the image. In FITS and Fortran standards, this is 1. In Numpy and C standards this is 0.

For a transformation that is not two-dimensional, the two-argument form must be used.

ra_dec_orderbool, optional

When True will ensure that world coordinates are always given and returned in as (ra, dec) pairs, regardless of the order of the axes specified by the in the CTYPE keywords. Default is False.

tolerancefloat, optional (default = 1.0e-4)

Tolerance of solution. Iteration terminates when the iterative solver estimates that the “true solution” is within this many pixels current estimate, more specifically, when the correction to the solution found during the previous iteration is smaller (in the sense of the L2 norm) than tolerance.

maxiterint, optional (default = 20)

Maximum number of iterations allowed to reach a solution.

quietbool, optional (default = False)

Do not throw NoConvergence exceptions when the method does not converge to a solution with the required accuracy within a specified number of maximum iterations set by maxiter parameter. Instead, simply return the found solution.

Returns:
resultarray

Returns the pixel coordinates. If the input was a single array and origin, a single array is returned, otherwise a tuple of arrays is returned.

Other Parameters:
adaptivebool, optional (default = False)

Specifies whether to adaptively select only points that did not converge to a solution within the required accuracy for the next iteration. Default is recommended for HST as well as most other instruments.

Note

The all_world2pix() uses a vectorized implementation of the method of consecutive approximations (see Notes section below) in which it iterates over all input points regardless until the required accuracy has been reached for all input points. In some cases it may be possible that almost all points have reached the required accuracy but there are only a few of input data points for which additional iterations may be needed (this depends mostly on the characteristics of the geometric distortions for a given instrument). In this situation it may be advantageous to set adaptive = True in which case all_world2pix() will continue iterating only over the points that have not yet converged to the required accuracy. However, for the HST’s ACS/WFC detector, which has the strongest distortions of all HST instruments, testing has shown that enabling this option would lead to a about 50-100% penalty in computational time (depending on specifics of the image, geometric distortions, and number of input points to be converted). Therefore, for HST and possibly instruments, it is recommended to set adaptive = False. The only danger in getting this setting wrong will be a performance penalty.

Note

When detect_divergence is True, all_world2pix() will automatically switch to the adaptive algorithm once divergence has been detected.

detect_divergencebool, optional (default = True)

Specifies whether to perform a more detailed analysis of the convergence to a solution. Normally all_world2pix() may not achieve the required accuracy if either the tolerance or maxiter arguments are too low. However, it may happen that for some geometric distortions the conditions of convergence for the method of consecutive approximations used by all_world2pix() may not be satisfied, in which case consecutive approximations to the solution will diverge regardless of the tolerance or maxiter settings.

When detect_divergence is False, these divergent points will be detected as not having achieved the required accuracy (without further details). In addition, if adaptive is False then the algorithm will not know that the solution (for specific points) is diverging and will continue iterating and trying to “improve” diverging solutions. This may result in NaN or Inf values in the return results (in addition to a performance penalties). Even when detect_divergence is False, all_world2pix(), at the end of the iterative process, will identify invalid results (NaN or Inf) as “diverging” solutions and will raise NoConvergence unless the quiet parameter is set to True.

When detect_divergence is True, all_world2pix() will detect points for which current correction to the coordinates is larger than the correction applied during the previous iteration if the requested accuracy has not yet been achieved. In this case, if adaptive is True, these points will be excluded from further iterations and if adaptive is False, all_world2pix() will automatically switch to the adaptive algorithm. Thus, the reported divergent solution will be the latest converging solution computed immediately before divergence has been detected.

Note

When accuracy has been achieved, small increases in current corrections may be possible due to rounding errors (when adaptive is False) and such increases will be ignored.

Note

Based on our testing using HST ACS/WFC images, setting detect_divergence to True will incur about 5-20% performance penalty with the larger penalty corresponding to adaptive set to True. Because the benefits of enabling this feature outweigh the small performance penalty, especially when adaptive = False, it is recommended to set detect_divergence to True, unless extensive testing of the distortion models for images from specific instruments show a good stability of the numerical method for a wide range of coordinates (even outside the image itself).

Note

Indices of the diverging inverse solutions will be reported in the divergent attribute of the raised NoConvergence exception object.

Raises:
NoConvergence

The method did not converge to a solution to the required accuracy within a specified number of maximum iterations set by the maxiter parameter. To turn off this exception, set quiet to True. Indices of the points for which the requested accuracy was not achieved (if any) will be listed in the slow_conv attribute of the raised NoConvergence exception object.

See NoConvergence documentation for more details.

MemoryError

Memory allocation failed.

SingularMatrixError

Linear transformation matrix is singular.

InconsistentAxisTypesError

Inconsistent or unrecognized coordinate axis types.

ValueError

Invalid parameter value.

ValueError

Invalid coordinate transformation parameters.

ValueError

x- and y-coordinate arrays are not the same size.

InvalidTransformError

Invalid coordinate transformation parameters.

InvalidTransformError

Ill-conditioned coordinate transformation parameters.

Notes

The order of the axes for the input world array is determined by the CTYPEia keywords in the FITS header, therefore it may not always be of the form (ra, dec). The lat, lng, lattyp, and lngtyp members can be used to determine the order of the axes.

Using the method of fixed-point iterations approximations we iterate starting with the initial approximation, which is computed using the non-distortion-aware wcs_world2pix() (or equivalent).

The all_world2pix() function uses a vectorized implementation of the method of consecutive approximations and therefore it is highly efficient (>30x) when all data points that need to be converted from sky coordinates to image coordinates are passed at once. Therefore, it is advisable, whenever possible, to pass as input a long array of all points that need to be converted to all_world2pix() instead of calling all_world2pix() for each data point. Also see the note to the adaptive parameter.

Examples

>>> import astropy.io.fits as fits
>>> import astropy.wcs as wcs
>>> import numpy as np
>>> import os
>>> filename = os.path.join(wcs.__path__[0], 'tests/data/j94f05bgq_flt.fits')
>>> hdulist = fits.open(filename)
>>> w = wcs.WCS(hdulist[('sci',1)].header, hdulist)
>>> hdulist.close()
>>> ra, dec = w.all_pix2world([1,2,3], [1,1,1], 1)
>>> print(ra)  
[ 5.52645627  5.52649663  5.52653698]
>>> print(dec)  
[-72.05171757 -72.05171276 -72.05170795]
>>> radec = w.all_pix2world([[1,1], [2,1], [3,1]], 1)
>>> print(radec)  
[[  5.52645627 -72.05171757]
 [  5.52649663 -72.05171276]
 [  5.52653698 -72.05170795]]
>>> x, y = w.all_world2pix(ra, dec, 1)
>>> print(x)  
[ 1.00000238  2.00000237  3.00000236]
>>> print(y)  
[ 0.99999996  0.99999997  0.99999997]
>>> xy = w.all_world2pix(radec, 1)
>>> print(xy)  
[[ 1.00000238  0.99999996]
 [ 2.00000237  0.99999997]
 [ 3.00000236  0.99999997]]
>>> xy = w.all_world2pix(radec, 1, maxiter=3,
...                      tolerance=1.0e-10, quiet=False)
Traceback (most recent call last):
...
NoConvergence: 'WCS.all_world2pix' failed to converge to the
requested accuracy. After 3 iterations, the solution is
diverging at least for one input point.
>>> # Now try to use some diverging data:
>>> divradec = w.all_pix2world([[1.0, 1.0],
...                             [10000.0, 50000.0],
...                             [3.0, 1.0]], 1)
>>> print(divradec)  
[[  5.52645627 -72.05171757]
 [  7.15976932 -70.8140779 ]
 [  5.52653698 -72.05170795]]
>>> # First, turn detect_divergence on:
>>> try:  
...   xy = w.all_world2pix(divradec, 1, maxiter=20,
...                        tolerance=1.0e-4, adaptive=False,
...                        detect_divergence=True,
...                        quiet=False)
... except wcs.wcs.NoConvergence as e:
...   print("Indices of diverging points: {0}"
...         .format(e.divergent))
...   print("Indices of poorly converging points: {0}"
...         .format(e.slow_conv))
...   print("Best solution:\n{0}".format(e.best_solution))
...   print("Achieved accuracy:\n{0}".format(e.accuracy))
Indices of diverging points: [1]
Indices of poorly converging points: None
Best solution:
[[  1.00000238e+00   9.99999965e-01]
 [ -1.99441636e+06   1.44309097e+06]
 [  3.00000236e+00   9.99999966e-01]]
Achieved accuracy:
[[  6.13968380e-05   8.59638593e-07]
 [  8.59526812e+11   6.61713548e+11]
 [  6.09398446e-05   8.38759724e-07]]
>>> raise e
Traceback (most recent call last):
...
NoConvergence: 'WCS.all_world2pix' failed to converge to the
requested accuracy.  After 5 iterations, the solution is
diverging at least for one input point.
>>> # This time turn detect_divergence off:
>>> try:  
...   xy = w.all_world2pix(divradec, 1, maxiter=20,
...                        tolerance=1.0e-4, adaptive=False,
...                        detect_divergence=False,
...                        quiet=False)
... except wcs.wcs.NoConvergence as e:
...   print("Indices of diverging points: {0}"
...         .format(e.divergent))
...   print("Indices of poorly converging points: {0}"
...         .format(e.slow_conv))
...   print("Best solution:\n{0}".format(e.best_solution))
...   print("Achieved accuracy:\n{0}".format(e.accuracy))
Indices of diverging points: [1]
Indices of poorly converging points: None
Best solution:
[[ 1.00000009  1.        ]
 [        nan         nan]
 [ 3.00000009  1.        ]]
Achieved accuracy:
[[  2.29417358e-06   3.21222995e-08]
 [             nan              nan]
 [  2.27407877e-06   3.13005639e-08]]
>>> raise e
Traceback (most recent call last):
...
NoConvergence: 'WCS.all_world2pix' failed to converge to the
requested accuracy.  After 6 iterations, the solution is
diverging at least for one input point.
array_index_to_world(*index_arrays)#

Convert array indices to world coordinates (represented by Astropy objects).

If a single high-level object is used to represent the world coordinates (i.e., if len(wcs.world_axis_object_classes) == 1), it is returned as-is (not in a tuple/list), otherwise a tuple of high-level objects is returned. See array_index_to_world_values for pixel indexing and ordering conventions.

array_index_to_world_values(*index_arrays)#

Convert array indices to world coordinates.

This is the same as pixel_to_world_values except that the indices should be given in (i, j) order, where for an image i is the row and j is the column (i.e. the opposite order to pixel_to_world_values).

If world_n_dim is 1, this method returns a single scalar or array, otherwise a tuple of scalars or arrays is returned.

calc_footprint(header=None, undistort=True, axes=None, center=True)[source]#

Calculates the footprint of the image on the sky.

A footprint is defined as the positions of the corners of the image on the sky after all available distortions have been applied.

Parameters:
headerHeader object, optional

Used to get NAXIS1 and NAXIS2 header and axes are mutually exclusive, alternative ways to provide the same information.

undistortbool, optional

If True, take SIP and distortion lookup table into account

axes(int, int), optional

If provided, use the given sequence as the shape of the image. Otherwise, use the NAXIS1 and NAXIS2 keywords from the header that was used to create this WCS object.

centerbool, optional

If True use the center of the pixel, otherwise use the corner.

Returns:
coord(4, 2) array of (x, y) coordinates.

The order is clockwise starting with the bottom left corner.

copy()[source]#

Return a shallow copy of the object.

Convenience method so user doesn’t have to import the copy stdlib module.

Warning

Use deepcopy instead of copy unless you know why you need a shallow copy.

deepcopy()[source]#

Return a deep copy of the object.

Convenience method so user doesn’t have to import the copy stdlib module.

det2im(*args)[source]#

Convert detector coordinates to image plane coordinates using distortion paper table-lookup correction.

The output is in absolute pixel coordinates, not relative to CRPIX.

Parameters:
*args

There are two accepted forms for the positional arguments:

  • 2 arguments: An N x 2 array of coordinates, and an origin.

  • more than 2 arguments: An array for each axis, followed by an origin. These arrays must be broadcastable to one another.

Here, origin is the coordinate in the upper left corner of the image. In FITS and Fortran standards, this is 1. In Numpy and C standards this is 0.

Returns:
resultarray

Returns the pixel coordinates. If the input was a single array and origin, a single array is returned, otherwise a tuple of arrays is returned.

Raises:
MemoryError

Memory allocation failed.

ValueError

Invalid coordinate transformation parameters.

dropaxis(dropax)[source]#

Remove an axis from the WCS.

Parameters:
wcsWCS

The WCS with naxis to be chopped to naxis-1

dropaxint

The index of the WCS to drop, counting from 0 (i.e., python convention, not FITS convention)

Returns:
WCS

A new WCS instance with one axis fewer

fix(translate_units='', naxis=None)[source]#

Perform the fix operations from wcslib, and warn about any changes it has made.

Parameters:
translate_unitsstr, optional

Specify which potentially unsafe translations of non-standard unit strings to perform. By default, performs none.

Although "S" is commonly used to represent seconds, its translation to "s" is potentially unsafe since the standard recognizes "S" formally as Siemens, however rarely that may be used. The same applies to "H" for hours (Henry), and "D" for days (Debye).

This string controls what to do in such cases, and is case-insensitive.

  • If the string contains "s", translate "S" to "s".

  • If the string contains "h", translate "H" to "h".

  • If the string contains "d", translate "D" to "d".

Thus '' doesn’t do any unsafe translations, whereas 'shd' does all of them.

naxisint array, optional

Image axis lengths. If this array is set to zero or None, then cylfix will not be invoked.

footprint_contains(coord, **kwargs)[source]#

Determines if a given SkyCoord is contained in the wcs footprint.

Parameters:
coordSkyCoord

The coordinate to check if it is within the wcs coordinate.

**kwargs

Additional arguments to pass to to_pixel

Returns:
responsebool

True means the WCS footprint contains the coordinate, False means it does not.

footprint_to_file(filename='footprint.reg', color='green', width=2, coordsys=None)[source]#

Writes out a ds9 style regions file. It can be loaded directly by ds9.

Parameters:
filenamestr, optional

Output file name - default is 'footprint.reg'

colorstr, optional

Color to use when plotting the line.

widthint, optional

Width of the region line.

coordsysstr, optional

Coordinate system. If not specified (default), the radesys value is used. For all possible values, see http://ds9.si.edu/doc/ref/region.html#RegionFileFormat

get_axis_types()[source]#

Similar to self.wcsprm.axis_types but provides the information in a more Python-friendly format.

Returns:
resultlist of dict

Returns a list of dictionaries, one for each axis, each containing attributes about the type of that axis.

Each dictionary has the following keys:

  • ‘coordinate_type’:

    • None: Non-specific coordinate type.

    • ‘stokes’: Stokes coordinate.

    • ‘celestial’: Celestial coordinate (including CUBEFACE).

    • ‘spectral’: Spectral coordinate.

  • ‘scale’:

    • ‘linear’: Linear axis.

    • ‘quantized’: Quantized axis (STOKES, CUBEFACE).

    • ‘non-linear celestial’: Non-linear celestial axis.

    • ‘non-linear spectral’: Non-linear spectral axis.

    • ‘logarithmic’: Logarithmic axis.

    • ‘tabular’: Tabular axis.

  • ‘group’

    • Group number, e.g. lookup table number

  • ‘number’

    • For celestial axes:

      • 0: Longitude coordinate.

      • 1: Latitude coordinate.

      • 2: CUBEFACE number.

    • For lookup tables:

      • the axis number in a multidimensional table.

CTYPEia in "4-3" form with unrecognized algorithm code will generate an error.

p4_pix2foc(*args)[source]#

Convert pixel coordinates to focal plane coordinates using distortion paper table-lookup correction.

The output is in absolute pixel coordinates, not relative to CRPIX.

Parameters:
*args

There are two accepted forms for the positional arguments:

  • 2 arguments: An N x 2 array of coordinates, and an origin.

  • more than 2 arguments: An array for each axis, followed by an origin. These arrays must be broadcastable to one another.

Here, origin is the coordinate in the upper left corner of the image. In FITS and Fortran standards, this is 1. In Numpy and C standards this is 0.

Returns:
resultarray

Returns the focal coordinates. If the input was a single array and origin, a single array is returned, otherwise a tuple of arrays is returned.

Raises:
MemoryError

Memory allocation failed.

ValueError

Invalid coordinate transformation parameters.

pix2foc(*args)[source]#

Convert pixel coordinates to focal plane coordinates using the SIP polynomial distortion convention and distortion paper table-lookup correction.

The output is in absolute pixel coordinates, not relative to CRPIX.

Parameters:
*args

There are two accepted forms for the positional arguments:

  • 2 arguments: An N x 2 array of coordinates, and an origin.

  • more than 2 arguments: An array for each axis, followed by an origin. These arrays must be broadcastable to one another.

Here, origin is the coordinate in the upper left corner of the image. In FITS and Fortran standards, this is 1. In Numpy and C standards this is 0.

Returns:
resultarray

Returns the focal coordinates. If the input was a single array and origin, a single array is returned, otherwise a tuple of arrays is returned.

Raises:
MemoryError

Memory allocation failed.

ValueError

Invalid coordinate transformation parameters.

pixel_to_world(*pixel_arrays)#

Convert pixel coordinates to world coordinates (represented by high-level objects).

If a single high-level object is used to represent the world coordinates (i.e., if len(wcs.world_axis_object_classes) == 1), it is returned as-is (not in a tuple/list), otherwise a tuple of high-level objects is returned. See pixel_to_world_values for pixel indexing and ordering conventions.

pixel_to_world_values(*pixel_arrays)#

Convert pixel coordinates to world coordinates.

This method takes pixel_n_dim scalars or arrays as input, and pixel coordinates should be zero-based. Returns world_n_dim scalars or arrays in units given by world_axis_units. Note that pixel coordinates are assumed to be 0 at the center of the first pixel in each dimension. If a pixel is in a region where the WCS is not defined, NaN can be returned. The coordinates should be specified in the (x, y) order, where for an image, x is the horizontal coordinate and y is the vertical coordinate.

If world_n_dim is 1, this method returns a single scalar or array, otherwise a tuple of scalars or arrays is returned.

printwcs()[source]#
proj_plane_pixel_area()[source]#

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.

Note

This function is concerned only about the transformation “image plane”->”projection plane” and not about the transformation “celestial sphere”->”projection plane”->”image plane”. Therefore, this function ignores distortions arising due to non-linear nature of most projections.

Note

This method only returns sensible answers if the WCS contains celestial axes, i.e., the celestial WCS object.

Returns:
areaQuantity

Area (in the projection plane) of the pixel at CRPIX location.

Raises:
ValueError

Pixel area is defined only for 2D pixels. Most likely the cd matrix of the celestial WCS is not a square matrix of second order.

Notes

Depending on the application, square root of the pixel area can be used to represent a single pixel scale of an equivalent square pixel whose area is equal to the area of a generally non-square pixel.

proj_plane_pixel_scales()[source]#

Calculate 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.

Note

This method is concerned only about the transformation “image plane”->”projection plane” and not about the transformation “celestial sphere”->”projection plane”->”image plane”. Therefore, this function ignores distortions arising due to non-linear nature of most projections.

Note

This method only returns sensible answers if the WCS contains celestial axes, i.e., the celestial WCS object.

Returns:
scalelist of Quantity

A vector of projection plane increments corresponding to each pixel side (axis).

reorient_celestial_first()[source]#

Reorient the WCS such that the celestial axes are first, followed by the spectral axis, followed by any others. Assumes at least celestial axes are present.

sip_foc2pix(*args)[source]#

Convert focal plane coordinates to pixel coordinates using the SIP polynomial distortion convention.

FITS WCS distortion paper table lookup distortion correction is not applied, even if that information existed in the FITS file that initialized this WCS object.

Parameters:
*args

There are two accepted forms for the positional arguments:

  • 2 arguments: An N x 2 array of coordinates, and an origin.

  • more than 2 arguments: An array for each axis, followed by an origin. These arrays must be broadcastable to one another.

Here, origin is the coordinate in the upper left corner of the image. In FITS and Fortran standards, this is 1. In Numpy and C standards this is 0.

Returns:
resultarray

Returns the pixel coordinates. If the input was a single array and origin, a single array is returned, otherwise a tuple of arrays is returned.

Raises:
MemoryError

Memory allocation failed.

ValueError

Invalid coordinate transformation parameters.

sip_pix2foc(*args)[source]#

Convert pixel coordinates to focal plane coordinates using the SIP polynomial distortion convention.

The output is in pixel coordinates, relative to CRPIX.

FITS WCS distortion paper table lookup correction is not applied, even if that information existed in the FITS file that initialized this WCS object. To correct for that, use pix2foc or p4_pix2foc.

Parameters:
*args

There are two accepted forms for the positional arguments:

  • 2 arguments: An N x 2 array of coordinates, and an origin.

  • more than 2 arguments: An array for each axis, followed by an origin. These arrays must be broadcastable to one another.

Here, origin is the coordinate in the upper left corner of the image. In FITS and Fortran standards, this is 1. In Numpy and C standards this is 0.

Returns:
resultarray

Returns the focal coordinates. If the input was a single array and origin, a single array is returned, otherwise a tuple of arrays is returned.

Raises:
MemoryError

Memory allocation failed.

ValueError

Invalid coordinate transformation parameters.

slice(view, numpy_order=True)[source]#

Slice a WCS instance using a Numpy slice. The order of the slice should be reversed (as for the data) compared to the natural WCS order.

Parameters:
viewtuple

A tuple containing the same number of slices as the WCS system. The step method, the third argument to a slice, is not presently supported.

numpy_orderbool

Use numpy order, i.e. slice the WCS so that an identical slice applied to a numpy array will slice the array and WCS in the same way. If set to False, the WCS will be sliced in FITS order, meaning the first slice will be applied to the last numpy index but the first WCS axis.

Returns:
wcs_newWCS

A new resampled WCS axis

sub(axes)[source]#

Extracts the coordinate description for a subimage from a WCS object.

The world coordinate system of the subimage must be separable in the sense that the world coordinates at any point in the subimage must depend only on the pixel coordinates of the axes extracted. In practice, this means that the PCi_ja matrix of the original image must not contain non-zero off-diagonal terms that associate any of the subimage axes with any of the non-subimage axes.

sub can also add axes to a wcsprm object. The new axes will be created using the defaults set by the Wcsprm constructor which produce a simple, unnamed, linear axis with world coordinates equal to the pixel coordinate. These default values can be changed before invoking set.

Parameters:
axesint or a sequence.
  • If an int, include the first N axes in their original order.

  • If a sequence, may contain a combination of image axis numbers (1-relative) or special axis identifiers (see below). Order is significant; axes[0] is the axis number of the input image that corresponds to the first axis in the subimage, etc. Use an axis number of 0 to create a new axis using the defaults.

  • If 0, [] or None, do a deep copy.

Coordinate axes types may be specified using either strings or special integer constants. The available types are:

  • 'longitude' / WCSSUB_LONGITUDE: Celestial longitude

  • 'latitude' / WCSSUB_LATITUDE: Celestial latitude

  • 'cubeface' / WCSSUB_CUBEFACE: Quadcube CUBEFACE axis

  • 'spectral' / WCSSUB_SPECTRAL: Spectral axis

  • 'stokes' / WCSSUB_STOKES: Stokes axis

  • 'temporal' / WCSSUB_TIME: Time axis (requires WCSLIB version 7.8 or greater)

  • 'celestial' / WCSSUB_CELESTIAL: An alias for the combination of 'longitude', 'latitude' and 'cubeface'.

Returns:
new_wcsWCS object
Raises:
MemoryError

Memory allocation failed.

InvalidSubimageSpecificationError

Invalid subimage specification (no spectral axis).

NonseparableSubimageCoordinateSystemError

Non-separable subimage coordinate system.

Notes

Combinations of subimage axes of particular types may be extracted in the same order as they occur in the input image by combining the integer constants with the ‘binary or’ (|) operator. For example:

wcs.sub([WCSSUB_LONGITUDE | WCSSUB_LATITUDE | WCSSUB_SPECTRAL])

would extract the longitude, latitude, and spectral axes in the same order as the input image. If one of each were present, the resulting object would have three dimensions.

For convenience, WCSSUB_CELESTIAL is defined as the combination WCSSUB_LONGITUDE | WCSSUB_LATITUDE | WCSSUB_CUBEFACE.

The codes may also be negated to extract all but the types specified, for example:

wcs.sub([
  WCSSUB_LONGITUDE,
  WCSSUB_LATITUDE,
  WCSSUB_CUBEFACE,
  -(WCSSUB_SPECTRAL | WCSSUB_STOKES)])

The last of these specifies all axis types other than spectral or Stokes. Extraction is done in the order specified by axes, i.e. a longitude axis (if present) would be extracted first (via axes[0]) and not subsequently (via axes[3]). Likewise for the latitude and cubeface axes in this example.

The number of dimensions in the returned object may be less than or greater than the length of axes. However, it will never exceed the number of axes in the input image.

swapaxes(ax0, ax1)[source]#

Swap axes in a WCS.

Parameters:
wcsWCS

The WCS to have its axes swapped

ax0int
ax1int

The indices of the WCS to be swapped, counting from 0 (i.e., python convention, not FITS convention)

Returns:
WCS

A new WCS instance with the same number of axes, but two swapped

to_fits(relax=False, key=None)[source]#

Generate an HDUList object with all of the information stored in this object. This should be logically identical to the input FITS file, but it will be normalized in a number of ways.

See to_header for some warnings about the output produced.

Parameters:
relaxbool or int, optional

Degree of permissiveness:

  • False (default): Write all extensions that are considered to be safe and recommended.

  • True: Write all recognized informal extensions of the WCS standard.

  • int: a bit field selecting specific extensions to write. See Header-writing relaxation constants for details.

keystr

The name of a particular WCS transform to use. This may be either ' ' or 'A'-'Z' and corresponds to the "a" part of the CTYPEia cards.

Returns:
hdulistHDUList
to_header(relax=None, key=None)[source]#

Generate an astropy.io.fits.Header object with the basic WCS and SIP information stored in this object. This should be logically identical to the input FITS file, but it will be normalized in a number of ways.

Warning

This function does not write out FITS WCS distortion paper information, since that requires multiple FITS header data units. To get a full representation of everything in this object, use to_fits.

Parameters:
relaxbool or int, optional

Degree of permissiveness:

  • False (default): Write all extensions that are considered to be safe and recommended.

  • True: Write all recognized informal extensions of the WCS standard.

  • int: a bit field selecting specific extensions to write. See Header-writing relaxation constants for details.

If the relax keyword argument is not given and any keywords were omitted from the output, an AstropyWarning is displayed. To override this, explicitly pass a value to relax.

keystr

The name of a particular WCS transform to use. This may be either ' ' or 'A'-'Z' and corresponds to the "a" part of the CTYPEia cards.

Returns:
headerastropy.io.fits.Header

Notes

The output header will almost certainly differ from the input in a number of respects:

  1. The output header only contains WCS-related keywords. In particular, it does not contain syntactically-required keywords such as SIMPLE, NAXIS, BITPIX, or END.

  2. Deprecated (e.g. CROTAn) or non-standard usage will be translated to standard (this is partially dependent on whether fix was applied).

  3. Quantities will be converted to the units used internally, basically SI with the addition of degrees.

  4. Floating-point quantities may be given to a different decimal precision.

  5. Elements of the PCi_j matrix will be written if and only if they differ from the unit matrix. Thus, if the matrix is unity then no elements will be written.

  6. Additional keywords such as WCSAXES, CUNITia, LONPOLEa and LATPOLEa may appear.

  7. The original keycomments will be lost, although to_header tries hard to write meaningful comments.

  8. Keyword order may be changed.

to_header_string(relax=None)[source]#

Identical to to_header, but returns a string containing the header cards.

wcs_pix2world(*args, **kwargs)[source]#

Transforms pixel coordinates to world coordinates by doing only the basic wcslib transformation.

No SIP or distortion paper table lookup correction is applied. To perform distortion correction, see all_pix2world, sip_pix2foc, p4_pix2foc, or pix2foc.

Parameters:
*args

There are two accepted forms for the positional arguments:

  • 2 arguments: An N x naxis array of coordinates, and an origin.

  • more than 2 arguments: An array for each axis, followed by an origin. These arrays must be broadcastable to one another.

Here, origin is the coordinate in the upper left corner of the image. In FITS and Fortran standards, this is 1. In Numpy and C standards this is 0.

For a transformation that is not two-dimensional, the two-argument form must be used.

ra_dec_orderbool, optional

When True will ensure that world coordinates are always given and returned in as (ra, dec) pairs, regardless of the order of the axes specified by the in the CTYPE keywords. Default is False.

Returns:
resultarray

Returns the world coordinates, in degrees. If the input was a single array and origin, a single array is returned, otherwise a tuple of arrays is returned.

Raises:
MemoryError

Memory allocation failed.

SingularMatrixError

Linear transformation matrix is singular.

InconsistentAxisTypesError

Inconsistent or unrecognized coordinate axis types.

ValueError

Invalid parameter value.

ValueError

Invalid coordinate transformation parameters.

ValueError

x- and y-coordinate arrays are not the same size.

InvalidTransformError

Invalid coordinate transformation parameters.

InvalidTransformError

Ill-conditioned coordinate transformation parameters.

Notes

The order of the axes for the result is determined by the CTYPEia keywords in the FITS header, therefore it may not always be of the form (ra, dec). The lat, lng, lattyp and lngtyp members can be used to determine the order of the axes.

wcs_world2pix(*args, **kwargs)[source]#

Transforms world coordinates to pixel coordinates, using only the basic wcslib WCS transformation. No SIP or distortion paper table lookup transformation is applied.

Parameters:
*args

There are two accepted forms for the positional arguments:

  • 2 arguments: An N x naxis array of coordinates, and an origin.

  • more than 2 arguments: An array for each axis, followed by an origin. These arrays must be broadcastable to one another.

Here, origin is the coordinate in the upper left corner of the image. In FITS and Fortran standards, this is 1. In Numpy and C standards this is 0.

For a transformation that is not two-dimensional, the two-argument form must be used.

ra_dec_orderbool, optional

When True will ensure that world coordinates are always given and returned in as (ra, dec) pairs, regardless of the order of the axes specified by the in the CTYPE keywords. Default is False.

Returns:
resultarray

Returns the pixel coordinates. If the input was a single array and origin, a single array is returned, otherwise a tuple of arrays is returned.

Raises:
MemoryError

Memory allocation failed.

SingularMatrixError

Linear transformation matrix is singular.

InconsistentAxisTypesError

Inconsistent or unrecognized coordinate axis types.

ValueError

Invalid parameter value.

ValueError

Invalid coordinate transformation parameters.

ValueError

x- and y-coordinate arrays are not the same size.

InvalidTransformError

Invalid coordinate transformation parameters.

InvalidTransformError

Ill-conditioned coordinate transformation parameters.

Notes

The order of the axes for the input world array is determined by the CTYPEia keywords in the FITS header, therefore it may not always be of the form (ra, dec). The lat, lng, lattyp and lngtyp members can be used to determine the order of the axes.

world_to_array_index(*world_objects)#

Convert world coordinates (represented by Astropy objects) to array indices.

If pixel_n_dim is 1, this method returns a single scalar or array, otherwise a tuple of scalars or arrays is returned. See world_to_array_index_values for pixel indexing and ordering conventions. The indices should be returned as rounded integers.

world_to_array_index_values(*world_arrays)#

Convert world coordinates to array indices.

This is the same as world_to_pixel_values except that the indices should be returned in (i, j) order, where for an image i is the row and j is the column (i.e. the opposite order to pixel_to_world_values). The indices should be returned as rounded integers.

If pixel_n_dim is 1, this method returns a single scalar or array, otherwise a tuple of scalars or arrays is returned.

world_to_pixel(*world_objects)#

Convert world coordinates (represented by Astropy objects) to pixel coordinates.

If pixel_n_dim is 1, this method returns a single scalar or array, otherwise a tuple of scalars or arrays is returned. See world_to_pixel_values for pixel indexing and ordering conventions.

world_to_pixel_values(*world_arrays)#

Convert world coordinates to pixel coordinates.

This method takes world_n_dim scalars or arrays as input in units given by world_axis_units. Returns pixel_n_dim scalars or arrays. Note that pixel coordinates are assumed to be 0 at the center of the first pixel in each dimension. If a world coordinate does not have a matching pixel coordinate, NaN can be returned. The coordinates should be returned in the (x, y) order, where for an image, x is the horizontal coordinate and y is the vertical coordinate.

If pixel_n_dim is 1, this method returns a single scalar or array, otherwise a tuple of scalars or arrays is returned.