DataArray

class mpdaf.obj.DataArray(filename=None, hdulist=None, data=None, mask=False, var=None, ext=None, unit=Unit(dimensionless), copy=True, dtype=None, primary_header=None, data_header=None, convert_float64=True, **kwargs)[source]

Bases: object

Parent class for Cube, Image and Spectrum.

Its primary purpose is to store pixel values and, optionally, also variances in masked Numpy arrays. For Cube objects these are 3D arrays indexed in the order [wavelength, image_y, image_x]. For Image objects they are 2D arrays indexed in the order [image_y, image_x]. For Spectrum objects they are 1D arrays.

Image arrays hold flat 2D map-projections of the sky. The X and Y axes of the image arrays are orthogonal on the sky at the tangent point of the projection. When the rotation angle of the projection on the sky is zero, the Y axis of the image arrays is along the declination axis, and the X axis is perpendicular to this, with the positive X axis pointing east.

Given a DataArray object obj, the data and variance arrays are accessible via the properties obj.data and obj.var. These two masked arrays share a single array of boolean masking elements, which is also accessible as a simple boolean array via the obj.mask property. The shared mask can be modified through any of the three properties:

obj.data[20:22] = numpy.ma.masked

Is equivalent to:

obj.var[20:22] = numpy.ma.masked

And is also equivalent to:

obj.mask[20:22] = True

All three of the above statements mask pixels 20 and 21 of the data and variance arrays, without changing their values.

Similarly, if one performs an operation on either obj.data or obj.var that modifies the mask, this is reflected in the shared mask of all three properties. For example, the following statement multiplies elements 20 and 21 of the data by 2.0, while changing the shared mask of these pixels to True. In this way the data and the variances of these pixels are consistently masked:

obj.data[20:22] *= numpy.ma.array([2.0,2.0], mask=[True,True])

The data and variance arrays can be completely replaced by assigning new arrays to the obj.data and obj.var properties, but these must have the same shape as before (ie. obj.shape). New arrays that are assigned to obj.data or obj.var can be simple numpy arrays, or a numpy masked arrays.

When a normal numpy array is assigned to obj.data, the obj.mask array is also assigned a mask array, whose elements are True wherever NaN or Inf values are found in the data array. An exception to this rule is if the mask has previously been disabled by assigning numpy.ma.nomask to obj.mask. In this case a mask array is not constructed.

When a numpy masked array is assigned to obj.data, then its mask is also assigned, unchanged, to obj.mask.

Assigning a normal numpy array to the obj.var attribute, has no effect on the contents of obj.mask. On the other hand, when a numpy masked array is assigned to obj.var the obj.mask array becomes the union of its current value and the mask of the provided array.

The ability to record variances for each element is optional. When no variances are stored, obj.var has the value None. To discard an unwanted variance array, None can be subsequently assigned to obj.var.

For cubes and spectra, the wavelengths of the spectral pixels are specified in the .wave member. For cubes and images, the world-coordinates of the image pixels are specified in the .wcs member.

When a DataArray object is constructed from a FITS file, the name of the file and the file’s primary header are recorded. If the data are read from a FITS extension, the header of this extension is also recorded. Alternatively, the primary header and data header can be passed to the DataArray constructor. Where FITS headers are neither provided, nor available in a provided FITS file, generic headers are substituted.

Methods are provided for masking and unmasking pixels, and performing basic arithmetic operations on pixels. Operations that are specific to cubes or spectra or images are provided elsewhere by derived classes.

Parameters:
filenamestr

FITS file name, default to None.

hdulistastropy.fits.HDUList

HDUList object, can be used instead of filename to avoid opening the FITS file multiple times.

dataarray_like

Data array, passed to numpy.ma.MaskedArray.

maskbool or numpy.ma.nomask or numpy.ndarray

Mask used for the creation of the .data MaskedArray. If mask is False (default value), a mask array of the same size of the data array is created. To avoid creating an array, it is possible to use numpy.ma.nomask, but in this case several methods will break if they use the mask.

vararray_like

Variance array, passed to numpy.ma.MaskedArray.

extint or tuple of int or str or tuple of str

Number/name of the data extension, or numbers/names of the data, variance, and optionally mask extensions.

unitastropy.units.Unit

Physical units of the data values, default to u.dimensionless_unscaled.

copybool

If True (default), then the data and variance arrays are copied. Passed to numpy.ma.MaskedArray.

dtypenumpy.dtype

Type of the data. Passed to numpy.ma.MaskedArray.

primary_headerastropy.io.fits.Header

FITS primary header instance.

data_headerastropy.io.fits.Header

FITS data header instance.

fits_kwargsdict

Additional arguments that can be passed to astropy.io.fits.open.

convert_float64bool

By default input arrays or FITS data are converted to float64, in order to increase precision to the detriment of memory usage.

Attributes:
filenamestr

FITS filename.

primary_headerastropy.io.fits.Header

FITS primary header instance.

data_headerastropy.io.fits.Header

FITS data header instance.

wcsmpdaf.obj.WCS

World coordinates.

wavempdaf.obj.WaveCoord

Wavelength coordinates

ndimint

The number of dimensions in the data and variance arrays : int

shapesequence

The lengths of each of the data axes.

datanumpy.ma.MaskedArray

Return data as a numpy.ma.MaskedArray.

varnumpy.ma.MaskedArray

Return variance as a numpy.ma.MaskedArray.

masknumpy.ndarray

The shared masking array of the data and variance arrays.

unitastropy.units.Unit

Physical units of the data values.

dtypenumpy.dtype

The type of the data.

Attributes Summary

data

Return data as a numpy.ma.MaskedArray.

dtype

The type of the data.

mask

The shared masking array of the data and variance arrays.

ndim

The number of dimensions in the data and variance arrays : int

shape

The lengths of each of the data axes.

var

Return variance as a numpy.ma.MaskedArray.

Methods Summary

abs([out])

Return a new object with the absolute value of the data.

clone([data_init, var_init])

Return a shallow copy with the same header and coordinates.

copy()

Return a copy of the object.

crop()

Reduce the size of the array to the smallest sub-array that keeps all unmasked pixels.

get_data_hdu([name, savemask, convert_float32])

Return an ImageHDU corresponding to the DATA extension.

get_dq_hdu([name, header])

Return an ImageHDU corresponding to the DQ (mask) extension.

get_stat_hdu([name, header, convert_float32])

Return an ImageHDU corresponding to the STAT extension.

get_wcs_header()

Return a FITS header containing coordinate descriptions.

info()

Print information.

mask_selection(ksel)

Mask selected pixels.

mask_variance(threshold)

Mask pixels with a variance above a threshold value.

new_from_obj(obj[, data, var, copy, unit])

Create a new object from another one, copying its attributes.

set_wcs([wcs, wave])

Set the world coordinates (spatial and/or spectral where pertinent).

sqrt([out])

Return a new object with positive data square-rooted, and negative data masked.

to_ds9([ds9id, newframe, zscale, cmap])

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

unmask()

Unmask the data (just invalid data (nan,inf) are masked).

write(filename[, savemask, checksum, ...])

Save the data to a FITS file.

Attributes Documentation

data

Return data as a numpy.ma.MaskedArray.

The DataArray constructor postpones reading data from FITS files until they are first used. Read the data array here if not already read.

Changes can be made to individual elements of the data property. When simple numeric values or Numpy array elements are assigned to elements of the data property, the values of these elements are updated and become unmasked.

When masked Numpy values or masked-array elements are assigned to elements of the data property, then these change both the values of the data property and the shared mask of the data and var properties.

Completely new arrays can also be assigned to the data property, provided that they have the same shape as before.

dtype

The type of the data.

mask

The shared masking array of the data and variance arrays.

This is a bool array which has the same shape as the data and variance arrays. Setting an element of this array to True, flags the corresponding pixels of the data and variance arrays, so that they don’t take part in subsequent calculations. Reverting this element to False, unmasks the pixel again.

This array can be modified either directly by assignments to elements of this property or by modifying the masks of the .data and .var arrays. An entirely new mask array can also be assigned to this property, provided that it has the same shape as the data array.

ndim

The number of dimensions in the data and variance arrays : int

shape

The lengths of each of the data axes.

var

Return variance as a numpy.ma.MaskedArray.

If variances have been provided for each data pixel, then this property can be used to record those variances. Normally this is a masked array which shares the mask of the data property. However if no variances have been provided, then this property is None.

Variances are typically provided along with the data values in the originating FITS file. Alternatively a variance array can be assigned to this property after the data have been read.

Note that any function that modifies the contents of the data array may need to update the array of variances accordingly. For example, after scaling pixel values by a constant factor c, the variances should be scaled by c**2.

When masked-array values are assigned to elements of the var property, the mask of the new values is assigned to the shared mask of the data and variance properties.

Completely new arrays can also be assigned to the var property. When a masked array is assigned to the var property, its mask is combined with the existing shared mask, rather than replacing it.

Methods Documentation

abs(out=None)[source]

Return a new object with the absolute value of the data.

Parameters:
outmpdaf.obj.DataArray, optional

Array of the same shape as input, into which the output is placed. By default, a new array is created.

clone(data_init=None, var_init=None)[source]

Return a shallow copy with the same header and coordinates.

Optionally fill the cloned array using values returned by provided functions.

Parameters:
data_initcallable(), optional

An optional function to use to create the data array (it takes the shape as parameter). For example np.zeros or np.empty can be used. It defaults to None, which results in the data attribute being None.

var_initcallable(), optional

An optional function to use to create the variance array, with the same specifics as data_init. This default to None, which results in the var attribute being assigned None.

copy()[source]

Return a copy of the object.

crop()[source]

Reduce the size of the array to the smallest sub-array that keeps all unmasked pixels.

This removes any margins around the array that only contain masked pixels. If all pixels are masked in the input cube, the data and variance arrays are deleted.

Returns:
itemlist of slice

The slices that were used to extract the sub-array.

get_data_hdu(name='DATA', savemask='dq', convert_float32=True)[source]

Return an ImageHDU corresponding to the DATA extension.

Parameters:
namestr

Extension name, DATA by default.

savemaskstr

If ‘dq’, the mask array is saved in a DQ extension. If ‘nan’, masked data are replaced by nan in a DATA extension. If ‘none’, masked array is not saved.

convert_float32bool

By default float64 arrays are converted to float32, in order to produce smaller files.

Returns:
astropy.io.fits.ImageHDU
get_dq_hdu(name='DQ', header=None)[source]

Return an ImageHDU corresponding to the DQ (mask) extension.

get_stat_hdu(name='STAT', header=None, convert_float32=True)[source]

Return an ImageHDU corresponding to the STAT extension.

Parameters:
namestr

Extension name, STAT by default.

headerastropy.io.fits.Header

Fits Header to put in the extension, typically to reuse the same as in the DATA extension. Otherwise it is created with the wcs.

convert_float32bool

By default float64 arrays are converted to float32, in order to produce smaller files.

Returns:
astropy.io.fits.ImageHDU
get_wcs_header()[source]

Return a FITS header containing coordinate descriptions.

info()[source]

Print information.

mask_selection(ksel)[source]

Mask selected pixels.

Parameters:
kseloutput of np.where

Elements depending on a condition

mask_variance(threshold)[source]

Mask pixels with a variance above a threshold value.

Parameters:
thresholdfloat

Threshold value.

classmethod new_from_obj(obj, data=None, var=None, copy=False, unit=None)[source]

Create a new object from another one, copying its attributes.

Parameters:
objmpdaf.obj.DataArray

The object to use as the template for the new object. This should be an object based on DataArray, such as an Image, Cube or Spectrum.

dataarray_like, optional

An optional data array, or None to indicate that obj.data should be used. The default is None.

vararray_like, optional

An optional variance array, or None to indicate that obj.var should be used, or False to indicate that the new object should not have any variances. The default is None.

copybool

Copy the data and variance arrays if True (default False).

set_wcs(wcs=None, wave=None)[source]

Set the world coordinates (spatial and/or spectral where pertinent).

Parameters:
wcsmpdaf.obj.WCS, optional

Spatial world coordinates. This argument is ignored when self is a Spectrum.

wavempdaf.obj.WaveCoord, optional

Spectral wavelength coordinates. This argument is ignored when self is an Image.

sqrt(out=None)[source]

Return a new object with positive data square-rooted, and negative data masked.

Parameters:
outmpdaf.obj.DataArray, optional

Array of the same shape as input, into which the output is placed. By default, a new array is created.

to_ds9(ds9id=None, newframe=False, zscale=False, cmap='grey')[source]

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

Parameters:
ds9idstr, optional

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

newframebool

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

unmask()[source]

Unmask the data (just invalid data (nan,inf) are masked).

write(filename, savemask='dq', checksum=False, convert_float32=True)[source]

Save the data to a FITS file.

Overwrite the file if it exists.

Parameters:
filenamestr

The FITS filename.

savemaskstr

If ‘dq’, the mask array is saved in a DQ extension If ‘nan’, masked data are replaced by nan in the DATA extension If ‘none’, masked array is not saved

checksumbool

If True, adds both DATASUM and CHECKSUM cards to the headers of all HDU’s written to the file.

convert_float32bool

By default float64 arrays are converted to float32, in order to produce smaller files.