datacube module
Provide the DataCube class for creating a data cube.
- class martini.datacube.DataCube(*, n_px_x: int, n_px_y: int, n_channels: int, px_size: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("arcsec")], channel_width: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("km / s")], spectral_centre: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("km / s")] = <Quantity 0. km / s>, ra: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("deg")] = <Quantity 0. deg>, dec: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("deg")] = <Quantity 0. deg>, stokes_axis: bool = False, coordinate_frame: BaseRADecFrame = <ICRS Frame>, specsys: str = 'icrs', velocity_centre: None = None)[source]
Bases:
objectHandles creation and management of the data cube itself.
Basic usage simply involves initializing with the parameters listed below. More advanced usage might arise if designing custom classes for other sub- modules, especially beams. To initialize a
DataCubefrom a saved state, seeload_state().- Parameters:
n_px_x (int) – Pixel count along the x (RA) axis. Even integers strongly preferred.
n_px_y (int) – Pixel count along the y (Dec) axis. Even integers strongly preferred.
n_channels (int) – Number of channels along the spectral axis.
px_size (Quantity) –
Quantity, with dimensions of angle. Angular scale of one pixel.channel_width (Quantity) –
Quantity, with dimensions of velocity or frequency. Step size along the spectral axis. Can be provided as a velocity or a frequency.spectral_centre (Quantity, optional) –
Quantitywith dimensions of velocity or frequency. Velocity (or frequency) of the centre along the spectral axis. Defaults to 0 km/s, normally set this to the systemic velocity of the source.ra (Quantity, optional) –
Quantity, with dimensions of angle. Right ascension of the cube centroid. Defaults to 0 degrees.dec (Quantity, optional) –
Quantitywith dimensions of angle. Declination of the cube centroid. Defaults to 0 degrees.stokes_axis (bool, optional) – Whether the datacube should be initialized with a Stokes’ axis.
coordinate_frame (BaseRADecFrame, optional) – The coordinate frame of the World Coordinate System (WCS) associated with the data cube. Recommended frames are
GCRS,ICRS,HCRS,LSRK,LSRDorLSR. The frame should be passed initialized, e.g.ICRS()(not justICRS).specsys (str, optional) – The spectral reference frame (standard of rest) of the World Coordinate System (WCS) associated with the data cube. Some common options include
"gcrs","icrs","hcrs","lsrk","lsrd","lsr". For a complete list, useastropy.coordinates.frame_transform_graph.get_names().velocity_centre (Quantity, deprecated) – Deprecated, use spectral centre instead.
Notes
The
channel_widthargument defines the channel spacing in either frequency or velocity units. If provided with units of frequency, the data cube channels will be evenly spaced in frequency. Conversely, if provided with units of velocity, the channels will be evenly spaced in velocity. Thespectral_centreis related and can also have units of frequency or velocity, but will be converted to have the same units as thechannel_widthin order to define the channels.- add_pad(pad: tuple[int, int]) None[source]
Resize the cube to add a padding region in the spatial direction.
Accurate convolution with a beam requires a cube padded according to the size of the beam kernel (its representation sampled on a grid with the same spacing). The beam class is required to handle defining the size of pad required.
- Parameters:
pad (tuple) – 2-tuple (or other sequence) containing the number of pixels to add in the x (RA) and y (Dec) directions.
See also
- property channel_edges: Annotated[Quantity, Unit('Hz')] | Annotated[Quantity, Unit('m / s')]
The edges of the channels from the coordinate system.
- property channel_maps: Iterator[Quantity]
An iterator over the channel maps.
An alias for
spatial_slices().- Returns:
The iterator over the spatial ‘slices’ of the cube.
- Return type:
iter
- property channel_mids: Annotated[Quantity, Unit('Hz')] | Annotated[Quantity, Unit('m / s')]
The centres of the channels from the coordinate system.
- copy() Self[source]
Produce a copy of the
DataCube.May be especially useful to create multiple datacubes with differing intermediate steps.
- drop_pad() None[source]
Remove the padding added using
add_pad().After convolution, the pad region contains meaningless information and can be discarded.
See also
- property frequency_channel_edges: Annotated[Quantity, Unit('Hz')]
The edges of the channels from the coordinate system in frequency units.
- property frequency_channel_mids: Annotated[Quantity, Unit('Hz')]
The centres of the channels from the coordinate system in frequency units.
- classmethod from_wcs(input_wcs: WCS, specsys: str | None = None) Self[source]
Create a DataCube from a World Coordinate System (WCS).
The WCS can be obtained, for instance, from a FITS header.
To create a MARTINI data cube with pixels and channels exactly matching an observed data cube (stored in a FITS file) would be a bit tedious using the usual
__init__()since the number of pixels, channel spacing and so on must be specified by hand. This function instead constructs aDataCubefrom aWCSinstance, that can be easily obtained from a FITS file or header (see example below). The resulting data cube has exactly the same dimensions (pixels and channels) and World Coordinate System (WCS) as the input WCS.- Parameters:
See also
Notes
MARTINI’s data cubes have a fixed axis ordering: first the RA axis, then the Dec axis, then the spectral axis, and finally the Stokes’ axis (if present). A
DataCubecreated withfrom_wcs()may therefore have its axes re-ordered (tranposed) relative to theinput_wcs. The FITS files output by MARTINI have the same axis ordering, and may therefore also be transposed relative to a data cube used to construct theWCSfor theinput_wcsargument.Examples
It is easy to initialize a
WCSfrom a FITS-format header. For example, given a FITS filemy_cube.fits, setting up aDataCubewith matching World Coordinate System (WCS) looks like:from astropy import wcs from astropy.io import fits from martini.datacube import DataCube with fits.open("my_cube.fits") as fitsfile: fits_hdr = fitsfile[0].header # header of the main HDU fits_wcs = wcs.WCS(fits_hdr) datacube = DataCube.from_wcs(fits_wcs)
- classmethod load_state(filename: str) Self[source]
Initialize a
DataCubefrom a state saved.Note that
h5pymust be installed for use. Note that ONLY theDataCubestate is restored, other modules and their configurations are not affected.- Parameters:
filename (str) – File to open.
- Returns:
A suitably initialized
DataCubeobject.- Return type:
See also
- save_state(filename: str, overwrite: bool = False) None[source]
Write the
DataCubestate to file for re-use.The state can be re-initialized (see
load_state()). Note thath5pymust be installed for use. NOT for outputting mock observations, for this seewrite_fits()andwrite_hdf5().- Parameters:
filename (str) – File to write.
overwrite (bool) – Whether to allow overwriting existing files. (default:
False).
See also
- property spatial_slices: Iterator[Quantity]
An iterator over the spatial ‘slices’ of the cube.
- Returns:
The iterator over the spatial ‘slices’ of the cube.
- Return type:
iter
- property spectra: Iterator[Quantity]
An iterator over the spectra (one in each spatial pixel).
- Returns:
The iterator over the spectra making up the cube.
- Return type:
iter
- property units: tuple[Annotated[Quantity, Unit('deg')], Annotated[Quantity, Unit('deg')], Annotated[Quantity, Unit('Hz')] | Annotated[Quantity, Unit('m / s')]] | tuple[Annotated[Quantity, Unit('deg')], Annotated[Quantity, Unit('deg')], Annotated[Quantity, Unit('Hz')] | Annotated[Quantity, Unit('m / s')], Annotated[Quantity, Unit(dimensionless)]]
The units of the DataCube’s World Coordinate System (WCS).
- Returns:
A tuple containing
Unitinstances describing the WCS units.- Return type:
tuple
- property velocity_channel_edges: Annotated[Quantity, Unit('m / s')]
The edges of the channels from the coordinate system in velocity units.