Cube

class mpdaf.obj.Cube(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: mpdaf.obj.arithmetic.ArithmeticMixin, mpdaf.obj.DataArray

This class manages Cube objects, which contain images at multiple wavelengths. The images are stored in 3D numpy.ndarray arrays. The axes of these arrays are, image wavelength, image y-axis, image x-axis. For MUSE images, the y-axis is typically aligned with the declination axis on the sky, but in general the y and x axes are not aligned with either declination or right-ascension. The actual orientation can be queried via the get_rot() method.

The 3D cube of images is contained in a property of the cube called .data. There is also a .var member. When variances are available, .var is a 3D array that contains the variances of each pixel in the .data array. When variances have not been provided, the .var member is given the value, None.

The .data and the .var are either numpy masked arrays or numpy ndarrays. When they are masked arrays, their shared mask can also be accessed separately via the .mask member, which holds a 3D array of bool elements. When .data and .var are normal numpy.ndarrays, the .mask member is given the value numpy.ma.nomask.

When a new Cube object is created, the data, variance and mask arrays can either be specified as arguments, or the name of a FITS file can be provided to load them from.

Parameters
filenamestr

An optional FITS file name from which to load the cube. None by default. This argument is ignored if the data argument is not None.

extint or (int,int) or str or (str,str)

The optional number/name of the data extension or the numbers/names of the data and variance extensions.

wcsmpdaf.obj.WCS

The world coordinates of the image pixels.

wavempdaf.obj.WaveCoord

The wavelength coordinates of the spectral pixels.

unitstr or astropy.units.Unit

The physical units of the data values. Defaults to astropy.units.dimensionless_unscaled.

dataarray_like

An optional array containing the values of each pixel in the cube (None by default). Where given, this array should be 3 dimensional, and the python ordering of its axes should be (wavelength,image_y,image_x).

vararray_like

An optional array containing the variances of each pixel in the cube (None by default). Where given, this array should be 3 dimensional, and the python ordering of its axes should be (wavelength,image_y,image_x).

copybool

If true (default), then the data and variance arrays are copied.

dtypenumpy.dtype

The type of the data (int, float)

Attributes
filenamestr

The name of the originating FITS file, if any. Otherwise None.

primary_headerastropy.io.fits.Header

The FITS primary header instance, if a FITS file was provided.

data_headerastropy.io.fits.Header

The FITS header of the DATA extension.

wcsmpdaf.obj.WCS

The world coordinates of the image pixels.

wavempdaf.obj.WaveCoord

The wavelength coordinates of the spectral pixels.

unitastropy.units.Unit

The 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(self[, out])

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

aperture(self, center, radius[, …])

Extract the spectrum of a circular aperture of given radius.

bandpass_image(self, wavelengths, sensitivities)

Given a cube of images versus wavelength and the bandpass filter-curve of a wide-band monochromatic instrument, extract an image from the cube that has the spectral response of the monochromatic instrument.

clone(self[, data_init, var_init])

Return a shallow copy with the same header and coordinates.

convolve(self, other[, inplace])

Convolve a Cube with a 3D array or another Cube, using the discrete convolution equation.

copy(self)

Return a copy of the object.

crop(self)

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

fftconvolve(self, other[, inplace])

Convolve a Cube with a 3D array or another Cube, using the Fourier convolution theorem.

get_band_image(self, name)

Generate an image using a known filter.

get_data_hdu(self[, name, savemask, …])

Return an ImageHDU corresponding to the DATA extension.

get_dq_hdu(self[, name, header])

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

get_end(self[, unit_wave, unit_wcs])

Return [lbda,y,x] at the center of pixel (-1,-1,-1).

get_image(self, wave[, is_sum, …])

Generate an image aggregating over a wavelength range.

get_range(self[, unit_wave, unit_wcs])

Return the range of wavelengths, declinations and right ascensions in the cube.

get_rot(self[, unit])

Return the rotation angle of the images of the cube, defined such that a rotation angle of zero aligns north along the positive Y axis, and a positive rotation angle rotates north away from the Y axis, in the sense of a rotation from north to east.

get_start(self[, unit_wave, unit_wcs])

Return [lbda,y,x] at the center of pixel (0,0,0).

get_stat_hdu(self[, name, header, …])

Return an ImageHDU corresponding to the STAT extension.

get_step(self[, unit_wave, unit_wcs])

Return the cube steps [dlbda,dy,dx].

get_wcs_header(self)

Return a FITS header containing coordinate descriptions.

info(self)

Print information.

loop_ima_multiprocessing(self, f[, cpu, verbose])

Use multiple processes to run a function on each image of a cube.

loop_spe_multiprocessing(self, f[, cpu, verbose])

Use multiple processes to run a function on each spectrum of a cube.

mask_ellipse(self, center, radius, posangle)

Mask values inside or outside an elliptical region.

mask_polygon(self, poly[, lmin, lmax, …])

Mask values inside or outside a polygonal region.

mask_region(self, center, radius[, lmin, …])

Mask values inside or outside a circular or rectangular region.

mask_selection(self, ksel)

Mask selected pixels.

mask_variance(self, threshold)

Mask pixels with a variance above a threshold value.

max(self[, axis])

Return the maximum over a given axis.

mean(self[, axis, weights])

Return a weighted or unweighted mean over a given axis or axes.

median(self[, axis])

Return the median over a given axis.

min(self[, axis])

Return the minimum over a given axis.

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

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

rebin(self, factor[, margin, inplace])

Combine neighboring pixels to reduce the size of a cube by integer factors along each axis.

select_lambda(self, lbda_min[, lbda_max, …])

Return the image or sub-cube corresponding to a wavelength range.

set_wcs(self[, wcs, wave])

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

spatial_erosion(self, npixels[, inplace])

Remove n pixels around the masked white image.

sqrt(self[, out])

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

subcube(self, center, size[, lbda, …])

Return a subcube view around a position and for a wavelength range.

subcube_circle_aperture(self, center, radius)

Extract a sub-cube that encloses a circular aperture of a specified radius and for a given wavelength range.

sum(self[, axis, weights])

Return a weighted or unweighted sum over a given axis or axes.

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

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

truncate(self, coord[, mask, unit_wave, …])

Return a sub-cube bounded by specified wavelength and spatial world-coordinates.

unmask(self)

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

write(self, 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(self, out=None)

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.

aperture(self, center, radius, unit_center=Unit("deg"), unit_radius=Unit("arcsec"), is_sum=True)[source]

Extract the spectrum of a circular aperture of given radius.

A spectrum is formed by summing/averaging the pixels within a specified circular region of each wavelength image. This yields a spectrum that has the same length as the wavelength axis of the cube.

Parameters
center(float,float)

The center of the aperture (y,x).

radiusfloat

The radius of the aperture. If the radius is None, or <= 0, the spectrum of the nearest image pixel to the specified center is returned.

unit_centerastropy.units.Unit

The units of the center coordinates (degrees by default) The special value, None, indicates that center is a 2D pixel index.

unit_radiusastropy.units.Unit

The units of the radius argument (arcseconds by default) The special value, None, indicates that the radius is specified in pixels.

is_sumbool

If True, compute the sum of the pixels, otherwise compute the arithmetic mean of the pixels.

Returns
Spectrum
bandpass_image(self, wavelengths, sensitivities, unit_wave=Unit("Angstrom"), interpolation='linear')[source]

Given a cube of images versus wavelength and the bandpass filter-curve of a wide-band monochromatic instrument, extract an image from the cube that has the spectral response of the monochromatic instrument.

For example, this can be used to create a MUSE image that has the same spectral characteristics as an HST image. The MUSE image can then be compared to the HST image without having to worry about any differences caused by different spectral sensitivities.

For each channel n of the cube, the filter-curve is integrated over the width of that channel to obtain a weight, w[n]. The output image is then given by the following weighted mean:

output_image = sum(w[n] * cube_image[n]) / sum(w[n])

In practice, to accomodate masked pixels, the w[n] array is expanded into a cube w[n,y,x], and the weights of individual masked pixels in the cube are zeroed before the above equation is applied.

If the wavelength axis of the cube only partly overlaps the bandpass of the filter-curve, the filter curve is truncated to fit within the bounds of the wavelength axis. A warning is printed to stderr if this occurs, because this results in an image that lacks flux from some of the wavelengths of the requested bandpass.

Parameters
wavelengthsnumpy.ndarray

An array of the wavelengths of the filter curve, listed in ascending order of wavelength. Outside the listed wavelengths the filter-curve is assumed to be zero.

sensitivitiesnumpy.ndarray

The relative flux sensitivities at the wavelengths in the wavelengths array. These sensititivies will be normalized, so only their relative values are important.

unit_waveastropy.units.Unit

The units used in the array of wavelengths. The default is angstroms. To specify pixel units, pass None.

interpolationstr

The form of interpolation to use to integrate over the filter curve. This should be one of:

"linear"     : Linear interpolation
"cubic"      : Cubic spline interpolation (very slow)

The default is linear interpolation. If the filter curve is well sampled and its sampling interval is narrower than the wavelength pixels of the cube, then this should be sufficient. Alternatively, if the sampling interval is significantly wider than the wavelength pixels of the cube, then cubic interpolation should be used instead. Beware that cubic interpolation is much slower than linear interpolation.

Returns
Image

An image formed from the filter-weighted mean of channels in the cube that overlap the bandpass of the filter curve.

clone(self, data_init=None, var_init=None)

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.

convolve(self, other, inplace=False)[source]

Convolve a Cube with a 3D array or another Cube, using the discrete convolution equation.

This function, which uses the discrete convolution equation, is usually slower than Cube.fftconvolve(). However it can be faster when other.data.size is small, and it always uses much less memory, so it is sometimes the only practical choice.

Masked values in self.data and self.var are replaced with zeros before the convolution is performed, but they are masked again after the convolution.

If self.var exists, the variances are propagated using the equation:

result.var = self.var (*) other**2

where (*) indicates convolution. This equation can be derived by applying the usual rules of error-propagation to the discrete convolution equation.

The speed of this function scales as O(Nd x No) where Nd=self.data.size and No=other.data.size.

Uses scipy.signal.convolve.

Parameters
otherCube or numpy.ndarray

The 3D array with which to convolve the cube in self.data. This can be an 3D array of the same size as self, or it can be a smaller array, such as a small 3D gaussian to use to smooth the larger cube.

When other contains a symmetric filtering function, such as a 3-dimensional gaussian, the center of the function should be placed at the center of pixel:

(other.shape - 1) // 2

If other is an MPDAF Cube object, note that only its data array is used. Masked values in this array are treated as zero, and any variances found in other.var are ignored.

inplacebool

If False (the default), return the results in a new Cube. If True, record the result in self and return that.

Returns
Cube
copy(self)

Return a copy of the object.

crop(self)

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.

fftconvolve(self, other, inplace=False)[source]

Convolve a Cube with a 3D array or another Cube, using the Fourier convolution theorem.

This function, which performs the convolution by multiplying the Fourier transforms of the two arrays, is usually much faster than Cube.convolve(), except when other.data.size is small. However it uses much more memory, so Cube.convolve() is sometimes a better choice.

Masked values in self.data and self.var are replaced with zeros before the convolution is performed, but they are masked again after the convolution.

If self.var exists, the variances are propagated using the equation:

result.var = self.var (*) other**2

where (*) indicates convolution. This equation can be derived by applying the usual rules of error-propagation to the discrete convolution equation.

The speed of this function scales as O(Nd x log(Nd)) where Nd=self.data.size. It temporarily allocates a pair of arrays that have the sum of the shapes of self.shape and other.shape, rounded up to a power of two along each axis. This can involve a lot of memory being allocated. For this reason, when other.shape is small, Cube.convolve() may be more efficient than Cube.fftconvolve().

Uses scipy.signal.fftconvolve.

Parameters
otherCube or numpy.ndarray

The 3D array with which to convolve the cube in self.data. This array can be the same size as self, or it can be a smaller array, such as a small 3D gaussian to use to smooth the larger cube.

When other contains a symmetric filtering function, such as a 3-dimensional gaussian, the center of the function should be placed at the center of pixel:

(other.shape - 1) // 2

If other is an MPDAF Cube object, note that only its data array is used. Masked values in this array are treated as zero, and any variances found in other.var are ignored.

inplacebool

If False (the default), return the results in a new Cube. If True, record the result in self and return that.

Returns
Cube
get_band_image(self, name)[source]

Generate an image using a known filter.

Parameters
namestr

Filter name. Must exist in the filter file taken from the MUSE DRS (lib/mpdaf/obj/filters/filter_list.fits). Available filters: Johnson_B, Johnson_V, Cousins_R, Cousins_I, SDSS_u, SDSS_g, SDSS_r, SDSS_i, SDSS_z, ACS_F475W, ACS_F550M, ACS_F555W, ACS_F606W, ACS_F625W, ACS_F775W, ACS_F814W, WFPC2_F555W, WFPC2_F675W, WFPC2_F814W, WFC3_F502N, WFC3_F555W, WFC3_F606W, WFC3_F625W, WFC3_F656N, WFC3_F775W, WFC3_F814W, Kron_V

get_data_hdu(self, name='DATA', savemask='dq', convert_float32=True)

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(self, name='DQ', header=None)

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

get_end(self, unit_wave=None, unit_wcs=None)[source]

Return [lbda,y,x] at the center of pixel (-1,-1,-1).

Parameters
unit_waveastropy.units.Unit

The wavelenth units to use for the starting wavelength. The default value, None, selects the intrinsic wavelength units of the cube.

unit_wcsastropy.units.Unit

The world coordinates units to use for the spatial starting position. The default value, None, selects the intrinsic world coordinates of the cube.

get_image(self, wave, is_sum=False, subtract_off=False, margin=10.0, fband=3.0, median_filter=0, unit_wave=Unit("Angstrom"), method='mean')[source]

Generate an image aggregating over a wavelength range.

This method creates an image aggregating all the slices between a wavelength range.

Parameters
wave(float, float)

The (lbda1,lbda2) interval of wavelength in angstrom.

unit_waveastropy.units.Unit

The wavelength units of the lbda1, lbda2, margin and median filter width arguments (angstrom by default). If None, lbda1, lbda2, margin and median filter width should be in pixels.

is_sumbool

If True, compute the sum of the images. Deprecated, use “sum” as aggregation method.

subtract_offbool

If True, subtract off a background image that is estimated from median filtering (if median_filter is not zero) or by combining some images from both below and above the chosen wavelength range. If the number of images between lbda1 and lbda2 is denoted, N, then the number of background images taken from below and above the wavelength range are:

nbelow = nabove = (fband * N) / 2   [rounded up to an integer]

where fband is a parameter of this function.

The wavelength ranges of the two groups of background images below and above the chosen wavelength range are separated from lower and upper edges of the chosen wavelength range by a value, margin, which is an argument of this function.

This scheme was developed by Jarle Brinchmann (jarle@strw.leidenuniv.nl)

The background is removed from the wavelength region of interest before aggregating it.

marginfloat

The wavelength or pixel offset of the centers of the ranges of background images below and above the chosen wavelength range. This has the units specified by the unit_wave argument. The zero points of the margin are one pixel below and above the chosen wavelength range. The default value is 10 angstroms.

fbandfloat

The ratio of the number of images used to form a background image and the number of images that are being combined. The default value is 3.0.

methodstr

Name of the Cube method used to aggregate the data. This method must accept the axis=0 parameter and return an image. Example: mean, sum, max.

median_filterfloat

Size of the median filter for background estimation. If set to 0, the classical off band images are used.

Returns
Image
get_range(self, unit_wave=None, unit_wcs=None)[source]

Return the range of wavelengths, declinations and right ascensions in the cube.

The minimum and maximum coordinates are returned as an array in the following order:

[lbda_min, dec_min, ra_min, lbda_max, dec_max, ra_max]

Note that when the rotation angle of the image on the sky is not zero, dec_min, ra_min, dec_max and ra_max are not at the corners of the image.

Parameters
unit_waveastropy.units.Unit

The wavelengths units.

unit_wcsastropy.units.Unit

The angular units of the returned sky coordinates.

Returns
numpy.ndarray

The range of right ascensions declinations and wavelengths, arranged as [lbda_min, dec_min, ra_min, lbda_max, dec_max, ra_max].

get_rot(self, unit=Unit("deg"))[source]

Return the rotation angle of the images of the cube, defined such that a rotation angle of zero aligns north along the positive Y axis, and a positive rotation angle rotates north away from the Y axis, in the sense of a rotation from north to east.

Note that the rotation angle is defined in a flat map-projection of the sky. It is what would be seen if the pixels of the image were drawn with their pixel widths scaled by the angular pixel increments returned by the get_axis_increments() method.

Parameters
unitastropy.units.Unit

The unit to give the returned angle (degrees by default).

Returns
outfloat

The angle between celestial north and the Y axis of the image, in the sense of an eastward rotation of celestial north from the Y-axis.

get_start(self, unit_wave=None, unit_wcs=None)[source]

Return [lbda,y,x] at the center of pixel (0,0,0).

Parameters
unit_waveastropy.units.Unit

The wavelenth units to use for the starting wavelength. The default value, None, selects the intrinsic wavelength units of the cube.

unit_wcsastropy.units.Unit

The world coordinates units to use for the spatial starting position. The default value, None, selects the intrinsic world coordinates of the cube.

get_stat_hdu(self, name='STAT', header=None, convert_float32=True)

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_step(self, unit_wave=None, unit_wcs=None)[source]

Return the cube steps [dlbda,dy,dx].

Parameters
unit_waveastropy.units.Unit, optional

The wavelength units of the returned wavelength step.

unit_wcsastropy.units.Unit, optional

The angular units of the returned spatial world-coordinate steps.

Returns
numpy.ndarray

Returns [dlbda, dy, dx] where, dlbda is the size of pixels along the wavelength axis, and dy and dx are the sizes of pixels along the Y and X axes of the image, respectively.

get_wcs_header(self)

Return a FITS header containing coordinate descriptions.

info(self)

Print information.

loop_ima_multiprocessing(self, f, cpu=None, verbose=True, **kargs)[source]

Use multiple processes to run a function on each image of a cube.

The provided function must accept an Image object as its first argument, such as a method function of the Image class. There are three options for the return value of the function.

  1. The return value of the provided function can be another Image, in which case loop_ima_multiprocessing() returns a Cube of the processed images.

  2. Alternatively if each call to the function returns a scalar number, then these numbers are assembled into a Spectrum that has the same spectral coordinates as the Cube. This Cube is returned by loop_ima_multiprocessing().

  3. Finally, if each call to the provided function returns a value that can be stored as an element of a numpy array, then a numpy array of these values is returned, ordered such that element k of the array contains the return value for spectral channel k of the cube.

Parameters
fcallable() or Image method

The function to be applied to each image of the cube. This function must either be a method of the Image class, or it must be a top-level function that accepts an Image object as its first argument. It should return either an Image, a number, or a value that can be recorded as an element of a numpy array. Note that the function must not be a lambda function, or a function that is defined within another function.

cpuint

The desired number of CPUs to use, or None to select the number of available CPUs. By default, the available number of CPUs is equal to the number of CPU cores on the computer, minus 1 CPU for the main process. However the variable, mpdaf.CPU can be assigned a smaller number by the user to limit the number that are available to MPDAF.

verbosebool

If True, a progress report is printed every 5 seconds.

kargskargs

An optional list of arguments to be passed to the function f(). The datatypes of all of the arguments in this list must support python pickling.

Returns
out :

Cube if f returns Image, Spectrum if f returns a float or int. np.array(dtype=object) for all others cases.

loop_spe_multiprocessing(self, f, cpu=None, verbose=True, **kargs)[source]

Use multiple processes to run a function on each spectrum of a cube.

The provided function must accept a Spectrum object as its first argument, such as a method function of the Spectrum class. There are three options for the return value of the function.

  1. The return value of the provided function can be another Spectrum, in which case loop_spe_multiprocessing() returns a Cube of the processed spectra.

  2. Alternatively if each call to the function returns a scalar number, then these numbers are assembled into an Image that has the same world coordinates as the Cube. This Cube is returned by loop_spe_multiprocessing().

  3. Finally, if each call to the provided function returns a value that can be stored as an element of a numpy array, then a numpy array of these values is returned. This array has the shape of the images in the cube, such that pixel of [i,j] of the returned array contains the result of processing the spectrum of pixel [i,j] in the cube.

Parameters
fcallable() or Spectrum method

The function to be applied to each spectrum of the cube. This function must either be a method of the Spectrum class, or it must be a top-level function that accepts an Spectrum object as its first argument. It should return either an Spectrum, a number, or a value that can be recorded as an element of a numpy array. Note that the function must not be a lambda function, or a function that is defined within another function.

cpuint

The desired number of CPUs to use, or None to select the number of available CPUs. By default, the available number of CPUs is equal to the number of CPU cores on the computer, minus 1 CPU for the main process. However the variable, mpdaf.CPU can be assigned a smaller number by the user to limit the number that are available to MPDAF.

verbosebool

If True, a progress report is printed every 5 seconds.

kargskargs

An optional list of arguments to be passed to the function f(). The datatypes of all of the arguments in this list must support python pickling.

Returns
out :

Cube if f returns Spectrum, Image if f returns a float or int, np.array(dtype=object) for all others cases.

mask_ellipse(self, center, radius, posangle, lmin=None, lmax=None, inside=True, unit_center=Unit("deg"), unit_radius=Unit("arcsec"), unit_wave=Unit("Angstrom"))[source]

Mask values inside or outside an elliptical region.

Parameters
center(float,float)

Center (y,x) of the region, where y,x are usually celestial coordinates along the Y and X axes of the image, but are interpretted as Y,X array-indexes if unit_center is changed to None.

radius(float,float)

The radii of the two orthogonal axes of the ellipse. When posangle is zero, radius[0] is the radius along the X axis of the image-array, and radius[1] is the radius along the Y axis of the image-array.

posanglefloat

The counter-clockwise position angle of the ellipse in degrees. When posangle is zero, the X and Y axes of the ellipse are along the X and Y axes of the image.

lminfloat

The minimum wavelength of the region.

lmaxfloat

The maximum wavelength of the region.

insidebool

If inside is True, pixels inside the region are masked. If inside is False, pixels outside the region are masked.

unit_waveastropy.units.Unit

The units of the lmin and lmax wavelength coordinates (Angstroms by default). If None, the units of the lmin and lmax arguments are assumed to be pixels.

unit_centerastropy.units.Unit

The units of the coordinates of the center argument (degrees by default). If None, the units of the center argument are assumed to be pixels.

unit_radiusastropy.units.Unit

The units of the radius argument (arcseconds by default). If None, the units are assumed to be pixels.

mask_polygon(self, poly, lmin=None, lmax=None, unit_poly=Unit("deg"), unit_wave=Unit("Angstrom"), inside=True)[source]

Mask values inside or outside a polygonal region.

Parameters
poly(float, float)

An array of (float,float) containing a set of (p,q) or (dec,ra) values for the polygon vertices.

lminfloat

The minimum wavelength of the region.

lmaxfloat

The maximum wavelength of the region.

unit_polyastropy.units.Unit

The units of the polygon coordinates (degrees by default). Use unit_poly=None for polygon coordinates in pixels.

unit_waveastropy.units.Unit

The units of the wavelengths lmin and lmax (angstrom by default).

insidebool

If inside is True, pixels inside the polygonal region are masked. If inside is False, pixels outside the polygonal region are masked.

mask_region(self, center, radius, lmin=None, lmax=None, inside=True, unit_center=Unit("deg"), unit_radius=Unit("arcsec"), unit_wave=Unit("Angstrom"), posangle=0.0)[source]

Mask values inside or outside a circular or rectangular region.

Parameters
center(float,float)

Center (y,x) of the region, where y,x are usually celestial coordinates along the Y and X axes of the image, but are interpretted as Y,X array-indexes if unit_center is changed to None.

radiusfloat or (float,float)

The radius of a circular region, or the half-width and half-height of a rectangular region, respectively.

lminfloat

The minimum wavelength of the region.

lmaxfloat

The maximum wavelength of the region.

insidebool

If inside is True, pixels inside the region are masked. If inside is False, pixels outside the region are masked.

unit_waveastropy.units.Unit

The units of the lmin and lmax wavelength coordinates (Angstroms by default). If None, the units of the lmin and lmax arguments are assumed to be pixels.

unit_centerastropy.units.Unit

The units of the coordinates of the center argument (degrees by default). If None, the units of the center argument are assumed to be pixels.

unit_radiusastropy.units.Unit

The units of the radius argument (arcseconds by default). If None, the units are assumed to be pixels.

posanglefloat

When the region is rectangular, this is the counter-clockwise rotation angle of the rectangle in degrees. When posangle is 0.0 (the default), the X and Y axes of the rectangle are along the X and Y axes of the image.

mask_selection(self, ksel)

Mask selected pixels.

Parameters
kseloutput of np.where

Elements depending on a condition

mask_variance(self, threshold)

Mask pixels with a variance above a threshold value.

Parameters
thresholdfloat

Threshold value.

max(self, axis=None)[source]

Return the maximum over a given axis.

Beware that if the pixels of the cube have associated variances, these are discarded by this function, because there is no way to estimate the effects of a maximum on the variance.

Parameters
axisint or tuple of int, optional

The axis or axes along which the maximum is computed.

The default (axis = None) computes the maximum over all the dimensions of the cube and returns a float.

axis = 0 computes the maximum over the wavelength dimension and returns an image.

axis = (1,2) computes the maximum over the (X,Y) axes and returns a spectrum.

mean(self, axis=None, weights=None)[source]

Return a weighted or unweighted mean over a given axis or axes.

The mean is computed as follows. Note that weights of 1.0 are implicitly used for each data-point if the weights option is None. Given N data points of values, d[i], with weights, w[i], the weighted mean of d[0..N-1] is given by:

mean = Sum(d[i] * w[i]) / Sum(w[i])  for i=0..N-1

If data point d[i] has a variance of v[i], then the variance of the mean is given by:

variance = Sum(v[i] * w[i]**2) / Sum(w[i])**2   for i=0..N-1

Note that if one substitutes 1/v[i] for w[i] in this equation, the result is a variance of 1/Sum(1/v[i]). If all the variances, v[i], happen to have the same value, v, then this further simplifies to v / N, which is equivalent to a standard deviation of sigma/sqrt(N). This is the familiar result that the signal-to-noise ratio of a mean of N samples increases as sqrt(N).

Parameters
axisint or tuple of int, optional

The axis or axes along which the mean is to be performed.

The default (axis = None) performs a mean over all the dimensions of the cube and returns a float.

axis = 0 performs a mean over the wavelength dimension and returns an image.

axis = (1,2) performs a mean over the (X,Y) axes and returns a spectrum.

weightsarray_like, optional

When an array of weights is provided via this argument, it used to perform a weighted mean, as described in the introductory documentation above.

The weights array can have the same shape as the cube, or they can be 1-D if axis=(1,2), or 2-D if axis=0. If weights is None then all non-masked data points are given a weight equal to one.

If the Cube provides an array of variances for each data-point, then a good choice for the array of weights is the reciprocal of this array, (ie. weights=1.0/cube.var). However beware that not all data-sets provide variance information.

median(self, axis=None)[source]

Return the median over a given axis.

Beware that if the pixels of the cube have associated variances, these are discarded by this function, because there is no way to estimate the effects of a median on the variance.

Parameters
axisint or tuple of int, optional

The axis or axes along which a median is performed.

The default (axis = None) performs a median over all the dimensions of the cube and returns a float.

axis = 0 performs a median over the wavelength dimension and returns an image.

axis = (1,2) performs a median over the (X,Y) axes and returns a spectrum.

min(self, axis=None)[source]

Return the minimum over a given axis.

Beware that if the pixels of the cube have associated variances, these are discarded by this function, because there is no way to estimate the effects of a minimum on the variance.

Parameters
axisint or tuple of int, optional

The axis or axes along which the minimum is computed.

The default (axis = None) computes the minimum over all the dimensions of the cube and returns a float.

axis = 0 computes the minimum over the wavelength dimension and returns an image.

axis = (1,2) computes the minimum over the (X,Y) axes and returns a spectrum.

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

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

rebin(self, factor, margin='center', inplace=False)[source]

Combine neighboring pixels to reduce the size of a cube by integer factors along each axis.

Each output pixel is the mean of n pixels, where n is the product of the reduction factors in the factor argument.

Parameters
factorint or (int,int,int)

The integer reduction factors along the wavelength, z array axis, and the image y and x array axes, respectively. Python notation: (nz,ny,nx).

margin{‘center’, ‘right’, ‘left’, ‘origin’}

When the dimensions of the input cube are not integer multiples of the reduction factor, the cube is truncated to remove just enough pixels that its dimensions are multiples of the reduction factor. This subcube is then rebinned in place of the original cube. The margin parameter determines which pixels of the input cube are truncated, and which remain.

The options are:
‘origin’ or ‘center’:

The starts of the axes of the output cube are coincident with the starts of the axes of the input cube.

‘center’:

The center of the output cube is aligned with the center of the input cube, within one pixel along each axis.

‘right’:

The ends of the axes of the output cube are coincident with the ends of the axes of the input cube.

inplacebool

If False, return a rebinned copy of the cube (the default). If True, rebin the original cube in-place, and return that.

select_lambda(self, lbda_min, lbda_max=None, unit_wave=Unit("Angstrom"))[source]

Return the image or sub-cube corresponding to a wavelength range.

Parameters
lbda_minfloat

The minimum wavelength to be selected.

lbda_maxfloat, optional

The maximum wavelength to be selected. If not given, the closest image to lbda_min is returned.

unit_waveastropy.units.Unit

The wavelength units of lbda_min and lbda_max. The default unit is angstrom.

Returns
outmpdaf.obj.Cube or mpdaf.obj.Image

If more than one spectral channel is selected, then a Cube object is returned that contains just the images of those channels. If a single channel is selected, then an Image object is returned, containing just the image of that channel.

set_wcs(self, wcs=None, wave=None)

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.

spatial_erosion(self, npixels, inplace=False)[source]

Remove n pixels around the masked white image.

Parameters
npixelsint

Erosion width in pixels

inplacebool

If False (the default), return the results in a new Cube. If True, record the result in self and return that.

sqrt(self, out=None)

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.

subcube(self, center, size, lbda=None, unit_center=Unit("deg"), unit_size=Unit("arcsec"), unit_wave=Unit("Angstrom"))[source]

Return a subcube view around a position and for a wavelength range.

Note: as this is a view on the original cube, both the cube and the sub-cube will be modified at the same time. If you need to make changes only to the sub-cube, copy it before.

Parameters
center(float,float)

Center (y, x) of the aperture.

sizefloat

The size to extract. It corresponds to the size along the delta axis and the image is square.

lbda(float, float), optional

If not None, tuple giving the wavelength range.

unit_centerastropy.units.Unit

Type of the center coordinates (degrees by default)

unit_sizeastropy.units.Unit

unit of the size value (arcseconds by default)

unit_waveastropy.units.Unit

Wavelengths unit (angstrom by default) If None, inputs are in pixels

Returns
Cube
subcube_circle_aperture(self, center, radius, lbda=None, unit_center=Unit("deg"), unit_radius=Unit("arcsec"), unit_wave=Unit("Angstrom"))[source]

Extract a sub-cube that encloses a circular aperture of a specified radius and for a given wavelength range.

Pixels outside the circle are masked.

Parameters
center(float,float)

The center of the aperture (y,x)

radiusfloat

The radius of the aperture.

lbda(float, float) or None

If not None, tuple giving the wavelength range.

unit_centerastropy.units.Unit

The units of the center coordinates (degrees by default) The special value, None, indicates that the center is a 2D array index.

unit_radiusastropy.units.Unit

The units of the radius argument (arcseconds by default) The special value, None, indicates that the radius is specified in pixels.

unit_waveastropy.units.Unit

Wavelengths unit (angstrom by default) If None, inputs are in pixels

Returns
Cube
sum(self, axis=None, weights=None)[source]

Return a weighted or unweighted sum over a given axis or axes.

If the weights option is None, the variance is the sum of variances. Otherwise the weighted sum is computed with the Cube.mean method, whose result is multiplied by the number of elements.

Parameters
axisint or tuple of int, optional

Axis or axes along which a sum is performed:

  • The default (axis = None) performs a sum over all the dimensions of the cube and returns a float.

  • axis = 0 performs a sum over the wavelength dimension and returns an image.

  • axis = (1,2) performs a sum over the (X,Y) axes and returns a spectrum.

weightsarray_like, optional

When an array of weights is provided via this argument, it used to perform a weighted sum. This involves obtaining a weighted mean using the Cube.mean() function, then scaling this by the number of points that were averaged along the specified axes. The number of points that is used to scale the mean to a sum, is the total number of points along the averaged axes, not the number of unmasked points that had finite weights. As a result, the sum behaves as though all pixels along the averaged axes had values equal to the mean, regardless of whether any of these were masked or had zero weight.

The weights array can have the same shape as the cube, or they can be 1-D if axis=(1,2), or 2-D if axis=0. If weights=None, then all non-masked data points are given a weight equal to one. Finally, if a scalar float is given, then the data are all weighted equally. This can be used to get an unweighted sum that behaves as though masked pixels in the input cube had been filled with the mean along the averaged axes before the sum was performed.

If the Cube provides an array of variances for each data-point, then a good choice for the array of weights is the reciprocal of this array, (ie. weights=1.0/cube.var). However note that not all data-sets provide variance information, so check that cube.var is not None before trying this.

Any weight elements that are masked, infinite or nan, are replaced with zero. As a result, if the weights are specified as 1.0/cube.var, then any zero-valued variances will not produce infinite weights.

to_ds9(self, ds9id=None, newframe=False, zscale=False, cmap='grey')

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?

truncate(self, coord, mask=True, unit_wave=Unit("Angstrom"), unit_wcs=Unit("deg"))[source]

Return a sub-cube bounded by specified wavelength and spatial world-coordinates.

Note that unless unit_wcs is None, the y-axis and x-axis boundary coordinates are along sky axes such as declination and right-ascension, which may not be parallel to the image array-axes. When they are not parallel, the returned image area will contain some pixels that are outside the requested range. These are masked by default. To prevent them from being masked, pass False to the mask argument.

Parameters
coorditerable

A list of coordinate boundaries:

[lbda_min, y_min, x_min, lbda_max, y_max, x_max]

Note that this is the order of the values returned by mpdaf.obj.cube.get_range(), so when these functions are used together, then can be used to extract a subcube whose bounds match those of an existing smaller cube.

maskbool

If True, pixels outside [y_min,y_max] and [x_min,x_max] are masked. This can be useful when the world-coordinate X and Y axis are not parallel with the image array-axes.

unit_waveastropy.units.Unit

The wavelength units of lbda_min and lbda_max elements of the coord array. If None, lbda_min and lbda_max are interpretted as pixel indexes along the wavelength axis.

unit_wcsastropy.units.Unit

The wavelength units of x_min,x_max,y_min and y_max elements of the coord array. If None, these values are interpretted as pixel indexes along the image axes.

Returns
outmpdaf.obj.Cube

A Cube object that contains the requested sub-cube.

unmask(self)

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

write(self, filename, savemask='dq', checksum=False, convert_float32=True)

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.