NDCube¶
- class ndcube.NDCube(data, wcs=None, uncertainty=None, mask=None, meta=None, unit=None, copy=False, **kwargs)[source]¶
Bases:
NDCubeBaseClass representing N-D data described by a single array and set of WCS transformations.
- Parameters:
data (
numpy.ndarray) – The array holding the actual data in this object.wcs (
astropy.wcs.wcsapi.BaseLowLevelWCS,astropy.wcs.wcsapi.BaseHighLevelWCS, optional) – The WCS object containing the axes’ information, optional only ifdatais anastropy.nddata.NDDataobject.uncertainty (any type, optional) – Uncertainty in the dataset. Should have an attribute uncertainty_type that defines what kind of uncertainty is stored, for example “std” for standard deviation or “var” for variance. A metaclass defining such an interface is NDUncertainty - but isn’t mandatory. If the uncertainty has no such attribute the uncertainty is stored as UnknownUncertainty. Defaults to None.
mask (any type, optional) – Mask for the dataset. Masks should follow the numpy convention that valid data points are marked by False and invalid ones with True. Defaults to None.
meta (dict-like object, optional) – Additional meta information about the dataset. If no meta is provided an empty collections.OrderedDict is created. Default is None.
unit (Unit-like or str, optional) – Unit for the dataset. Strings that can be converted to a Unit are allowed. Default is None.
extra_coords (iterable of
tuple, each with three entries) – (str,int,astropy.units.quantityor array-like) Gives the name, axis of data, and values of coordinates of a data axis not included in the WCS object.copy (bool, optional) – Indicates whether to save the arguments as copy. True copies every attribute before saving it while False tries to save every parameter as reference. Note however that it is not always possible to save the input as reference. Default is False.
Attributes Summary
Returns the physical types associated with each array axis.
A
BaseHighLevelWCSobject which combines.wcswith.extra_coords.The stored dataset.
The array dimensions of the cube.
An
ExtraCoordsobject holding extra coordinates aligned to array axes.A
GlobalCoordsobject holding coordinate metadata not aligned to an array axis.Mask for the dataset, if any.
Additional meta information about the dataset.
A
MatplotlibPlotterinstance providing visualization methods.Image representation of the PSF for the dataset.
Uncertainty in the dataset, if any.
Unit for the dataset, if any.
A world coordinate system (WCS) for the dataset, if any.
Methods Summary
axis_world_coords(*axes[, pixel_corners, wcs])Returns WCS coordinate values of all pixels for all axes.
axis_world_coords_values(*axes[, ...])Returns WCS coordinate values of all pixels for desired axes.
crop(*points[, wcs])Crop to the smallest cube in pixel space containing the world coordinate points.
crop_by_values(*points[, units, wcs])Crop to the smallest cube in pixel space containing the world coordinate points.
explode_along_axis(axis)Separates slices of NDCubes along a given axis into an NDCubeSequence of (N-1)DCubes.
plot([axes, plot_axes, axes_coordinates, ...])Visualize the
NDCube.reproject_to(target_wcs[, algorithm, ...])Reprojects this
NDCubeto the coordinates described by another WCS object.Attributes Documentation
- array_axis_physical_types¶
Returns the physical types associated with each array axis.
Returns an iterable of tuples where each tuple corresponds to an array axis and holds strings denoting the physical types associated with that array axis. Since multiple physical types can be associated with one array axis, tuples can be of different lengths. Likewise, as a single physical type can correspond to multiple array axes, the same physical type string can appear in multiple tuples.
The physical types are drawn from the WCS ExtraCoords objects.
- combined_wcs¶
A
BaseHighLevelWCSobject which combines.wcswith.extra_coords.
- dimensions¶
- extra_coords¶
An
ExtraCoordsobject holding extra coordinates aligned to array axes.
- global_coords¶
A
GlobalCoordsobject holding coordinate metadata not aligned to an array axis.
- mask¶
Mask for the dataset, if any.
Masks should follow the
numpyconvention that valid data points are marked byFalseand invalid ones withTrue.- Type:
any type
- plotter = None¶
A
MatplotlibPlotterinstance providing visualization methods.The type of this attribute can be changed to provide custom visualization functionality.
- psf¶
- uncertainty¶
Uncertainty in the dataset, if any.
Should have an attribute
uncertainty_typethat defines what kind of uncertainty is stored, such as'std'for standard deviation or'var'for variance. A metaclass defining such an interface isNDUncertaintybut isn’t mandatory.- Type:
any type
- unit¶
Unit for the dataset, if any.
- Type:
Unit
- wcs¶
A world coordinate system (WCS) for the dataset, if any.
- Type:
any type
Methods Documentation
- axis_world_coords(*axes, pixel_corners=False, wcs=None)¶
Returns WCS coordinate values of all pixels for all axes.
- Parameters:
axes (
intorstr, or multipleintorstr, optional) – Axis number in numpy ordering or unique substring ofworld_axis_physical_typesof axes for which real world coordinates are desired. axes=None implies all axes will be returned.pixel_corners (
bool, optional) – IfTruethen instead of returning the coordinates at the centers of the pixels, the coordinates at the pixel corners will be returned. This increases the size of the output by 1 in all dimensions as all corners are returned.wcs (
astropy.wcs.wcsapi.BaseHighLevelWCS, optional) – The WCS object to used to calculate the world coordinates. Although technically this can be any valid WCS, it will typically beself.wcs,self.extra_coords, orself.combined_wcswhich combines both the WCS and extra coords. Defaults to the.wcsproperty.
- Returns:
axes_coords – An iterable of “high level” objects giving the real world coords for the axes requested by user. For example, a tuple of
SkyCoordobjects. The types returned are determined by the WCS object. The dimensionality of these objects should match that of their corresponding array dimensions, unlesspixel_corners=Truein which case the length along each axis will be 1 greater than the number of pixels.- Return type:
Example
>>> NDCube.all_world_coords(('lat', 'lon')) >>> NDCube.all_world_coords(2)
- axis_world_coords_values(*axes, pixel_corners=False, wcs=None)¶
Returns WCS coordinate values of all pixels for desired axes.
- Parameters:
axes (
intorstr, or multipleintorstr, optional) – Axis number in numpy ordering or unique substring ofworld_axis_physical_typesof axes for which real world coordinates are desired. axes=None implies all axes will be returned.pixel_corners (
bool, optional) – IfTruethen instead of returning the coordinates of the pixel centers the coordinates of the pixel corners will be returned. This increases the size of the output along each dimension by 1 as all corners are returned.wcs (
astropy.wcs.wcsapi.BaseHighLevelWCS, optional) – The WCS object to used to calculate the world coordinates. Although technically this can be any valid WCS, it will typically beself.wcs,self.extra_coords, orself.combined_wcs, combing both the WCS and extra coords. Defaults to the.wcsproperty.
- Returns:
axes_coords – An iterable of “high level” objects giving the real world coords for the axes requested by user. For example, a tuple of
SkyCoordobjects. The types returned are determined by the WCS object. The dimensionality of these objects should match that of their corresponding array dimensions, unlesspixel_corners=Truein which case the length along each axis will be 1 greater than the number of pixels.- Return type:
Example
>>> NDCube.all_world_coords_values(('lat', 'lon')) >>> NDCube.all_world_coords_values(2)
- crop(*points, wcs=None)¶
Crop to the smallest cube in pixel space containing the world coordinate points.
- Parameters:
points (iterable of iterables) –
Tuples of high level coordinate objects e.g.
SkyCoord. The coordinates of the points must be specified in Cartesian (WCS) order as they are passed toworld_to_array_index. Therefore their number and order must be compatible with the API of that method.It is possible to not specify a coordinate for an axis by replacing any object with
None. Any coordinate replaced byNonewill not be used to calculate pixel coordinates, and therefore not affect the calculation of the final bounding box.wcs (
astropy.wcs.wcsapi.BaseLowLevelWCS) – The WCS to use to calculate the pixel coordinates based on the input. Will default to the.wcsproperty if not given. While any valid WCS could be used it is expected that either the.wcs,.combined_wcs, or.extra_coordsproperties will be used.
- Returns:
result
- Return type:
- crop_by_values(*points, units=None, wcs=None)¶
Crop to the smallest cube in pixel space containing the world coordinate points.
- Parameters:
points (iterable of iterables) –
Tuples of coordinates as
Quantityobjects. The coordinates of the points must be specified in Cartesian (WCS) order as they are passed toworld_to_array_index_values. Therefore their number and order must be compatible with the API of that method.It is possible to not specify a coordinate for an axis by replacing any coordinate with
None. Any coordinate replaced byNonewill not be used to calculate pixel coordinates, and therefore not affect the calculation of the final bounding box. Note that you must specify either none or all coordinates for any correlated axes, e.g. both spatial coordinates.units (iterable of
astropy.units.Unit) – The unit of the corresponding entries in each point. Must therefore be the same length as the number of world axes. Only used if the corresponding type is not aastropy.units.QuantityorNone.wcs (
astropy.wcs.wcsapi.BaseLowLevelWCS) – The WCS to use to calculate the pixel coordinates based on the input. Will default to the.wcsproperty if not given. While any valid WCS could be used it is expected that either the.wcs,.combined_wcs, or.extra_coordsproperties will be used.
- Returns:
result
- Return type:
- explode_along_axis(axis)¶
Separates slices of NDCubes along a given axis into an NDCubeSequence of (N-1)DCubes.
- Parameters:
axis (
int) – The array axis along which the data is to be changed.- Returns:
result
- Return type:
- plot(axes=None, plot_axes=None, axes_coordinates=None, axes_units=None, data_unit=None, wcs=None, **kwargs)¶
Visualize the
NDCube.- Parameters:
axes (
WCSAxesor None:, optional) – The axes to plot onto. If None the current axes will be used.plot_axes (
list, optional) – A list of length equal to the number of pixel dimensions in array axis order. This list selects which cube axes are displayed on which plot axes. For an image plot this list should contain'x'and'y'for the plot axes andNonefor all the other elements. For a line plot it should only contain'x'andNonefor all the other elements.axes_unit (
list, optional) – A list of length equal to the number of world dimensions specifying the units of each axis, orNoneto use the default unit for that axis.axes_coordinates (
list, optional) – A list of length equal to the number of pixel dimensions. For each axis the value of the list should either be a string giving the world axis type orNoneto use the default axis from the WCS.data_unit (
astropy.unit.Unit) – The data is changed to the unit given or theNDCube.unitif not given.wcs (
astropy.wcs.wcsapi.BaseHighLevelWCS) – The WCS object to define the coordinates of the plot axes.kwargs – Additional keyword arguments are given to the underlying plotting infrastructure which depends on the dimensionality of the data and whether 1 or 2 plot_axes are defined: - Animations:
mpl_animators.ArrayAnimatorWCS- Static 2-D images:matplotllib.pyplot.imshow- Static 1-D line plots:matplotllib.pyplot.plot
- reproject_to(target_wcs, algorithm='interpolation', shape_out=None, return_footprint=False, **reproject_args)¶
Reprojects this
NDCubeto the coordinates described by another WCS object.- Parameters:
target_wcs (
astropy.wcs.wcsapi.BaseHighLevelWCS,astropy.wcs.wcsapi.BaseLowLevelWCS,) – orastropy.io.fits.HeaderThe WCS object to which theNDCubeis to be reprojected.algorithm (
str{‘interpolation’, ‘adaptive’, ‘exact’}) – The algorithm to use for reprojecting. When set to ‘interpolation’ ~`reproject.reproject_interp` is used, when set to ‘adaptive’ ~`reproject.reproject_adaptive` is used and when set to ‘exact’ ~`reproject.reproject_exact` is used.shape_out (
tuple, optional) – The shape of the output data array. The ordering of the dimensions must follow NumPy ordering and not the WCS pixel shape. If not specified,array_shapeattribute (if available) from the low level API of thetarget_wcsis used.return_footprint (
bool) – IfTrue`the footprint is returned in addition to the newNDCube. Defaults toFalse.**reproject_args – All other arguments are passed through to the reproject function being called. The function being called depends on the
algorithm=keyword argument, see that for more details.
- Returns:
resampled_cube (
ndcube.NDCube) – A new resultant NDCube object, the suppliedtarget_wcswill be the.wcsattribute of the outputNDCube.footprint (
numpy.ndarray) – Footprint of the input array in the output array. Values of 0 indicate no coverage or valid values in the input image, while values of 1 indicate valid values.
See also
reproject.reproject_interpreproject.reproject_adaptivereproject.reproject_exact
Notes
This method doesn’t support handling of the
mask,extra_coords, anduncertaintyattributes yet. However,metaandglobal_coordsare copied to the outputndcube.NDCube.
