Bases: astropy.wcs.WCSBase
WCS objects perform standard WCS transformations, and correct for SIP and Paper IV table-lookup distortions, based on the WCS keywords and supplementary data read from a FITS file.
Parameters: | header : astropy.io.fits header object, string, dict-like, or None, optional
fobj : An astropy.io.fits file (hdulist) object, optional
key : str, optional
minerr : float, optional
relax : bool or int, optional
naxis : int or sequence, optional keysel : sequence of flags, optional
colsel : sequence of int, optional
fix : bool, optional
translate_units : str, optional |
---|---|
Raises: | MemoryError
ValueError
KeyError
AssertionError
|
Notes
astropy.wcs supports arbitrary n dimensions for the core WCS (the transformations handled by WCSLIB). However, the Paper IV 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 Paper IV 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.
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.
When the header includes duplicate keywords, in most cases the last encountered is used.
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
axis_type_names | World names for each coordinate axis |
Methods Summary
all_pix2world(*args, **kwargs) | Transforms pixel coordinates to world coordinates. |
all_world2pix(*args, **kwargs) | Transforms world coordinates to pixel coordinates, using numerical iteration to invert the method all_pix2world within a tolerance of 1e-6 pixels. |
calcFootprint(*args, **kwargs) | Deprecated since version 0.4. |
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 Paper IV table-lookup distortion 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_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 Paper IV table-lookup distortion correction. |
pix2foc(*args) | Convert pixel coordinates to focal plane coordinates using the SIP polynomial distortion convention and Paper IV table-lookup distortion correction. |
printwcs() | |
reorient_celestial_first() | Reorient the WCS such that the celestial axes are first, followed by the spectral axis, followed by any others. |
rotateCD(theta) | |
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 astropy.io.fits.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. |
Attributes Documentation
World names for each coordinate axis
Returns: | A list of names along each axis |
---|
Methods Documentation
Transforms pixel coordinates to world coordinates.
Performs all of the following in order:
Parameters: | args : flexible
ra_dec_order : bool, optional |
---|---|
Returns: | result : array
|
Raises: | MemoryError
SingularMatrixError
InconsistentAxisTypesError
ValueError
ValueError
ValueError
InvalidTransformError
InvalidTransformError
|
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.
Transforms world coordinates to pixel coordinates, using numerical iteration to invert the method all_pix2world within a tolerance of 1e-6 pixels.
Note that to use this function, you must have Scipy installed.
Parameters: | args : flexible
ra_dec_order : bool, optional tolerance : float, optional
|
---|---|
Returns: | result : array
|
Raises: | MemoryError
SingularMatrixError
InconsistentAxisTypesError
ValueError
ValueError
ValueError
InvalidTransformError
InvalidTransformError
|
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.
Deprecated since version 0.4: The calcFootprint function is deprecated and may be removed in a future version. Use calc_footprint instead.
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: | header : Header object, optional undistort : bool, optional
axes : length 2 sequence ints, optional
center : bool, optional
|
---|---|
Returns: | coord : (4, 2) array of (x, y) coordinates.
|
Return a shallow copy of the object.
Convenience method so user doesn’t have to import the copy stdlib module.
Return a deep copy of the object.
Convenience method so user doesn’t have to import the copy stdlib module.
Convert detector coordinates to image plane coordinates using Paper IV table-lookup distortion correction.
The output is in absolute pixel coordinates, not relative to CRPIX.
Parameters: | args : flexible
|
---|---|
Returns: | result : array
|
Raises: | MemoryError
ValueError
|
Remove an axis from the WCS.
Parameters: | wcs : WCS
dropax : int
|
---|---|
Returns: | A new WCS instance with one axis fewer |
Perform the fix operations from wcslib, and warn about any changes it has made.
Parameters: | translate_units : str, optional
naxis : int array[naxis], optional
|
---|
Writes out a ds9 style regions file. It can be loaded directly by ds9.
Parameters: | filename : str, optional
color : str, optional
width : int, optional
|
---|
Similar to self.wcsprm.axis_types but provides the information in a more Python-friendly format.
Returns: | result : list of dicts
|
---|
Convert pixel coordinates to focal plane coordinates using Paper IV table-lookup distortion correction.
The output is in absolute pixel coordinates, not relative to CRPIX.
Parameters: | args : flexible
|
---|---|
Returns: | result : array
|
Raises: | MemoryError
ValueError
|
Convert pixel coordinates to focal plane coordinates using the SIP polynomial distortion convention and Paper IV table-lookup distortion correction.
The output is in absolute pixel coordinates, not relative to CRPIX.
Parameters: | args : flexible
|
---|---|
Returns: | result : array
|
Raises: | MemoryError
ValueError
|
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.
Convert focal plane coordinates to pixel coordinates using the SIP polynomial distortion convention.
Paper IV table lookup distortion correction is not applied, even if that information existed in the FITS file that initialized this WCS object.
Parameters: | args : flexible
|
---|---|
Returns: | result : array
|
Raises: | MemoryError
ValueError
|
Convert pixel coordinates to focal plane coordinates using the SIP polynomial distortion convention.
The output is in pixel coordinates, relative to CRPIX.
Paper IV table lookup distortion 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 : flexible
|
---|---|
Returns: | result : array
|
Raises: | MemoryError
ValueError
|
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: | view : tuple
numpy_order : bool
|
---|---|
Returns: | wcs_new : WCS
|
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: | axes : int or a sequence.
|
---|---|
Returns: | new_wcs : WCS object |
Raises: | MemoryError
InvalidSubimageSpecificationError
NonseparableSubimageCoordinateSystem
|
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.
Swap axes in a WCS.
Parameters: | wcs: `~astropy.wcs.WCS`
ax0: int ax1: int
|
---|---|
Returns: | A new WCS instance with the same number of axes, but two swapped |
Generate an astropy.io.fits.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: | relax : bool or int, optional
key : str
|
---|---|
Returns: | hdulist : astropy.io.fits.HDUList |
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 Paper IV distortion information, since that requires multiple FITS header data units. To get a full representation of everything in this object, use to_fits.
Parameters: | relax : bool or int, optional
key : str
|
---|---|
Returns: | header : astropy.io.fits.Header |
Notes
The output header will almost certainly differ from the input in a number of respects:
- The output header only contains WCS-related keywords. In particular, it does not contain syntactically-required keywords such as SIMPLE, NAXIS, BITPIX, or END.
- Deprecated (e.g. CROTAn) or non-standard usage will be translated to standard (this is partially dependent on whether fix was applied).
- Quantities will be converted to the units used internally, basically SI with the addition of degrees.
- Floating-point quantities may be given to a different decimal precision.
- 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.
- Additional keywords such as WCSAXES, CUNITia, LONPOLEa and LATPOLEa may appear.
- The original keycomments will be lost, although to_header tries hard to write meaningful comments.
- Keyword order may be changed.
Identical to to_header, but returns a string containing the header cards.
Transforms pixel coordinates to world coordinates by doing only the basic wcslib transformation.
No SIP or Paper IV table lookup distortion correction is applied. To perform distortion correction, see all_pix2world, sip_pix2foc, p4_pix2foc, or pix2foc.
Parameters: | args : flexible
ra_dec_order : bool, optional |
---|---|
Returns: | result : array
|
Raises: | MemoryError
SingularMatrixError
InconsistentAxisTypesError
ValueError
ValueError
ValueError
InvalidTransformError
InvalidTransformError
|
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.
Transforms world coordinates to pixel coordinates, using only the basic wcslib WCS transformation. No SIP or Paper IV table lookup distortion is applied.
Parameters: | args : flexible
ra_dec_order : bool, optional |
---|---|
Returns: | result : array
|
Raises: | MemoryError
SingularMatrixError
InconsistentAxisTypesError
ValueError
ValueError
ValueError
InvalidTransformError
InvalidTransformError
|
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.