SpectralCube

class spectral_cube.SpectralCube(data, wcs, mask=None, meta=None, fill_value=nan, header=None, allow_huge_operations=False, beam=None, read_beam=True, wcs_tolerance=0.0, **kwargs)[source]

Bases: spectral_cube.spectral_cube.BaseSpectralCube

Attributes Summary

base The data type ‘base’ of the cube - useful for, e.g., joblib
fill_value The replacement value used by filled_data().
filled_data Return a portion of the data array, with excluded mask values replaced by fill_value.
hdu HDU version of self
hdulist
header
latitude_extrema
longitude_extrema
mask
meta
ndim Dimensionality of the data
shape Length of cube along each axis
size Number of elements in the cube
spatial_coordinate_map
spectral_axis A Quantity array containing the central values of each channel along the spectral axis.
spectral_extrema
unit The flux unit
unitless Return a copy of self with unit set to None
unmasked_data Return a view of the subset of the underlying data, ignoring the mask.
velocity_convention The equivalencies that describes the spectral axis
wcs
world Return a list of the world coordinates in a cube, projection, or a view of it.
world_extrema

Methods Summary

apply_function(function[, axis, weights, ...]) Apply a function to valid data along the specified axis or to the whole
apply_numpy_function(function[, fill, ...]) Apply a numpy function to the cube
argmax(*args, **kwargs) Return the index of the maximum data value.
argmin(*args, **kwargs) Return the index of the minimum data value.
chunked([chunksize]) Not Implemented.
closest_spectral_channel(value) Find the index of the closest spectral channel to the specified spectral coordinate.
convolve_to(beam[, convolve]) Convolve each channel in the cube to a specified beam
find_lines([velocity_offset, ...]) Using astroquery’s splatalogue interface, search for lines within the spectral band.
flattened([slice, weights]) Return a slice of the cube giving only the valid data (i.e., removing
flattened_world([view]) Retrieve the world coordinates corresponding to the extracted flattened
get_mask_array() Convert the mask to a boolean numpy array
linewidth_fwhm([how]) Compute a (FWHM) linewidth map along the spectral axis.
linewidth_sigma([how]) Compute a (sigma) linewidth map along the spectral axis.
max(*args, **kwargs) Return the maximum data value of the cube, optionally over an axis.
mean(*args, **kwargs) Return the mean of the cube, optionally over an axis.
median([axis, iterate_rays]) Compute the median of an array, optionally along an axis.
min(*args, **kwargs) Return the minimum data value of the cube, optionally over an axis.
minimal_subcube([spatial_only]) Return the minimum enclosing subcube where the mask is valid
moment([order, axis, how]) Compute moments along the spectral axis.
moment0([axis, how]) Compute the zeroth moment along an axis.
moment1([axis, how]) Compute the 1st moment along an axis.
moment2([axis, how]) Compute the 2nd moment along an axis.
percentile(q[, axis, iterate_rays]) Return percentiles of the data.
read(filename[, format, hdu]) Read a spectral cube from a file.
reproject(header[, order]) Reproject the cube into a new header.
spectral_interpolate(spectral_grid[, ...]) Resample the cube spectrally onto a specific grid
spectral_slab(lo, hi) Extract a new cube between two spectral coordinates
spectral_smooth(kernel[, convolve]) Smooth the cube along the spectral dimension
std(*args, **kwargs) Return the standard deviation of the cube, optionally over an axis.
subcube([xlo, xhi, ylo, yhi, zlo, zhi, ...]) Extract a sub-cube spatially and spectrally.
subcube_from_ds9region(ds9region[, allow_empty]) Extract a masked subcube from a ds9 region or a pyregion Region object
subcube_from_mask(region_mask) Given a mask, return the minimal subcube that encloses the mask
subcube_slices_from_mask(region_mask[, ...]) Given a mask, return the slices corresponding to the minimum subcube
sum(*args, **kwargs) Return the sum of the cube, optionally over an axis.
to(*args, **kwargs) Return the cube converted to the given unit (assuming it is equivalent).
to_ds9([ds9id, newframe]) Send the data to ds9 (this will create a copy in memory)
to_glue([name, glue_app, dataset, start_gui]) Send data to a new or existing Glue application
to_pvextractor() Open the cube in a quick viewer written in matplotlib that allows you
to_yt([spectral_factor, nprocs]) Convert a spectral cube to a yt object that can be further analyzed in yt.
unmasked_copy() Return a copy of the cube with no mask (i.e., all data included)
with_fill_value(fill_value) Create a new SpectralCube with a different fill_value.
with_mask(mask[, inherit_mask, wcs_tolerance]) Return a new SpectralCube instance that contains a composite mask of the current SpectralCube and the new mask.
with_spectral_unit(unit[, ...]) Returns a new Cube with a different Spectral Axis unit
world_spines() Returns a list of 1D arrays, for the world coordinates along each pixel axis.
write(filename[, overwrite, format]) Write the spectral cube to a file.

Attributes Documentation

base

The data type ‘base’ of the cube - useful for, e.g., joblib

fill_value

The replacement value used by filled_data().

fill_value is immutable; use with_fill_value() to create a new cube with a different fill value.

filled_data
Return a portion of the data array, with excluded mask values replaced by fill_value.
Returns:

data : Quantity

The masked data.

Notes

Supports efficient Numpy slice notation, like filled_data[0:3, :, 2:4]

hdu

HDU version of self

hdulist
header
latitude_extrema
longitude_extrema
mask
meta
ndim

Dimensionality of the data

shape

Length of cube along each axis

size

Number of elements in the cube

spatial_coordinate_map
spectral_axis

A Quantity array containing the central values of each channel along the spectral axis.

spectral_extrema
unit

The flux unit

unitless

Return a copy of self with unit set to None

unmasked_data
Return a view of the subset of the underlying data, ignoring the mask.
Returns:

data : Quantity instance

The unmasked data

Notes

Supports efficient Numpy slice notation, like unmasked_data[0:3, :, 2:4]

velocity_convention

The equivalencies that describes the spectral axis

wcs
world

Return a list of the world coordinates in a cube, projection, or a view of it.

SpatialCoordMixinClass.world is called with bracket notation, like a NumPy array:

c.world[0:3, :, :]
Returns:

[v, y, x] : list of NumPy arrays

The 3 world coordinates at each pixel in the view. For a 2D image, the output is [y, x].

Notes

Supports efficient Numpy slice notation, like world[0:3, :, 2:4]

Examples

Extract the first 3 velocity channels of the cube:

>>> v, y, x = c.world[0:3]

Extract all the world coordinates:

>>> v, y, x = c.world[:, :, :]

Extract every other pixel along all axes:

>>> v, y, x = c.world[::2, ::2, ::2]

Extract all the world coordinates for a 2D image:

>>> y, x = c.world[:, :]
world_extrema

Methods Documentation

apply_function(function, axis=None, weights=None, unit=None, projection=False, progressbar=False, **kwargs)

Apply a function to valid data along the specified axis or to the whole cube, optionally using a weight array that is the same shape (or at least can be sliced in the same way)

Parameters:

function : function

A function that can be applied to a numpy array. Does not need to be nan-aware

axis : 1, 2, 3, or None

The axis to operate along. If None, the return is scalar.

weights : (optional) np.ndarray

An array with the same shape (or slicing abilities/results) as the data cube

unit : (optional) Unit

The unit of the output projection or value. Not all functions should return quantities with units.

projection : bool

Return a projection if the resulting array is 2D?

progressbar : bool

Show a progressbar while iterating over the slices/rays through the cube?

Returns:

result : Projection or Quantity or float

The result depends on the value of axis, projection, and unit. If axis is None, the return will be a scalar with or without units. If axis is an integer, the return will be a Projection if projection is set

apply_numpy_function(function, fill=nan, reduce=True, how='auto', projection=False, unit=None, check_endian=False, progressbar=False, **kwargs)

Apply a numpy function to the cube

Parameters:

function : Numpy ufunc

A numpy ufunc to apply to the cube

fill : float

The fill value to use on the data

reduce : bool

reduce indicates whether this is a reduce-like operation, that can be accumulated one slice at a time. sum/max/min are like this. argmax/argmin/stddev are not

how : cube | slice | ray | auto

How to compute the moment. All strategies give the same result, but certain strategies are more efficient depending on data size and layout. Cube/slice/ray iterate over decreasing subsets of the data, to conserve memory. Default=’auto’

projection : bool

Return a Projection if the resulting array is 2D or a OneDProjection if the resulting array is 1D and the sum is over both spatial axes?

unit : None or astropy.units.Unit

The unit to include for the output array. For example, SpectralCube.max calls SpectralCube.apply_numpy_function(np.max, unit=self.unit), inheriting the unit from the original cube. However, for other numpy functions, e.g. numpy.argmax, the return is an index and therefore unitless.

check_endian : bool

A flag to check the endianness of the data before applying the function. This is only needed for optimized functions, e.g. those in the bottleneck package.

progressbar : bool

Show a progressbar while iterating over the slices through the cube?

kwargs : dict

Passed to the numpy function.

Returns:

result : Projection or Quantity or float

The result depends on the value of axis, projection, and unit. If axis is None, the return will be a scalar with or without units. If axis is an integer, the return will be a Projection if projection is set

argmax(*args, **kwargs)

Return the index of the maximum data value.

The return value is arbitrary if all pixels along axis are excluded from the mask.

Ignores excluded mask elements.

Parameters:

axis : int (optional)

The axis to collapse, or None to perform a global aggregation

how : cube | slice | ray | auto

How to compute the aggregation. All strategies give the same result, but certain strategies are more efficient depending on data size and layout. Cube/slice/ray iterate over decreasing subsets of the data, to conserve memory. Default=’auto’

argmin(*args, **kwargs)

Return the index of the minimum data value.

The return value is arbitrary if all pixels along axis are excluded from the mask

Ignores excluded mask elements.

Parameters:

axis : int (optional)

The axis to collapse, or None to perform a global aggregation

how : cube | slice | ray | auto

How to compute the aggregation. All strategies give the same result, but certain strategies are more efficient depending on data size and layout. Cube/slice/ray iterate over decreasing subsets of the data, to conserve memory. Default=’auto’

chunked(chunksize=1000)

Not Implemented.

Iterate over chunks of valid data

closest_spectral_channel(value)

Find the index of the closest spectral channel to the specified spectral coordinate.

Parameters:

value : Quantity

The value of the spectral coordinate to search for.

convolve_to(beam, convolve=<function convolve_fft>)[source]

Convolve each channel in the cube to a specified beam

Parameters:

beam : radio_beam.Beam

The beam to convolve to

convolve : function

The astropy convolution function to use, either astropy.convolution.convolve or astropy.convolution.convolve_fft

Returns:

cube : SpectralCube

A SpectralCube with a single beam

find_lines(velocity_offset=None, velocity_convention=None, rest_value=None, **kwargs)

Using astroquery’s splatalogue interface, search for lines within the spectral band. See astroquery.splatalogue.Splatalogue for information on keyword arguments

Parameters:

velocity_offset : u.km/u.s equivalent

An offset by which the spectral axis should be shifted before searching splatalogue. This value will be added to the velocity, so if you want to redshift a spectrum, make this value positive, and if you want to un-redshift it, make this value negative.

velocity_convention : ‘radio’, ‘optical’, ‘relativistic’

The doppler convention to pass to with_spectral_unit

rest_value : u.GHz equivalent

The rest frequency (or wavelength or energy) to be passed to with_spectral_unit

flattened(slice=(), weights=None)

Return a slice of the cube giving only the valid data (i.e., removing bad values)

Parameters:

slice: 3-tuple

A length-3 tuple of view (or any equivalent valid slice of a cube)

weights: (optional) np.ndarray

An array with the same shape (or slicing abilities/results) as the data cube

flattened_world(view=())

Retrieve the world coordinates corresponding to the extracted flattened version of the cube

get_mask_array()

Convert the mask to a boolean numpy array

linewidth_fwhm(how='auto')

Compute a (FWHM) linewidth map along the spectral axis.

For an explanation of the how parameter, see moment().

linewidth_sigma(how='auto')

Compute a (sigma) linewidth map along the spectral axis.

For an explanation of the how parameter, see moment().

max(*args, **kwargs)

Return the maximum data value of the cube, optionally over an axis.

Ignores excluded mask elements.

Parameters:

axis : int (optional)

The axis to collapse, or None to perform a global aggregation

how : cube | slice | ray | auto

How to compute the aggregation. All strategies give the same result, but certain strategies are more efficient depending on data size and layout. Cube/slice/ray iterate over decreasing subsets of the data, to conserve memory. Default=’auto’

mean(*args, **kwargs)

Return the mean of the cube, optionally over an axis.

Ignores excluded mask elements.

Parameters:

axis : int (optional)

The axis to collapse, or None to perform a global aggregation

how : cube | slice | ray | auto

How to compute the aggregation. All strategies give the same result, but certain strategies are more efficient depending on data size and layout. Cube/slice/ray iterate over decreasing subsets of the data, to conserve memory. Default=’auto’

median(axis=None, iterate_rays=False, **kwargs)

Compute the median of an array, optionally along an axis.

Ignores excluded mask elements.

Parameters:

axis : int (optional)

The axis to collapse

iterate_rays : bool

Iterate over individual rays? This mode is slower but can save RAM costs, which may be extreme for large cubes

Returns:

med : ndarray

The median

min(*args, **kwargs)

Return the minimum data value of the cube, optionally over an axis.

Ignores excluded mask elements.

Parameters:

axis : int (optional)

The axis to collapse, or None to perform a global aggregation

how : cube | slice | ray | auto

How to compute the aggregation. All strategies give the same result, but certain strategies are more efficient depending on data size and layout. Cube/slice/ray iterate over decreasing subsets of the data, to conserve memory. Default=’auto’

minimal_subcube(spatial_only=False)

Return the minimum enclosing subcube where the mask is valid

Parameters:

spatial_only: bool

Only compute the minimal subcube in the spatial dimensions

moment(order=0, axis=0, how='auto')

Compute moments along the spectral axis.

Moments are defined as follows:

Moment 0:

\[M_0 \int I dl\]

Moment 1:

\[M_1 = \frac{\int I l dl}{M_0}\]

Moment N:

\[M_N = \frac{\int I (l - M_1)^N dl}{M_0}\]

Warning

Note that these follow the mathematical definitions of moments, and therefore the second moment will return a variance map. To get linewidth maps, you can instead use the linewidth_fwhm() or linewidth_sigma() methods.

Parameters:

order : int

The order of the moment to take. Default=0

axis : int

The axis along which to compute the moment. Default=0

how : cube | slice | ray | auto

How to compute the moment. All strategies give the same result, but certain strategies are more efficient depending on data size and layout. Cube/slice/ray iterate over decreasing subsets of the data, to conserve memory. Default=’auto’

Returns:

map [, wcs]

The moment map (numpy array) and, if wcs=True, the WCS object describing the map

Notes

Generally, how=’cube’ is fastest for small cubes that easily fit into memory. how=’slice’ is best for most larger datasets. how=’ray’ is probably only a good idea for very large cubes whose data are contiguous over the axis of the moment map.

For the first moment, the result for axis=1, 2 is the angular offset relative to the cube face. For axis=0, it is the absolute velocity/frequency of the first moment.

moment0(axis=0, how='auto')

Compute the zeroth moment along an axis.

See moment().

moment1(axis=0, how='auto')

Compute the 1st moment along an axis.

For an explanation of the axis and how parameters, see moment().

moment2(axis=0, how='auto')

Compute the 2nd moment along an axis.

For an explanation of the axis and how parameters, see moment().

percentile(q, axis=None, iterate_rays=False, **kwargs)

Return percentiles of the data.

Parameters:

q : float

The percentile to compute

axis : int, or None

Which axis to compute percentiles over

iterate_rays : bool

Iterate over individual rays? This mode is slower but can save RAM costs, which may be extreme for large cubes

read(filename, format=None, hdu=None, **kwargs)

Read a spectral cube from a file.

If the file contains Stokes axes, they will automatically be dropped. If you want to read in all Stokes informtion, use read() instead.

Parameters:

filename : str

The file to read the cube from

format : str

The format of the file to read. (Currently limited to ‘fits’ and ‘casa_image’)

hdu : int or str

For FITS files, the HDU to read in (can be the ID or name of an HDU).

kwargs : dict

If the format is ‘fits’, the kwargs are passed to open().

reproject(header, order='bilinear')

Reproject the cube into a new header. Fills the data with the cube’s fill_value to replace bad values before reprojection.

Parameters:

header : astropy.io.fits.Header

A header specifying a cube in valid WCS

order : int or str, optional

The order of the interpolation (if mode is set to 'interpolation'). This can be either one of the following strings:

  • ‘nearest-neighbor’
  • ‘bilinear’
  • ‘biquadratic’
  • ‘bicubic’

or an integer. A value of 0 indicates nearest neighbor interpolation.

spectral_interpolate(spectral_grid, suppress_smooth_warning=False, fill_value=None)[source]

Resample the cube spectrally onto a specific grid

Parameters:

spectral_grid : array

An array of the spectral positions to regrid onto

suppress_smooth_warning : bool

If disabled, a warning will be raised when interpolating onto a grid that does not nyquist sample the existing grid. Disable this if you have already appropriately smoothed the data.

fill_value : float

Value for extrapolated spectral values that lie outside of the spectral range defined in the original data. The default is to use the nearest spectral channel in the cube.

Returns:

cube : SpectralCube

spectral_slab(lo, hi)

Extract a new cube between two spectral coordinates

Parameters:

lo, hi : Quantity

The lower and upper spectral coordinate for the slab range. The units should be compatible with the units of the spectral axis. If the spectral axis is in frequency-equivalent units and you want to select a range in velocity, or vice-versa, you should first use with_spectral_unit() to convert the units of the spectral axis.

spectral_smooth(kernel, convolve=<function convolve>, **kwargs)[source]

Smooth the cube along the spectral dimension

Parameters:

kernel : Kernel1D

A 1D kernel from astropy

convolve : function

The astropy convolution function to use, either astropy.convolution.convolve or astropy.convolution.convolve_fft

kwargs : dict

Passed to the convolve function

std(*args, **kwargs)

Return the standard deviation of the cube, optionally over an axis.

Ignores excluded mask elements.

Parameters:

axis : int (optional)

The axis to collapse, or None to perform a global aggregation

how : cube | slice | ray | auto

How to compute the aggregation. All strategies give the same result, but certain strategies are more efficient depending on data size and layout. Cube/slice/ray iterate over decreasing subsets of the data, to conserve memory. Default=’auto’

subcube(xlo='min', xhi='max', ylo='min', yhi='max', zlo='min', zhi='max', rest_value=None)

Extract a sub-cube spatially and spectrally.

Parameters:

[xyz]lo/[xyz]hi : int or Quantity or min/max

The endpoints to extract. If given as a quantity, will be interpreted as World coordinates. If given as a string or int, will be interpreted as pixel coordinates.

subcube_from_ds9region(ds9region, allow_empty=False)

Extract a masked subcube from a ds9 region or a pyregion Region object (only functions on celestial dimensions)

Parameters:

ds9region: str or `pyregion.Shape`

The region to extract

allow_empty: bool

If this is False, an exception will be raised if the region contains no overlap with the cube

subcube_from_mask(region_mask)

Given a mask, return the minimal subcube that encloses the mask

Parameters:

region_mask: `masks.MaskBase` or boolean `numpy.ndarray`

The mask with appropraite WCS or an ndarray with matched coordinates

subcube_slices_from_mask(region_mask, spatial_only=False)

Given a mask, return the slices corresponding to the minimum subcube that encloses the mask

Parameters:

region_mask: `masks.MaskBase` or boolean `numpy.ndarray`

The mask with appropriate WCS or an ndarray with matched coordinates

spatial_only: bool

Return only slices that affect the spatial dimensions; the spectral dimension will be left unchanged

sum(*args, **kwargs)

Return the sum of the cube, optionally over an axis.

Ignores excluded mask elements.

Parameters:

axis : int (optional)

The axis to collapse, or None to perform a global aggregation

how : cube | slice | ray | auto

How to compute the aggregation. All strategies give the same result, but certain strategies are more efficient depending on data size and layout. Cube/slice/ray iterate over decreasing subsets of the data, to conserve memory. Default=’auto’

to(*args, **kwargs)

Return the cube converted to the given unit (assuming it is equivalent). If conversion was required, this will be a copy, otherwise it will

to_ds9(ds9id=None, newframe=False)

Send the data to ds9 (this will create a copy in memory)

Parameters:

ds9id: None or string

The DS9 session ID. If ‘None’, a new one will be created. To find your ds9 session ID, open the ds9 menu option File:XPA:Information and look for the XPA_METHOD string, e.g. XPA_METHOD:  86ab2314:60063. You would then calll this function as cube.to_ds9('86ab2314:60063')

newframe: bool

Send the cube to a new frame or to the current frame?

to_glue(name=None, glue_app=None, dataset=None, start_gui=True)

Send data to a new or existing Glue application

Parameters:

name : str or None

The name of the dataset within Glue. If None, defaults to ‘SpectralCube’. If a dataset with the given name already exists, a new dataset with “_” appended will be added instead.

glue_app : GlueApplication or None

A glue application to send the data to. If this is not specified, a new glue application will be started if one does not already exist for this cube. Otherwise, the data will be sent to the existing glue application, self._glue_app.

dataset : glue.core.Data or None

An existing Data object to add the cube to. This is a good way to compare cubes with the same dimensions. Supercedes glue_app

start_gui : bool

Start the GUI when this is run. Set to False for testing.

to_pvextractor()

Open the cube in a quick viewer written in matplotlib that allows you to create PV extractions within the GUI

to_yt(spectral_factor=1.0, nprocs=None, **kwargs)

Convert a spectral cube to a yt object that can be further analyzed in yt.

Parameters:

spectral_factor : float, optional

Factor by which to stretch the spectral axis. If set to 1, one pixel in spectral coordinates is equivalent to one pixel in spatial coordinates.

If using yt 3.0 or later, additional keyword arguments will be passed

onto yt’s ``FITSDataset`` constructor. See the yt documentation

(http://yt-project.org/docs/3.0/examining/loading_data.html?#fits-data)

for details on options for reading FITS data.

unmasked_copy()

Return a copy of the cube with no mask (i.e., all data included)

with_fill_value(fill_value)

Create a new SpectralCube with a different fill_value.

Notes

This method is fast (it does not copy any data)

with_mask(mask, inherit_mask=True, wcs_tolerance=None)

Return a new SpectralCube instance that contains a composite mask of the current SpectralCube and the new mask. Values of the mask that are True will be included (masks are analogous to numpy boolean index arrays, they are the inverse of the .mask attribute of a numpy masked array).

Parameters:

mask : MaskBase instance, or boolean numpy array

The mask to apply. If a boolean array is supplied, it will be converted into a mask, assuming that True values indicate included elements.

inherit_mask : bool (optional, default=True)

If True, combines the provided mask with the mask currently attached to the cube

wcs_tolerance : None or float

The tolerance of difference in WCS parameters between the cube and the mask. Defaults to self._wcs_tolerance (which itself defaults to 0.0) if unspecified

Returns:

new_cube : SpectralCube

A cube with the new mask applied.

Notes

This operation returns a view into the data, and not a copy.

with_spectral_unit(unit, velocity_convention=None, rest_value=None)

Returns a new Cube with a different Spectral Axis unit

Parameters:

unit : Unit

Any valid spectral unit: velocity, (wave)length, or frequency. Only vacuum units are supported.

velocity_convention : ‘relativistic’, ‘radio’, or ‘optical’

The velocity convention to use for the output velocity axis. Required if the output type is velocity. This can be either one of the above strings, or an astropy.units equivalency.

rest_value : Quantity

A rest wavelength or frequency with appropriate units. Required if output type is velocity. The cube’s WCS should include this already if the input type is velocity, but the WCS’s rest wavelength/frequency can be overridden with this parameter.

world_spines()

Returns a list of 1D arrays, for the world coordinates along each pixel axis.

Raises error if this operation is ill-posed (e.g. rotated world coordinates, strong distortions)

This method is not currently implemented. Use world() instead.

write(filename, overwrite=False, format=None)

Write the spectral cube to a file.

Parameters:

filename : str

The path to write the file to

format : str

The format of the file to write. (Currently limited to ‘fits’)

overwrite : bool

If True, overwrite filename if it exists