SpectrogramCubeSequence#

class irispy.spectrograph.SpectrogramCubeSequence(data_list, meta=None, common_axis=0, **kwargs)[source]#

Bases: SpectrogramSequence

Class representing spectrogram data described by a collection of separate WCSes.

So each individual SpectrogramCube within represents a single complete raster scan. The sequence contains multiple such cubes till the end of the observation.

Parameters:
  • data_list (list) – List of SpectrogramCube objects from the same spectral window and OBS ID.

  • meta (dict or header object, optional) – Metadata associated with the sequence.

  • common_axis (int, optional) – The axis of the NDCubes corresponding to time.

Attributes Summary

array_axis_physical_types

The physical types associated with each array axis, including the sequence axis.

celestial

Return the celestial coordinates for each pixel.

common_axis_coords

The coordinate values at each location along the common axis across all cubes.

cube_like_array_axis_physical_types

The physical types associated with each array axis, omitting the sequence axis.

cube_like_dimensions

The length of each array axis as if all cubes were concatenated along the common axis.

cube_like_shape

The length of each array axis as if all cubes were concatenated along the common axis.

dimensions

The length of each axis including the sequence axis.

exposure_time

Return the exposure time for each exposure.

index_as_cube

Slice the NDCubesequence instance as a single cube concatenated along the common axis.

plotter

sequence_axis_coords

Return the coordinate values along the sequence axis.

shape

spectral_axis

Return the spectral coordinates for each pixel.

time

Return the time coordinates for each pixel.

Methods Summary

apply_exposure_time_correction([undo, copy, ...])

Applies or undoes exposure time correction to data and uncertainty and adjusts unit.

crop(*points[, wcses])

Crop cubes in sequence to smallest pixel-space bounding box containing the input points.

crop_by_values(*points[, units, wcses])

Crop cubes in sequence to smallest pixel-space bounding box containing the input points.

explode_along_axis(axis)

Separates slices of N-D cubes along a given cube axis into (N-1)D cubes.

plot(*args, **kwargs)

Visualize the NDCubeSequence.

plot_as_cube(*args, **kwargs)

Attributes Documentation

array_axis_physical_types#

The physical types associated with each array axis, including the sequence axis.

celestial#
common_axis_coords#

The coordinate values at each location along the common axis across all cubes.

Only coordinates associated with the common axis in all cubes in the sequence are returned. Coordinates from different cubes are concatenated along the common axis. They thus represent the coordinate values at each location as if all cubes in the sequence were concatenated along the common axis.

cube_like_array_axis_physical_types#

The physical types associated with each array axis, omitting the sequence axis.

cube_like_dimensions#

The length of each array axis as if all cubes were concatenated along the common axis.

cube_like_shape#

The length of each array axis as if all cubes were concatenated along the common axis.

dimensions#

The length of each axis including the sequence axis.

exposure_time#
index_as_cube#

Slice the NDCubesequence instance as a single cube concatenated along the common axis.

Example

>>> # Say we have three Cubes each cube has common_axis=0 is time and shape=(3,3,3)
>>> data_list = [cubeA, cubeB, cubeC]
>>> cs = NDCubeSequence(data_list, meta=None, common_axis=0)
>>> # return zeroth time slice of cubeB in via normal NDCubeSequence indexing.
>>> cs[1,:,0,:]
>>> # Return same slice using this function
>>> cs.index_as_cube[3:6, 0, :]
plotter = None#
sequence_axis_coords#

Return the coordinate values along the sequence axis.

These are compiled from the GlobalCoords objects attached to each NDCube where each cube represents a location along the sequence axis. Only coordinates that are common to all cubes are returned.

shape#
spectral_axis#
time#

Methods Documentation

apply_exposure_time_correction(
undo=False,
copy=False,
force=False,
)#

Applies or undoes exposure time correction to data and uncertainty and adjusts unit.

Correction is only applied (undone) if the object’s unit doesn’t (does) already include inverse time. This can be overridden so that correction is applied (undone) regardless of unit by setting force=True.

Parameters:
  • undo (bool) – If False, exposure time correction is applied. If True, exposure time correction is removed. Default=False

  • copy (bool) – If True a new instance with the converted data values is returned. If False, the current instance is overwritten. Default=False

  • force (bool) – If not True, applies (undoes) exposure time correction only if unit doesn’t (does) already include inverse time. If True, correction is applied (undone) regardless of unit. Unit is still adjusted accordingly.

Returns:

If copy=False, the original input is modified with the exposure time correction applied (undone). If copy=True, a new sunraster.SpectrogramSequence is returned with the correction applied (undone).

Return type:

None or sunraster.SpectrogramSequence

crop(*points, wcses=None)#

Crop cubes in sequence to smallest pixel-space bounding box containing the input points.

Each input point is given as a tuple of high-level world coordinate objects. This method does not crop the sequence axis. Instead input points are passed to the crop method of cubes in sequence. Note, therefore, that the input points do not include an entry for world coords only associated with the sequence axis. For a description of how the bounding box is defined, see the docstring of the ndcube.NDCube.crop() method. In cases where the cubes are not aligned, all cubes are cropped to the same region in pixel space. This region will be the smallest that encompasses the input points in all cubes while maintaining consistent array shape between the cubes.

Parameters:
  • points – Passed to ndcube.NDCube.crop() as the points arg without checking.

  • wcses (iterable of WCS objects or str, optional) – The WCS objects to be used to crop the cubes. There must by one WCS per cube. Alternatively, can be a string giving the name of the cube wcs attribute to be used, namely, ‘wcs’, ‘combined_wcs’, or ‘extra_coords’. Default=None is equivalent to ‘wcs’.

Returns:

The cropped sequence.

Return type:

NDCubeSequence

crop_by_values(*points, units=None, wcses=None)#

Crop cubes in sequence to smallest pixel-space bounding box containing the input points.

Each input point is given as a tuple of low-level world coordinate objects. This method does not crop the sequence axis. Instead input points are passed to the ndcube.NDCube.crop_by_values() method of cubes in sequence. Note, therefore, that the input points do not include an entry for world coords only associated with the sequence axis. For a description of how the bounding box is defined, see the docstring of the NDCube.crop_by_values method. In cases where the cubes are not aligned, all cubes are cropped to the same region in pixel space. This region will be the smallest that encompasses the input points in all cubes while maintaining consistent array shape between the cubes.

Parameters:
  • points – Passed to ndcube.NDCube.crop() as the points arg without checking.

  • units – Passed to ndcube.NDCube.crop_by_values() as the units kwarg.

  • wcses (iterable of WCS objects or str, optional) – The WCS objects to be used to crop the cubes. There must by one WCS per cube. Alternatively, can be a string giving the name of the cube wcs attribute to be used, namely, ‘wcs’, ‘combined_wcs’, or ‘extra_coords’. Default=None is equivalent to ‘wcs’.

Returns:

The cropped sequence.

Return type:

NDCubeSequence

explode_along_axis(axis)#

Separates slices of N-D cubes along a given cube axis into (N-1)D cubes.

Parameters:

axis (int) – The axis along which the data is to be changed.

Returns:

New sequence of (N-1)D cubes broken up along given axis.

Return type:

ndcube.NDCubeSequence

plot(*args, **kwargs)[source]#

Visualize the NDCubeSequence.

Parameters:
plot_as_cube(*args, **kwargs)#