SpectralCube

class spectral_cube.SpectralCube(*args, **kwargs)[source]

Bases: BaseSpectralCube, BeamMixinClass

Attributes Summary

base

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

beam

fill_value

The replacement value used by filled_data.

filled_data

Return a portion of the data array, with excluded mask values

hdu

HDU version of self

hdulist

header

latitude_extrema

longitude_extrema

mask

meta

ndim

Dimensionality of the data

pixels_per_beam

read

Read and parse a dataset and return as a SpectralCube

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

unitless_filled_data

Return a portion of the data array, with excluded mask values

unmasked_data

Return a view of the subset of the underlying data,

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

world_extrema

write

Write this SpectralCube object out in the specified format.

Methods Summary

apply_function(function[, axis, weights, ...])

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)

apply_function_parallel_spatial(function[, ...])

Apply a function in parallel along the spatial dimension.

apply_function_parallel_spectral(function[, ...])

Apply a function in parallel along the spectral dimension.

apply_numpy_function(function[, fill, ...])

Apply a numpy function to the cube

argmax([axis, how])

Return the index of the maximum data value.

argmax_world(axis, **kwargs)

Return the spatial or spectral index of the maximum value along a line of sight.

argmin([axis, how])

Return the index of the minimum data value.

argmin_world(axis, **kwargs)

Return the spatial or spectral index of the minimum value along a line of sight.

check_jybeam_smoothing([raise_error_jybm])

This runs for spatial resolution operations (e.g. spatial_smooth) and either an error or warning when smoothing will affect brightness in Jy/beam operations.

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, update_function])

Convolve each channel in the cube to a specified beam

downsample_axis(factor, axis[, estimator, ...])

Downsample the cube by averaging over factor pixels along an axis.

filled([fill_value])

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 bad values)

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])

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

linewidth_sigma([how])

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

mad_std([axis, how])

Use astropy's mad_std to computer the standard deviation

mask_channels(goodchannels)

Helper function to mask out channels.

max([axis, how])

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

mean([axis, how])

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([axis, how])

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.

plot_channel_maps(nx, ny, channels[, ...])

Make channel maps from a spectral cube

reproject(header[, order, use_memmap, filled])

Spatially reproject the cube into a new header.

sigma_clip_spectrally(threshold[, verbose, ...])

Run astropy's sigma clipper along the spectral axis, converting all bad (excluded) values to NaN.

spatial_filter(ksize, filter[, ...])

Smooth the image in each spatial-spatial plane of the cube using a scipy.ndimage filter.

spatial_smooth(kernel[, convolve, ...])

Smooth the image in each spatial-spatial plane of the cube.

spatial_smooth_median(ksize[, ...])

Smooth the image in each spatial-spatial plane of the cube using a median filter.

spectral_filter(ksize, filter[, use_memmap, ...])

Smooth the cube along the spectral dimension using a scipy.ndimage filter.

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, verbose, ...])

Smooth the cube along the spectral dimension

spectral_smooth_median(ksize[, use_memmap, ...])

Smooth the cube along the spectral dimension

std([axis, how, ddof])

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_crtfregion(crtf_region[, ...])

Extract a masked subcube from a CRTF region.

subcube_from_ds9region(ds9_region[, allow_empty])

Extract a masked subcube from a ds9 region (only functions on celestial dimensions)

subcube_from_mask(region_mask)

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

subcube_from_regions(region_list[, ...])

Extract a masked subcube from a list of regions.Region object (only functions on celestial dimensions)

subcube_slices_from_mask(region_mask[, ...])

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

sum([axis, how])

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

to(unit[, equivalencies])

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 create PV extractions within the GUI

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_beam(beam[, raise_error_jybm])

Attach a beam object to the SpectralCube.

with_fill_value(fill_value)

Create a new object 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.

Attributes Documentation

base

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

beam
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:
dataQuantity

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

pixels_per_beam
read

Read and parse a dataset and return as a SpectralCube

This allows easily reading a dataset in several supported data formats using syntax such as:

>>> from spectral_cube import SpectralCube
>>> cube1 = SpectralCube.read('cube.fits', format='fits')
>>> cube2 = SpectralCube.read('cube.image', format='casa')

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

Get help on the available readers for SpectralCube using the``help()`` method:

>>> SpectralCube.read.help()  # Get help reading SpectralCube and list supported formats
>>> SpectralCube.read.help('fits')  # Get detailed help on SpectralCube FITS reader
>>> SpectralCube.read.list_formats()  # Print list of available formats

See also: http://docs.astropy.org/en/stable/io/unified.html

Parameters:
*argstuple, optional

Positional arguments passed through to data reader. If supplied the first argument is typically the input filename.

formatstr

File format specifier.

**kwargsdict, optional

Keyword arguments passed through to data reader.

Returns:
cubeSpectralCube

SpectralCube corresponding to dataset

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

unitless_filled_data

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

Returns:
datanumpy.array

The masked data.

Notes

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

unmasked_data

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

Returns:
dataQuantity 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
write

Write this SpectralCube object out in the specified format.

This allows easily writing a dataset in many supported data formats using syntax such as:

>>> data.write('data.fits', format='fits')

Get help on the available writers for SpectralCube using the``help()`` method:

>>> SpectralCube.write.help()  # Get help writing SpectralCube and list supported formats
>>> SpectralCube.write.help('fits')  # Get detailed help on SpectralCube FITS writer
>>> SpectralCube.write.list_formats()  # Print list of available formats

See also: http://docs.astropy.org/en/stable/io/unified.html

Parameters:
*argstuple, optional

Positional arguments passed through to data writer. If supplied the first argument is the output filename.

formatstr

File format specifier.

**kwargsdict, optional

Keyword arguments passed through to data writer.

Methods Documentation

apply_function(function, axis=None, weights=None, unit=None, projection=False, progressbar=False, update_function=None, keep_shape=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:
functionfunction

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

axis1, 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.

projectionbool

Return a projection if the resulting array is 2D?

progressbarbool

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

keep_shapebool

If True, the returned object will be the same dimensionality as the cube.

update_functionfunction

An alternative tracker for the progress of applying the function to the cube data. If progressbar is True, this argument is ignored.

Returns:
resultProjection 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_function_parallel_spatial(function, num_cores=None, verbose=0, use_memmap=True, parallel=True, **kwargs)

Apply a function in parallel along the spatial dimension. The function will be performed on data with masked values replaced with the cube’s fill value.

Parameters:
functionfunction

The function to apply in the spatial dimension. It must take two arguments: an array representing an image and a boolean array representing the mask. It may also accept **kwargs. The function must return an object with the same shape as the input spectrum.

num_coresint or None

The number of cores to use if running in parallel

verboseint

Verbosity level to pass to joblib

use_memmapbool

If specified, a memory mapped temporary file on disk will be written to rather than storing the intermediate spectra in memory.

parallelbool

If set to False, will force the use of a single core without using joblib.

kwargsdict

Passed to function

apply_function_parallel_spectral(function, num_cores=None, verbose=0, use_memmap=True, parallel=True, **kwargs)

Apply a function in parallel along the spectral dimension. The function will be performed on data with masked values replaced with the cube’s fill value.

Parameters:
functionfunction

The function to apply in the spectral dimension. It must take two arguments: an array representing a spectrum and a boolean array representing the mask. It may also accept **kwargs. The function must return an object with the same shape as the input spectrum.

num_coresint or None

The number of cores to use if running in parallel

verboseint

Verbosity level to pass to joblib

use_memmapbool

If specified, a memory mapped temporary file on disk will be written to rather than storing the intermediate spectra in memory.

parallelbool

If set to False, will force the use of a single core without using joblib.

kwargsdict

Passed to function

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

Apply a numpy function to the cube

Parameters:
functionNumpy ufunc

A numpy ufunc to apply to the cube

fillfloat

The fill value to use on the data

reducebool

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

howcube | 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’

projectionbool

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?

unitNone 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_endianbool

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.

progressbarbool

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

kwargsdict

Passed to the numpy function.

Returns:
resultProjection 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(axis=None, how='auto', **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:
axisint (optional)

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

howcube | 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’

argmax_world(axis, **kwargs)

Return the spatial or spectral index of the maximum value along a line of sight.

Parameters:
axisint

The axis to return the peak location along. e.g., axis=0 will return the value of the spectral axis at the peak value.

kwargsdict

Passed to argmax.

argmin(axis=None, how='auto', **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:
axisint (optional)

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

howcube | 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_world(axis, **kwargs)

Return the spatial or spectral index of the minimum value along a line of sight.

Parameters:
axisint

The axis to return the peak location along. e.g., axis=0 will return the value of the spectral axis at the peak value.

kwargsdict

Passed to argmin.

check_jybeam_smoothing(raise_error_jybm=True)

This runs for spatial resolution operations (e.g. spatial_smooth) and either an error or warning when smoothing will affect brightness in Jy/beam operations.

This is also true for using the with_beam and with_beams methods, including 1D spectra with Jy/beam units.

Parameters:
raise_error_jybmbool, optional

Raises a BeamUnitsError when True (default). When False, it triggers a BeamWarning.

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:
valueQuantity

The value of the spectral coordinate to search for.

convolve_to(beam, convolve=<function convolve_fft>, update_function=None, **kwargs)

Convolve each channel in the cube to a specified beam

Warning

The current implementation of convolve_to creates an in-memory copy of the whole cube to store the convolved data. Issue #506 notes that this is a problem, and it is on our to-do list to fix.

Parameters:
beamradio_beam.Beam

The beam to convolve to

convolvefunction

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

update_functionmethod

Method that is called to update an external progressbar If provided, it disables the default astropy.utils.console.ProgressBar

kwargsdict

Keyword arguments to pass to the convolution function

Returns:
cubeSpectralCube

A SpectralCube with a single beam

downsample_axis(factor, axis, estimator=<function nanmean>, truncate=False, use_memmap=True, progressbar=True)

Downsample the cube by averaging over factor pixels along an axis. Crops right side if the shape is not a multiple of factor.

The WCS will be ‘downsampled’ by the specified factor as well. If the downsample factor is odd, there will be an offset in the WCS.

There is both an in-memory and a memory-mapped implementation; the default is to use the memory-mapped version. Technically, the ‘large data’ warning doesn’t apply when using the memory-mapped version, but the warning is still there anyway.

Parameters:
myarrndarray

The array to downsample

factorint

The factor to downsample by

axisint

The axis to downsample along

estimatorfunction

defaults to mean. You can downsample by summing or something else if you want a different estimator (e.g., downsampling error: you want to sum & divide by sqrt(n))

truncatebool

Whether to truncate the last chunk or average over a smaller number. e.g., if you downsample [1,2,3,4] by a factor of 3, you could get either [2] or [2,4] if truncate is True or False, respectively.

use_memmapbool

Use a memory map on disk to avoid loading the whole cube into memory (several times)? If set, the warning about large cubes can be ignored (though you still have to override the warning)

progressbarbool

Include a progress bar? Only works with use_memmap=True

filled(fill_value=None)
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_offsetu.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_valueu.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().

mad_std(axis=None, how='cube', **kwargs)

Use astropy’s mad_std to computer the standard deviation

Ignores excluded mask elements.

Parameters:
axisint (optional)

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

howcube | 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’

mask_channels(goodchannels)

Helper function to mask out channels. This function is equivalent to adding a mask with cube[view] where view is broadcastable to the cube shape, but it accepts 1D arrays that are not normally broadcastable.

Parameters:
goodchannelsarray

A 1D boolean array declaring which channels should be kept.

Returns:
cubeSpectralCube

A cube with the specified channels masked

max(axis=None, how='auto', **kwargs)

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

Ignores excluded mask elements.

Parameters:
axisint (optional)

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

howcube | 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(axis=None, how='cube', **kwargs)

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

Ignores excluded mask elements.

Parameters:
axisint (optional)

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

howcube | 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:
axisint (optional)

The axis to collapse

iterate_raysbool

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

Returns:
medndarray

The median

min(axis=None, how='auto', **kwargs)

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

Ignores excluded mask elements.

Parameters:
axisint (optional)

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

howcube | 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, where \(I\) is the intensity in a channel and \(x\) is the spectral coordinate:

Moment 0:

\[M_0 \int I dx\]

Moment 1:

\[M_1 = \frac{\int I x dx}{M_0}\]

Moment N:

\[M_N = \frac{\int I (x - M_1)^N dx}{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:
orderint

The order of the moment to take. Default=0

axisint

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

howcube | 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:
qfloat

The percentile to compute

axisint, or None

Which axis to compute percentiles over

iterate_raysbool

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

plot_channel_maps(nx, ny, channels, contourkwargs={}, output_file=None, fig=None, fig_smallest_dim_inches=8, decimals=3, zoom=1, textcolor=None, cmap='gray_r', tighten=False, textxloc=0.5, textyloc=0.9, savefig_kwargs={}, **kwargs)

Make channel maps from a spectral cube

Parameters:
input_filestr

Name of the input spectral cube

nx, nyint

Number of sub-plots in the x and y direction

channelslist

List of channels to show

cmapstr

The name of a colormap to use for the imshow colors

contourkwargsdict

Keyword arguments passed to contour

textcolorNone or str

Color of the label text to overlay. If None, will be determined automatically. If 'notext', no text will be added.

textxlocfloat
textylocfloat

Text label X,Y-location in axis fraction units

output_filestr

Name of the matplotlib plot

figmatplotlib figure

The figure object to plot onto. Will be overridden to enforce a specific aspect ratio.

fig_smallest_dim_inchesfloat

The size of the smallest dimension (either width or height) of the figure in inches. The other dimension will be selected based on the aspect ratio of the data: it cannot be a free parameter.

decimalsint, optional

Number of decimal places to show in spectral value

zoomint, optional

How much to zoom in. In future versions of this function, the pointing center will be customizable.

tightenbool

Call plt.tight_layout() after plotting?

savefig_kwargsdict

Keyword arguments to pass to savefig (e.g., bbox_inches='tight')

kwargsdict

Passed to imshow

reproject(header, order='bilinear', use_memmap=False, filled=True, **kwargs)

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

If you want to reproject a cube both spatially and spectrally, you need to use spectral_interpolate as well.

Warning

The current implementation of reproject requires that the whole cube be loaded into memory. Issue #506 notes that this is a problem, and it is on our to-do list to fix.

Parameters:
headerastropy.io.fits.Header

A header specifying a cube in valid WCS

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

use_memmapbool

If specified, a memory mapped temporary file on disk will be written to rather than storing the intermediate spectra in memory.

filledbool

Fill the masked values with the cube’s fill value before reprojection? Note that setting filled=False will use the raw data array, which can be a workaround that prevents loading large data into memory.

kwargsdict

Passed to reproject.reproject_interp.

sigma_clip_spectrally(threshold, verbose=0, use_memmap=True, num_cores=None, **kwargs)

Run astropy’s sigma clipper along the spectral axis, converting all bad (excluded) values to NaN.

Parameters:
thresholdfloat

The sigma parameter in astropy.stats.sigma_clip, which refers to the number of sigma above which to cut.

verboseint

Verbosity level to pass to joblib

Other Parameters:
parallelbool

Use joblib to parallelize the operation. If set to False, will force the use of a single core without using joblib.

num_coresint or None

The number of cores to use when applying this function in parallel across the cube.

use_memmapbool

If specified, a memory mapped temporary file on disk will be written to rather than storing the intermediate spectra in memory.

spatial_filter(ksize, filter, update_function=None, raise_error_jybm=True, **kwargs)

Smooth the image in each spatial-spatial plane of the cube using a scipy.ndimage filter.

Parameters:
ksizeint

Size of the filter in pixels (scipy.ndimage.*_filter).

filterfunction

A filter from scipy.ndimage.

update_functionmethod

Method that is called to update an external progressbar If provided, it disables the default astropy.utils.console.ProgressBar

raise_error_jybmbool, optional

Raises a BeamUnitsError when smoothing a cube in Jy/beam units, since the brightness is dependent on the spatial resolution.

kwargsdict

Passed to the convolve function

Other Parameters:
parallelbool

Use joblib to parallelize the operation. If set to False, will force the use of a single core without using joblib.

num_coresint or None

The number of cores to use when applying this function in parallel across the cube.

use_memmapbool

If specified, a memory mapped temporary file on disk will be written to rather than storing the intermediate spectra in memory.

spatial_smooth(kernel, convolve=<function convolve>, raise_error_jybm=True, **kwargs)

Smooth the image in each spatial-spatial plane of the cube.

Parameters:
kernelKernel2D

A 2D kernel from astropy

convolvefunction

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

raise_error_jybmbool, optional

Raises a BeamUnitsError when smoothing a cube in Jy/beam units, since the brightness is dependent on the spatial resolution.

kwargsdict

Passed to the convolve function

Other Parameters:
parallelbool

Use joblib to parallelize the operation. If set to False, will force the use of a single core without using joblib.

num_coresint or None

The number of cores to use when applying this function in parallel across the cube.

use_memmapbool

If specified, a memory mapped temporary file on disk will be written to rather than storing the intermediate spectra in memory.

spatial_smooth_median(ksize, update_function=None, raise_error_jybm=True, filter=<function median_filter>, **kwargs)

Smooth the image in each spatial-spatial plane of the cube using a median filter.

Parameters:
ksizeint

Size of the median filter in pixels (scipy.ndimage.median_filter)

filterfunction

A filter from scipy.ndimage. The default is the median filter.

update_functionmethod

Method that is called to update an external progressbar If provided, it disables the default astropy.utils.console.ProgressBar

raise_error_jybmbool, optional

Raises a BeamUnitsError when smoothing a cube in Jy/beam units, since the brightness is dependent on the spatial resolution.

kwargsdict

Passed to the convolve function

Other Parameters:
parallelbool

Use joblib to parallelize the operation. If set to False, will force the use of a single core without using joblib.

num_coresint or None

The number of cores to use when applying this function in parallel across the cube.

use_memmapbool

If specified, a memory mapped temporary file on disk will be written to rather than storing the intermediate spectra in memory.

spectral_filter(ksize, filter, use_memmap=True, verbose=0, num_cores=None, **kwargs)

Smooth the cube along the spectral dimension using a scipy.ndimage filter.

Parameters:
ksizeint

Size of the filter in spectral channels.

filterfunction

A filter from scipy.ndimage.

Other Parameters:
parallelbool

Use joblib to parallelize the operation. If set to False, will force the use of a single core without using joblib.

num_coresint or None

The number of cores to use when applying this function in parallel across the cube.

use_memmapbool

If specified, a memory mapped temporary file on disk will be written to rather than storing the intermediate spectra in memory.

spectral_interpolate(spectral_grid, suppress_smooth_warning=False, fill_value=None, update_function=None)

Resample the cube spectrally onto a specific grid

Parameters:
spectral_gridarray

An array of the spectral positions to regrid onto

suppress_smooth_warningbool

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_valuefloat

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.

update_functionmethod

Method that is called to update an external progressbar If provided, it disables the default astropy.utils.console.ProgressBar

Returns:
cubeSpectralCube
spectral_slab(lo, hi)

Extract a new cube between two spectral coordinates

Parameters:
lo, hiQuantity

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>, verbose=0, use_memmap=True, num_cores=None, **kwargs)

Smooth the cube along the spectral dimension

Note that the mask is left unchanged in this operation.

Parameters:
kernelKernel1D

A 1D kernel from astropy

convolvefunction

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

verboseint

Verbosity level to pass to joblib

kwargsdict

Passed to the convolve function

Other Parameters:
parallelbool

Use joblib to parallelize the operation. If set to False, will force the use of a single core without using joblib.

num_coresint or None

The number of cores to use when applying this function in parallel across the cube.

use_memmapbool

If specified, a memory mapped temporary file on disk will be written to rather than storing the intermediate spectra in memory.

spectral_smooth_median(ksize, use_memmap=True, verbose=0, num_cores=None, filter=<function median_filter>, **kwargs)

Smooth the cube along the spectral dimension

Parameters:
ksizeint

Size of the median filter (scipy.ndimage.median_filter)

verboseint

Verbosity level to pass to joblib

kwargsdict

Not used at the moment.

Other Parameters:
parallelbool

Use joblib to parallelize the operation. If set to False, will force the use of a single core without using joblib.

num_coresint or None

The number of cores to use when applying this function in parallel across the cube.

use_memmapbool

If specified, a memory mapped temporary file on disk will be written to rather than storing the intermediate spectra in memory.

std(axis=None, how='cube', ddof=0, **kwargs)

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

Parameters:
axisint (optional)

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

howcube | 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’

Other Parameters:
ddofint

Means Delta Degrees of Freedom. The divisor used in calculations is N - ddof, where N represents the number of elements. By default ddof is zero.

Ignores excluded mask elements.
subcube(xlo='min', xhi='max', ylo='min', yhi='max', zlo='min', zhi='max', rest_value=None)

Extract a sub-cube spatially and spectrally.

When spatial WCS dimensions are given as an Quantity, the spatial coordinates of the ‘lo’ and ‘hi’ corners are solved together. This minimizes WCS variations due to the sky curvature when slicing from a large (>1 deg) image.

Parameters:
[xyz]lo/[xyz]hiint 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_crtfregion(crtf_region, allow_empty=False)

Extract a masked subcube from a CRTF region.

Parameters:
crtf_region: str

The CRTF region(s) string 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_ds9region(ds9_region, allow_empty=False)

Extract a masked subcube from a ds9 region (only functions on celestial dimensions)

Parameters:
ds9_region: str

The DS9 region(s) 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: `~spectral_cube.masks.MaskBase` or boolean `numpy.ndarray`

The mask with appropraite WCS or an ndarray with matched coordinates

subcube_from_regions(region_list, allow_empty=False, minimize=True)

Extract a masked subcube from a list of regions.Region object (only functions on celestial dimensions)

Parameters:
region_list: ``regions.Region`` list

The region(s) to extract

allow_empty: bool, optional

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

minimizebool

Run minimal_subcube(). This is mostly redundant, since the bounding box of the region is already used, but it will sometimes slice off a one-pixel rind depending on the details of the region shape. If minimize is disabled, there will potentially be a ring of NaN values around the outside.

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: `~spectral_cube.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(axis=None, how='auto', **kwargs)

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

Ignores excluded mask elements.

Parameters:
axisint (optional)

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

howcube | 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(unit, equivalencies=())

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:
namestr 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_appGlueApplication 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.

datasetglue.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_guibool

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_factorfloat, 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/doc/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_beam(beam, raise_error_jybm=True)[source]

Attach a beam object to the SpectralCube.

Parameters:
beamBeam

Beam object defining the resolution element of the SpectralCube.

with_fill_value(fill_value)

Create a new object 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:
maskMaskBase 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_maskbool (optional, default=True)

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

wcs_toleranceNone 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_cubeSpectralCube

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:
unitUnit

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_valueQuantity

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.