WCS

class gwcs.wcs.WCS(forward_transform=None, input_frame='detector', output_frame=None, name='')[source]

Bases: GWCSAPIMixin

Basic WCS class.

Parameters:
forward_transformModel or a list

The transform between input_frame and output_frame. A list of (frame, transform) tuples where frame is the starting frame and transform is the transform from this frame to the next one or output_frame. The last tuple is (transform, None), where None indicates the end of the pipeline.

input_framestr, CoordinateFrame

A coordinates object or a string name.

output_framestr, CoordinateFrame

A coordinates object or a string name.

namestr

a name for this WCS

Attributes Summary

available_frames

List all frames in this WCS object.

backward_transform

Return the total backward transform if available - from output to input coordinate system.

bounding_box

Return the range of acceptable values for each input axis.

forward_transform

Return the total forward transform - from input to output coordinate frame.

input_frame

Return the input coordinate frame.

name

Return the name for this WCS.

output_frame

Return the output coordinate frame.

pipeline

Return the pipeline structure.

unit

The unit of the coordinates in the output coordinate system.

Methods Summary

__call__(*args, **kwargs)

Executes the forward transform.

attach_compound_bounding_box(cbbox, ...)

fix_inputs(fixed)

Return a new unique WCS by fixing inputs to constant values.

footprint([bounding_box, center, axis_type])

Return the footprint in world coordinates.

get_transform(from_frame, to_frame)

Return a transform between two coordinate frames.

in_image(*args, **kwargs)

This method tests if one or more of the input world coordinates are contained within forward transformation's image and that it maps to the domain of definition of the forward transformation.

insert_frame(input_frame, transform, ...)

Insert a new frame into an existing pipeline.

insert_transform(frame, transform[, after])

Insert a transform before (default) or after a coordinate frame.

invert(*args, **kwargs)

Invert coordinates from output frame to input frame using analytical or user-supplied inverse.

numerical_inverse(*args[, tolerance, ...])

Invert coordinates from output frame to input frame using numerical inverse.

set_transform(from_frame, to_frame, transform)

Set/replace the transform between two coordinate frames.

to_fits([bounding_box, max_pix_error, ...])

Construct a FITS WCS -TAB-based approximation to the WCS in the form of a FITS header and a binary table extension.

to_fits_sip([bounding_box, max_pix_error, ...])

Construct a SIP-based approximation to the WCS for the axes corresponding to the CelestialFrame in the form of a FITS header.

to_fits_tab([bounding_box, bin_ext_name, ...])

Construct a FITS WCS -TAB-based approximation to the WCS in the form of a FITS header and a binary table extension.

transform(from_frame, to_frame, *args, **kwargs)

Transform positions between two frames.

Attributes Documentation

available_frames

List all frames in this WCS object.

Returns:
available_framesdict

{frame_name: frame_object or None}

backward_transform

Return the total backward transform if available - from output to input coordinate system.

Raises:
NotImplementedError

An analytical inverse does not exist.

bounding_box

Return the range of acceptable values for each input axis. The order of the axes is axes_order.

forward_transform

Return the total forward transform - from input to output coordinate frame.

input_frame

Return the input coordinate frame.

name

Return the name for this WCS.

output_frame

Return the output coordinate frame.

pipeline

Return the pipeline structure.

unit

The unit of the coordinates in the output coordinate system.

Methods Documentation

__call__(*args, **kwargs)[source]

Executes the forward transform.

argsfloat or array-like

Inputs in the input coordinate system, separate inputs for each dimension.

with_unitsbool

If True returns a SkyCoord or SpectralCoord object, by using the units of the output cooridnate frame. Optional, default=False.

with_bounding_boxbool, optional

If True(default) values in the result which correspond to any of the inputs being outside the bounding_box are set to fill_value.

fill_valuefloat, optional

Output value for inputs outside the bounding_box (default is np.nan).

attach_compound_bounding_box(cbbox, selector_args)[source]
fix_inputs(fixed)[source]

Return a new unique WCS by fixing inputs to constant values.

Parameters:
fixeddict

Keyword arguments with fixed values corresponding to self.selector.

Returns:
new_wcsWCS

A new unique WCS corresponding to the values in fixed.

Examples

>>> w = WCS(pipeline, selector={"spectral_order": [1, 2]}) 
>>> new_wcs = w.set_inputs(spectral_order=2) 
>>> new_wcs.inputs 
    ("x", "y")
footprint(bounding_box=None, center=False, axis_type='all')[source]

Return the footprint in world coordinates.

Parameters:
bounding_boxtuple of floats: (start, stop)

prop: bounding_box

centerbool

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

axis_typestr

A supported output_frame.axes_type or "all" (default). One of ['spatial', 'spectral', 'temporal'] or a custom type.

Returns:
coordndarray

Array of coordinates in the output_frame mapping corners to the output frame. For spatial coordinates the order is clockwise, starting from the bottom left corner.

get_transform(from_frame, to_frame)[source]

Return a transform between two coordinate frames.

Parameters:
from_framestr or CoordinateFrame

Initial coordinate frame name of object.

to_framestr, or instance of CoordinateFrame

End coordinate frame name or object.

Returns:
transformModel

Transform between two frames.

in_image(*args, **kwargs)[source]

This method tests if one or more of the input world coordinates are contained within forward transformation’s image and that it maps to the domain of definition of the forward transformation. In practical terms, this function tests that input world coordinate(s) can be converted to input frame and that it is within the forward transformation’s bounding_box when defined.

Parameters:
argsfloat, array like, SkyCoord or

Unit coordinates to be inverted

kwargsdict

keyword arguments to be passed either to backward_transform (when defined) or to the iterative invert method.

Returns:
resultbool, numpy.ndarray

A single boolean value or an array of boolean values with True indicating that the WCS footprint contains the coordinate and False if input is outside the footprint.

insert_frame(input_frame, transform, output_frame)[source]

Insert a new frame into an existing pipeline. This frame must be anchored to a frame already in the pipeline by a transform. This existing frame is identified solely by its name, although an entire CoordinateFrame can be passed (e.g., the input_frame or output_frame attribute). This frame is never modified.

Parameters:
input_framestr or CoordinateFrame

Coordinate frame at start of new transform

transformModel

New transform to be inserted in the pipeline

output_frame: str or `~gwcs.coordinate_frames.CoordinateFrame`

Coordinate frame at end of new transform

insert_transform(frame, transform, after=False)[source]

Insert a transform before (default) or after a coordinate frame.

Append (or prepend) a transform to the transform connected to frame.

Parameters:
framestr or CoordinateFrame

Coordinate frame which sets the point of insertion.

transformModel

New transform to be inserted in the pipeline

afterbool

If True, the new transform is inserted in the pipeline immediately after frame.

invert(*args, **kwargs)[source]

Invert coordinates from output frame to input frame using analytical or user-supplied inverse. When neither analytical nor user-supplied inverses are defined, a numerical solution will be attempted using numerical_inverse().

Note

Currently numerical inverse is implemented only for 2D imaging WCS.

Parameters:
argsfloat, array like, SkyCoord or Unit

Coordinates to be inverted. The number of arguments must be equal to the number of world coordinates given by world_n_dim.

with_bounding_boxbool, optional

If True (default) values in the result which correspond to any of the inputs being outside the bounding_box are set to fill_value.

fill_valuefloat, optional

Output value for inputs outside the bounding_box (default is np.nan).

with_unitsbool, optional

If True returns a SkyCoord or SpectralCoord object, by using the units of the output cooridnate frame. Default is False.

Returns:
resulttuple or value

Returns a tuple of scalar or array values for each axis. Unless input_frame.naxes == 1 when it shall return the value.

Other Parameters:
kwargsdict

Keyword arguments to be passed to numerical_inverse() (when defined) or to the iterative invert method.

numerical_inverse(*args, tolerance=1e-05, maxiter=50, adaptive=True, detect_divergence=True, quiet=True, with_bounding_box=True, fill_value=nan, with_units=False, **kwargs)[source]

Invert coordinates from output frame to input frame using numerical inverse.

Note

Currently numerical inverse is implemented only for 2D imaging WCS.

Note

This method uses a combination of vectorized fixed-point iterations algorithm and scipy.optimize.root. The later is used for input coordinates for which vectorized algorithm diverges.

Parameters:
argsfloat, array like, SkyCoord or Unit

Coordinates to be inverted. The number of arguments must be equal to the number of world coordinates given by world_n_dim.

with_bounding_boxbool, optional

If True (default) values in the result which correspond to any of the inputs being outside the bounding_box are set to fill_value.

fill_valuefloat, optional

Output value for inputs outside the bounding_box (default is np.nan).

with_unitsbool, optional

If True returns a SkyCoord or SpectralCoord object, by using the units of the output cooridnate frame. Default is False.

tolerancefloat, optional

Absolute 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. Default tolerance is 1.0e-5.

maxiterint, optional

Maximum number of iterations allowed to reach a solution. Default is 50.

quietbool, optional

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. Default is True.

adaptivebool, optional

Specifies whether to adaptively select only points that did not converge to a solution within the required accuracy for the next iteration. Default (True) is recommended.

Note

The numerical_inverse() 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 numerical_inverse() will continue iterating only over the points that have not yet converged to the required accuracy.

Note

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

detect_divergencebool, optional

Specifies whether to perform a more detailed analysis of the convergence to a solution. Normally numerical_inverse() 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 the method of consecutive approximations used by numerical_inverse() 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, numerical_inverse(), 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 (default), numerical_inverse() 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, numerical_inverse() 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 JWST NIRCAM images, setting detect_divergence to True will incur about 5-10% 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.

Returns:
resulttuple

Returns a tuple of scalar or array values for each axis.

Raises:
NoConvergence

The iterative 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.

NotImplementedError

Numerical inverse has not been implemented for this WCS.

ValueError

Invalid argument values.

Examples

>>> from astropy.utils.data import get_pkg_data_filename
>>> from gwcs import NoConvergence
>>> import asdf
>>> import numpy as np
>>> filename = get_pkg_data_filename('data/nircamwcs.asdf', package='gwcs.tests')
>>> with asdf.open(filename, copy_arrays=True, lazy_load=False, ignore_missing_extensions=True) as af:
...    w = af.tree['wcs']
>>> ra, dec = w([1,2,3], [1,1,1])
>>> assert np.allclose(ra, [5.927628, 5.92757069, 5.92751337]);
>>> assert np.allclose(dec, [-72.01341247, -72.01341273, -72.013413])
>>> x, y = w.numerical_inverse(ra, dec)
>>> assert np.allclose(x, [1.00000005, 2.00000005, 3.00000006]);
>>> assert np.allclose(y, [1.00000004, 0.99999979, 1.00000015]);
>>> x, y = w.numerical_inverse(ra, dec, maxiter=3, tolerance=1.0e-10, quiet=False)
Traceback (most recent call last):
...
gwcs.wcs.NoConvergence: 'WCS.numerical_inverse' failed to converge to the
requested accuracy after 3 iterations.
>>> w.numerical_inverse(
...     *w([1, 300000, 3], [2, 1000000, 5], with_bounding_box=False),
...     adaptive=False,
...     detect_divergence=True,
...     quiet=False,
...     with_bounding_box=False
... )
Traceback (most recent call last):
...
gwcs.wcs.NoConvergence: 'WCS.numerical_inverse' failed to converge to the
requested accuracy. After 4 iterations, the solution is diverging at
least for one input point.
>>> # Now try to use some diverging data:
>>> divra, divdec = w([1, 300000, 3], [2, 1000000, 5], with_bounding_box=False)
>>> assert np.allclose(divra, [5.92762673, 148.21600848, 5.92750827])
>>> assert np.allclose(divdec, [-72.01339464, -7.80968079, -72.01334172])
>>> try:  
...     x, y = w.numerical_inverse(divra, divdec, maxiter=20,
...                                tolerance=1.0e-4, adaptive=True,
...                                detect_divergence=True,
...                                quiet=False)
... except NoConvergence as e:
...     print(f"Indices of diverging points: {e.divergent}")
...     print(f"Indices of poorly converging points: {e.slow_conv}")
...     print(f"Best solution:\n{e.best_solution}")
...     print(f"Achieved accuracy:\n{e.accuracy}")
Indices of diverging points: None
Indices of poorly converging points: [1]
Best solution:
[[1.00000040e+00 1.99999841e+00]
 [6.33507833e+17 3.40118820e+17]
 [3.00000038e+00 4.99999841e+00]]
Achieved accuracy:
[[2.75925982e-05 1.18471543e-05]
 [3.65405005e+04 1.31364188e+04]
 [2.76552923e-05 1.14789013e-05]]
set_transform(from_frame, to_frame, transform)[source]

Set/replace the transform between two coordinate frames.

Parameters:
from_framestr or CoordinateFrame

Initial coordinate frame.

to_framestr, or instance of CoordinateFrame

End coordinate frame.

transformModel

Transform between from_frame and to_frame.

to_fits(bounding_box=None, max_pix_error=0.25, degree=None, max_inv_pix_error=0.25, inv_degree=None, npoints=32, crpix=None, projection='TAN', bin_ext_name='WCS-TABLE', coord_col_name='coordinates', sampling=1, verbose=False)[source]

Construct a FITS WCS -TAB-based approximation to the WCS in the form of a FITS header and a binary table extension. For the description of the FITS WCS -TAB convention, see “Representations of spectral coordinates in FITS” in Greisen, E. W. et al. A&A 446 (2) 747-771 (2006) . If WCS contains celestial frame, PC/CD formalism will be used for the celestial axes.

Note

SIP distortion fitting requires that the WCS object has only two celestial axes. When WCS does not contain celestial axes, SIP fitting parameters (max_pix_error, degree, max_inv_pix_error, inv_degree, and projection) are ignored. When a WCS, in addition to celestial frame, contains other types of axes, SIP distortion fitting is disabled (ony linear terms are fitted for celestial frame).

Parameters:
bounding_boxtuple, optional

Specifies the range of acceptable values for each input axis. The order of the axes is axes_order. For two image axes bounding_box is of the form ((xmin, xmax), (ymin, ymax)).

max_pix_errorfloat, optional

Maximum allowed error over the domain of the pixel array. This error is the equivalent pixel error that corresponds to the maximum error in the output coordinate resulting from the fit based on a nominal plate scale.

degreeint, iterable, None, optional

Degree of the SIP polynomial. Default value None indicates that all allowed degree values ([1...9]) will be considered and the lowest degree that meets accuracy requerements set by max_pix_error will be returned. Alternatively, degree can be an iterable containing allowed values for the SIP polynomial degree. This option is similar to default None but it allows caller to restrict the range of allowed SIP degrees used for fitting. Finally, degree can be an integer indicating the exact SIP degree to be fit to the WCS transformation. In this case max_pixel_error is ignored.

Note

When WCS object has When degree is None and the WCS object has

max_inv_pix_errorfloat, optional

Maximum allowed inverse error over the domain of the pixel array in pixel units. If None, no inverse is generated.

inv_degreeint, iterable, None, optional

Degree of the SIP polynomial. Default value None indicates that all allowed degree values ([1...9]) will be considered and the lowest degree that meets accuracy requerements set by max_pix_error will be returned. Alternatively, degree can be an iterable containing allowed values for the SIP polynomial degree. This option is similar to default None but it allows caller to restrict the range of allowed SIP degrees used for fitting. Finally, degree can be an integer indicating the exact SIP degree to be fit to the WCS transformation. In this case max_inv_pixel_error is ignored.

npointsint, optional

The number of points in each dimension to sample the bounding box for use in the SIP fit. Minimum number of points is 3.

crpixlist of float, None, optional

Coordinates (1-based) of the reference point for the new FITS WCS. When not provided, i.e., when set to None (default) the reference pixel will be chosen near the center of the bounding box for axes corresponding to the celestial frame.

projectionstr, Pix2SkyProjection, optional

Projection to be used for the created FITS WCS. It can be specified as a string of three characters specifying a FITS projection code from Table 13 in Representations of World Coordinates in FITS (Paper I), Greisen, E. W., and Calabretta, M. R., A & A, 395, 1061-1075, 2002. Alternatively, it can be an instance of one of the astropy’s Pix2Sky_* projection models inherited from Pix2SkyProjection.

bin_ext_namestr, optional

Extension name for the BinTableHDU HDU for those axes groups that will be converted using FITW WCS’ -TAB algorith. Extension version will be determined automatically based on the number of separable group of axes.

coord_col_namestr, optional

Field name of the coordinate array in the structured array stored in BinTableHDU data. This corresponds to TTYPEi field in the FITS header of the binary table extension.

samplingfloat, tuple, optional

The target “density” of grid nodes per pixel to be used when creating the coordinate array for the -TAB FITS WCS convention. It is equal to 1/step where step is the distance between grid nodes in pixels. sampling can be specified as a single number to be used for all axes or as a tuple of numbers that specify the sampling for each image axis.

verbosebool, optional

Print progress of fits.

Returns:
hdrHeader

Header with WCS-TAB information associated (to be used) with image data.

hdulista list of BinTableHDU

A Python list of binary table extensions containing the coordinate array for TAB extensions; one extension per separable axes group.

Raises:
ValueError

When bounding_box is not defined either through the input bounding_box parameter or this object’s bounding_box property.

ValueError

When sampling is a tuple of length larger than 1 that does not match the number of image axes.

RuntimeError

If the number of image axes (~gwcs.WCS.pixel_n_dim) is larger than the number of world axes (~gwcs.WCS.world_n_dim).

to_fits_sip(bounding_box=None, max_pix_error=0.25, degree=None, max_inv_pix_error=0.25, inv_degree=None, npoints=32, crpix=None, projection='TAN', verbose=False)[source]

Construct a SIP-based approximation to the WCS for the axes corresponding to the CelestialFrame in the form of a FITS header.

The default mode in using this attempts to achieve roughly 0.25 pixel accuracy over the whole image.

Parameters:
bounding_boxtuple, optional

A pair of tuples, each consisting of two numbers Represents the range of pixel values in both dimensions ((xmin, xmax), (ymin, ymax))

max_pix_errorfloat, optional

Maximum allowed error over the domain of the pixel array. This error is the equivalent pixel error that corresponds to the maximum error in the output coordinate resulting from the fit based on a nominal plate scale. Ignored when degree is an integer or a list with a single degree.

degreeint, iterable, None, optional

Degree of the SIP polynomial. Default value None indicates that all allowed degree values ([1...9]) will be considered and the lowest degree that meets accuracy requerements set by max_pix_error will be returned. Alternatively, degree can be an iterable containing allowed values for the SIP polynomial degree. This option is similar to default None but it allows caller to restrict the range of allowed SIP degrees used for fitting. Finally, degree can be an integer indicating the exact SIP degree to be fit to the WCS transformation. In this case max_pixel_error is ignored.

max_inv_pix_errorfloat, optional

Maximum allowed inverse error over the domain of the pixel array in pixel units. If None, no inverse is generated. Ignored when degree is an integer or a list with a single degree.

inv_degreeint, iterable, None, optional

Degree of the SIP polynomial. Default value None indicates that all allowed degree values ([1...9]) will be considered and the lowest degree that meets accuracy requerements set by max_pix_error will be returned. Alternatively, degree can be an iterable containing allowed values for the SIP polynomial degree. This option is similar to default None but it allows caller to restrict the range of allowed SIP degrees used for fitting. Finally, degree can be an integer indicating the exact SIP degree to be fit to the WCS transformation. In this case max_inv_pixel_error is ignored.

npointsint, optional

The number of points in each dimension to sample the bounding box for use in the SIP fit. Minimum number of points is 3.

crpixlist of float, None, optional

Coordinates (1-based) of the reference point for the new FITS WCS. When not provided, i.e., when set to None (default) the reference pixel will be chosen near the center of the bounding box for axes corresponding to the celestial frame.

projectionstr, Pix2SkyProjection, optional

Projection to be used for the created FITS WCS. It can be specified as a string of three characters specifying a FITS projection code from Table 13 in Representations of World Coordinates in FITS (Paper I), Greisen, E. W., and Calabretta, M. R., A & A, 395, 1061-1075, 2002. Alternatively, it can be an instance of one of the astropy’s Pix2Sky_* projection models inherited from Pix2SkyProjection.

verbosebool, optional

Print progress of fits.

Returns:
FITS header with all SIP WCS keywords
Raises:
ValueError

If the WCS is not at least 2D, an exception will be raised. If the specified accuracy (both forward and inverse, both rms and maximum) is not achieved an exception will be raised.

Notes

Use of this requires a judicious choice of required accuracies. Attempts to use higher degrees (~7 or higher) will typically fail due to floating point problems that arise with high powers.

to_fits_tab(bounding_box=None, bin_ext_name='WCS-TABLE', coord_col_name='coordinates', sampling=1)[source]

Construct a FITS WCS -TAB-based approximation to the WCS in the form of a FITS header and a binary table extension. For the description of the FITS WCS -TAB convention, see “Representations of spectral coordinates in FITS” in Greisen, E. W. et al. A&A 446 (2) 747-771 (2006) .

Parameters:
bounding_boxtuple, optional

Specifies the range of acceptable values for each input axis. The order of the axes is axes_order. For two image axes bounding_box is of the form ((xmin, xmax), (ymin, ymax)).

bin_ext_namestr, optional

Extension name for the BinTableHDU HDU for those axes groups that will be converted using FITW WCS’ -TAB algorith. Extension version will be determined automatically based on the number of separable group of axes.

coord_col_namestr, optional

Field name of the coordinate array in the structured array stored in BinTableHDU data. This corresponds to TTYPEi field in the FITS header of the binary table extension.

samplingfloat, tuple, optional

The target “density” of grid nodes per pixel to be used when creating the coordinate array for the -TAB FITS WCS convention. It is equal to 1/step where step is the distance between grid nodes in pixels. sampling can be specified as a single number to be used for all axes or as a tuple of numbers that specify the sampling for each image axis.

Returns:
hdrHeader

Header with WCS-TAB information associated (to be used) with image data.

bin_table_hduBinTableHDU

Binary table extension containing the coordinate array.

Raises:
ValueError

When bounding_box is not defined either through the input bounding_box parameter or this object’s bounding_box property.

ValueError

When sampling is a tuple of length larger than 1 that does not match the number of image axes.

RuntimeError

If the number of image axes (~gwcs.WCS.pixel_n_dim) is larger than the number of world axes (~gwcs.WCS.world_n_dim).

transform(from_frame, to_frame, *args, **kwargs)[source]

Transform positions between two frames.

Parameters:
from_framestr or CoordinateFrame

Initial coordinate frame.

to_framestr, or instance of CoordinateFrame

Coordinate frame into which to transform.

argsfloat or array-like

Inputs in from_frame, separate inputs for each dimension.

output_with_unitsbool

If True - returns a SkyCoord or SpectralCoord object.

with_bounding_boxbool, optional

If True(default) values in the result which correspond to any of the inputs being outside the bounding_box are set to fill_value.

fill_valuefloat, optional

Output value for inputs outside the bounding_box (default is np.nan).