Image

class mpdaf.obj.Image(filename=None, ext=None, wcs=None, data=None, var=None, unit=Unit(dimensionless), copy=True, dtype=None, **kwargs)[source]

Bases: ArithmeticMixin, DataArray

Manage image, optionally including a variance and a bad pixel mask.

Parameters:
filenamestr

Possible filename (.fits, .png or .bmp).

extint or (int,int) or string or (string,string)

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

wcsmpdaf.obj.WCS

World coordinates.

unitstr or astropy.units.Unit

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

datafloat array

Array containing the pixel values of the image. None by default.

varfloat array

Array containing the variance. None by default.

copybool

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

dtypenumpy.dtype

Type of the data (int, float)

Attributes:
filenamestr

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

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.

add_gaussian_noise(sigma[, interp])

Add Gaussian noise to image in place.

adjust_coordinates(ref[, nsigma, inplace])

Given a reference image of the sky that is expected to overlap with the current image, attempt to fit for any offset between the sky coordinate system of the current image and that of the reference image.

align_with_image(other[, flux, inplace, ...])

Resample the image to give it the same orientation, position, resolution and size as a given image.

background([niter, sigma])

Compute the image background with sigma-clipping.

clone([data_init, var_init])

Return a shallow copy with the same header and coordinates.

convolve(other[, inplace])

Convolve an Image with a 2D array or another Image, using the discrete convolution equation.

copy()

Return a new copy of an Image object.

correlate2d(other[, interp])

Return the cross-correlation of the image with an array/image

crop()

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

ee([center, radius, unit_center, ...])

Compute ensquared/encircled energy.

ee_size([center, unit_center, etot, frac, ...])

Compute the size of the square centered on (y,x) containing the fraction of the energy.

eer_curve([center, unit_center, ...])

Return containing enclosed energy as function of radius.

estimate_coordinate_offset(ref[, nsigma])

Given a reference image of the sky that is expected to overlap with the current image, attempt to fit for any offset between the sky coordinate system of the current image and that of the reference image.

fftconvolve(other[, inplace])

Convolve an Image with a 2D array or another Image, using the Fourier convolution theorem.

fftconvolve_gauss([center, flux, fwhm, ...])

Return the convolution of the image with a 2D gaussian.

fftconvolve_moffat([center, flux, a, q, n, ...])

Return the convolution of the image with a 2D moffat.

fwhm([center, radius, unit_center, unit_radius])

fwhm_gauss([center, radius, unit_center, ...])

Compute the Gaussian fwhm.

gauss_fit([pos_min, pos_max, center, flux, ...])

Perform Gaussian fit on image.

gaussian_filter([sigma, interp, inplace])

Return an image containing Gaussian filter applied to the current image.

get_axis_increments([unit])

Return the displacements on the sky that result from incrementing the array indexes of the image by one along the Y and X axes, respectively.

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_end([unit])

Return [y,x] corresponding to pixel (-1,-1).

get_range([unit])

Return the minimum and maximum right-ascensions and declinations in the image array.

get_rot([unit])

Return the rotation angle of the image, 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_spatial_fmax([rot])

Return the spatial-frequency band-limits of the image along the Y and X axes.

get_start([unit])

Return [y,x] corresponding to pixel (0,0).

get_stat_hdu([name, header, convert_float32])

Return an ImageHDU corresponding to the STAT extension.

get_step([unit])

Return the angular height and width of a pixel along the Y and X axes of the image array.

get_wcs_header()

Return a FITS header containing coordinate descriptions.

info()

Print information.

inside(coord[, unit])

Return True if coord is inside image.

mask_ellipse(center, radius, posangle[, ...])

Mask values inside or outside an elliptical region.

mask_polygon(poly[, unit, inside])

Mask values inside or outside a polygonal region.

mask_region(center, radius[, unit_center, ...])

Mask values inside or outside a circular or rectangular region.

mask_selection(ksel)

Mask selected pixels.

mask_variance(threshold)

Mask pixels with a variance above a threshold value.

moffat_fit([pos_min, pos_max, center, fwhm, ...])

Perform moffat fit on image.

moments([unit])

Return [width_y, width_x] first moments of the 2D gaussian.

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

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

norm([typ, value])

Normalize in place total flux to value (default 1).

peak([center, radius, unit_center, ...])

Find image peak location.

peak_detection(nstruct, niter[, threshold])

Return a list of peak locations.

plot([title, scale, vmin, vmax, zscale, ...])

Plot the image with axes labeled in pixels.

rebin(factor[, margin, inplace])

Combine neighboring pixels to reduce the size of an image by integer factors along each axis.

regrid(newdim, refpos, refpix, newinc[, ...])

Resample an image of the sky to select its angular resolution, to specify the position of the sky in the image array, and optionally to reflect one or more of its axes.

resample(newdim, newstart, newstep[, flux, ...])

Resample an image of the sky to select its angular resolution and to specify which sky position appears at the center of pixel [0,0].

rotate([theta, interp, reshape, order, ...])

Rotate the sky within an image in the sense of a rotation from north to east.

segment([shape, minsize, minpts, ...])

Segment the image in a number of smaller images.

set_spatial_fmax([newfmax, rot])

Specify the spatial-frequency band-limits of the image along the Y and X axis.

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.

subimage(center, size[, unit_center, ...])

Return a view on a square or rectangular part.

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

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

truncate(y_min, y_max, x_min, x_max[, mask, ...])

Return a sub-image that contains a specified area of the sky.

unmask()

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

update_spatial_fmax(newfmax[, rot])

Update the spatial-frequency band-limits recorded for the current image.

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)

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.

add_gaussian_noise(sigma, interp='no')[source]

Add Gaussian noise to image in place.

Parameters:
sigmafloat

Standard deviation.

interp‘no’ | ‘linear’ | ‘spline’

if ‘no’, data median value replaced masked values. if ‘linear’, linear interpolation of the masked values. if ‘spline’, spline interpolation of the masked values.

adjust_coordinates(ref, nsigma=1.0, inplace=False)[source]

Given a reference image of the sky that is expected to overlap with the current image, attempt to fit for any offset between the sky coordinate system of the current image and that of the reference image. Apply this offset to the coordinates of the current image, to bring it into line with the reference image.

This function calls self.estimate_coordinate_offset() to fit for the offset between the coordinate systems of the two images, then adjusts the coordinate reference pixel of the current image to bring its coordinates into line with those of the reference image.

Parameters:
refImage

The image of the sky that is to be used as the coordinate reference. The sky coverage of this image should overlap with that of self. Ideally the resolution of this image should be at least as good as the resolution of self.

nsigmafloat

Only values that exceed this many standard deviations above the mean of each image will be used.

inplacebool

If False, return a shifted copy of the image (the default). If True, shift the original image in-place, and return that.

Returns:
outImage

A version of self in which the sky coordinates have been shifted to match those of the reference image.

align_with_image(other, flux=False, inplace=False, cutoff=0.25, antialias=True, window='blackman')[source]

Resample the image to give it the same orientation, position, resolution and size as a given image.

The image is first rotated to give it the same orientation on the sky as the other image. The resampling process also eliminates any shear terms from the original image, so that its pixels can be correctly drawn on a rectangular grid.

Secondly the image is resampled. This changes its resolution, shifts the image such that the same points on the sky appear in the same pixels as in the other image, and changes the dimensions of the image array to match that of the other image.

The rotation and resampling processes are performed as separate steps because the anti-aliasing filter that needs to be applied in the resampling step reduces the resolution, is difficult to implement before the axes have been rotated to the final orientation.

Parameters:
otherImage

The image to be aligned with.

fluxbool

This tells the function whether the pixel units of the image are flux densities (flux=True), such as erg/s/cm2/Hz, or whether they are per-steradian brightness units (flux=False), such as erg/s/cm2/Hz/steradian. It needs to know this when it changes the pixel size, because when pixel sizes change, resampled flux densities need to be corrected for the change in the area per pixel, where resampled brightnesses don’t.

inplacebool

If False, return an aligned copy of the image (the default). If True, align the original image in-place, and return that.

cutofffloat

Mask each output pixel where at least this fraction of the pixel was interpolated from dummy values given to masked input pixels.

antialiasbool

By default, when the resolution of an image axis is about to be reduced, a low pass filter is first applied to suppress high spatial frequencies that can not be represented by the reduced sampling interval. If this is not done, high-frequency noise and sharp edges get folded back to lower frequencies, where they increase the noise level of the image and introduce ringing artefacts next to sharp edges, such as CCD saturation spikes and bright unresolved stars. This filtering can be disabled by passing False to the antialias argument.

windowstr

The type of window function to use for antialiasing in the Fourier plane. The following windows are supported:

blackman

This window suppresses ringing better than any other window, at the expense of lowered image resolution. In the image plane, the PSF of this window is approximately gaussian, with a standard deviation of around 0.96*newstep, and a FWHM of about 2.3*newstep.

gaussian

A truncated gaussian window. This has a smaller PSF than the blackman window, however gaussians never fall to zero, so either significant ringing will be seen due to truncation of the gaussian, or low-level aliasing will occur, depending on the spatial frequency coverage of the image beyond the folding frequency. It can be a good choice for images that only contain smoothly varying features. It is equivalent to a convolution of the image with both an airy profile and a gaussian of standard deviation 0.724*newstep (FWHM 1.704*newstep).

rectangle

This window simply zeros all spatial frequencies above the highest that can be correctly sampled by the new pixel size. This gives the best resolution of any of the windows, but this is marred by the strong sidelobes of the resulting airy-profile, especially near bright point sources and CCD saturation lines.

background(niter=3, sigma=3.0)[source]

Compute the image background with sigma-clipping.

Returns the background value and its standard deviation.

Parameters:
niterint

Number of iterations.

sigmafloat

Number of sigma used for the clipping.

Returns:
out2-dim float array
clone(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(other, inplace=False)[source]

Convolve an Image with a 2D array or another Image, using the discrete convolution equation.

This function, which uses the discrete convolution equation, is usually slower than Image.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:
otherImage or numpy.ndarray

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

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

(other.shape - 1) // 2

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

inplacebool

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

Returns:
outImage
copy()[source]

Return a new copy of an Image object.

correlate2d(other, interp='no')[source]

Return the cross-correlation of the image with an array/image

Uses scipy.signal.correlate2d.

Parameters:
other2d-array or Image

Second Image or 2d-array.

interp‘no’ | ‘linear’ | ‘spline’

if ‘no’, data median value replaced masked values. if ‘linear’, linear interpolation of the masked values. if ‘spline’, spline interpolation of the masked values.

crop()

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.

ee(center=None, radius=0, unit_center=Unit('deg'), unit_radius=Unit('arcsec'), frac=False, cont=0)[source]

Compute ensquared/encircled energy.

Parameters:
center(float,float)

Center of the explored region. If center is None, the full image is explored.

radiusfloat or (float,float)

Radius defined the explored region. If float, it defined a circular region (encircled energy). If (float,float), it defined a rectangular region (ensquared energy).

unit_centerastropy.units.Unit

Type of the center coordinates. Degrees by default (use None for coordinates in pixels).

unit_radiusastropy.units.Unit

Radius unit. Arcseconds by default (use None for radius in pixels)

fracbool

If frac is True, result is given relative to the total energy of the full image.

contfloat

Continuum value.

Returns:
outfloat

Ensquared/encircled flux.

ee_size(center=None, unit_center=Unit('deg'), etot=None, frac=0.9, cont=0, unit_size=Unit('arcsec'))[source]

Compute the size of the square centered on (y,x) containing the fraction of the energy.

Parameters:
center(float,float)

Center (y,x) of the explored region. If center is None, center of the image is used.

unitastropy.units.Unit

Type of the center coordinates. Degrees by default (use None for coordinates in pixels).

etotfloat
Total energy used to comute the ratio.

If etot is not set, it is computed from the full image.

fracfloat in ]0,1]

Fraction of energy.

contfloat

continuum value

unit_centerastropy.units.Unit

Type of the center coordinates. Degrees by default (use None for coordinates in pixels).

unit_sizeastropy.units.Unit

Size unit. Arcseconds by default (use None for sier in pixels).

Returns:
outfloat array
eer_curve(center=None, unit_center=Unit('deg'), unit_radius=Unit('arcsec'), etot=None, cont=0)[source]

Return containing enclosed energy as function of radius.

The enclosed energy ratio (EER) shows how much light is concentrated within a certain radius around the image-center.

Parameters:
center(float,float)

Center of the explored region. If center is None, center of the image is used.

unit_centerastropy.units.Unit

Type of the center coordinates. Degrees by default (use None for coordinates in pixels).

unit_radiusastropy.units.Unit

Radius units (arcseconds by default)/

etotfloat

Total energy used to comute the ratio. If etot is not set, it is computed from the full image.

contfloat

Continuum value.

Returns:
out(float array, float array)

Radius array, EER array

estimate_coordinate_offset(ref, nsigma=1.0)[source]

Given a reference image of the sky that is expected to overlap with the current image, attempt to fit for any offset between the sky coordinate system of the current image and that of the reference image. The returned value is designed to be added to the coordinate reference pixel values of self.wcs.

This function performs the following steps:

  1. The align_with_image() method is called to resample the reference image onto the same coordinate grid as the current image.

  2. The two images are then cross-correlated, after zeroing all background values in the images below nsigma standard deviations above the mean.

  3. The peak in the auto-correlation image is found and its sub-pixel position is estimated by a simple quadratic interpolation. This position, relative to the center of the auto-correlation image, gives the average position offset between similar features in the two images.

Parameters:
refImage

The image of the sky that is to be used as the coordinate reference. The sky coverage of this image should overlap with that of self. Ideally the resolution of this image should be at least as good as the resolution of self.

nsigmafloat

Only values that exceed this many standard deviations above the mean of each image will be used.

Returns:
outfloat,float

The pixel offsets that would need to be added to the coordinate reference pixel values, crpix2 and crpix1, of self.wcs to make the features in self line up with those in the reference image.

fftconvolve(other, inplace=False)[source]

Convolve an Image with a 2D array or another Image, using the Fourier convolution theorem.

This function, which performs the convolution by multiplying the Fourier transforms of the two images, is usually much faster than Image.convolve(), except when other.data.size is small. However it uses much more memory, so Image.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, Image.convolve() may be more efficient than Image.fftconvolve().

Uses scipy.signal.fftconvolve.

Parameters:
otherImage or numpy.ndarray

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

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

(other.shape - 1) // 2

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

inplacebool

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

Returns:
outImage
fftconvolve_gauss(center=None, flux=1.0, fwhm=(1.0, 1.0), peak=False, rot=0.0, factor=1, unit_center=Unit('deg'), unit_fwhm=Unit('arcsec'), inplace=False)[source]

Return the convolution of the image with a 2D gaussian.

Parameters:
center(float,float)

Gaussian center (y_peak, x_peak). If None the center of the image is used. The unit is given by the unit_center parameter (degrees by default).

fluxfloat

Integrated gaussian flux or gaussian peak value if peak is True.

fwhm(float,float)

Gaussian fwhm (fwhm_y,fwhm_x). The unit is given by the unit_fwhm parameter (arcseconds by default).

peakbool

If true, flux contains a gaussian peak value.

rotfloat

Angle position in degree.

factorint

If factor<=1, gaussian value is computed in the center of each pixel. If factor>1, for each pixel, gaussian value is the sum of the gaussian values on the factor*factor pixels divided by the pixel area.

unit_centerastropy.units.Unit

type of the center and position coordinates. Degrees by default (use None for coordinates in pixels).

unit_fwhmastropy.units.Unit

FWHM unit. Arcseconds by default (use None for radius in pixels)

inplacebool

If False, return a convolved copy of the image (default value). If True, convolve the original image in-place, and return that.

Returns:
outImage
fftconvolve_moffat(center=None, flux=1.0, a=1.0, q=1.0, n=2, peak=False, rot=0.0, factor=1, unit_center=Unit('deg'), unit_a=Unit('arcsec'), inplace=False)[source]

Return the convolution of the image with a 2D moffat.

Parameters:
center(float,float)

Gaussian center (y_peak, x_peak). If None the center of the image is used. The unit is given by the unit_center parameter (degrees by default).

fluxfloat

Integrated gaussian flux or gaussian peak value if peak is True.

afloat

Half width at half maximum of the image in the absence of atmospheric scattering. 1 by default. The unit is given by the unit_a parameter (arcseconds by default).

qfloat

Axis ratio, 1 by default.

nint

Atmospheric scattering coefficient. 2 by default.

rotfloat

Angle position in degree.

factorint

If factor<=1, moffat value is computed in the center of each pixel. If factor>1, for each pixel, moffat value is the sum of the moffat values on the factor*factor pixels divided by the pixel area.

peakbool

If true, flux contains a gaussian peak value.

unit_centerastropy.units.Unit

type of the center and position coordinates. Degrees by default (use None for coordinates in pixels).

unit_aastropy.units.Unit

a unit. Arcseconds by default (use None for radius in pixels)

inplacebool

If False, return a convolved copy of the image (default value). If True, convolve the original image in-place, and return that.

Returns:
outImage
fwhm(center=None, radius=0, unit_center=Unit('deg'), unit_radius=Unit('arcsec'))[source]
fwhm_gauss(center=None, radius=0, unit_center=Unit('deg'), unit_radius=Unit('arcsec'))[source]

Compute the Gaussian fwhm.

Parameters:
center(float,float)

Center of the explored region. If center is None, the full image is explored.

radiusfloat or (float,float)

Radius defined the explored region.

unit_centerastropy.units.Unit

type of the center coordinates. Degrees by default (use None for coordinates in pixels).

unit_radiusastropy.units.Unit

Radius unit. Arcseconds by default (use None for radius in pixels)

Returns:
outarray of float

[fwhm_y,fwhm_x], returned in unit_radius (arcseconds by default).

gauss_fit(pos_min=None, pos_max=None, center=None, flux=None, fwhm=None, circular=False, cont=0, fit_back=True, rot=0, peak=False, factor=1, weight=True, plot=False, unit_center=Unit('deg'), unit_fwhm=Unit('arcsec'), maxiter=0, verbose=True, full_output=0)[source]

Perform Gaussian fit on image.

Parameters:
pos_min(float,float)

Minimum y and x values. Their unit is given by the unit_center parameter (degrees by default).

pos_max(float,float)

Maximum y and x values. Their unit is given by the unit_center parameter (degrees by default).

center(float,float)

Initial gaussian center (y_peak,x_peak) If None it is estimated. The unit is given by the unit_center parameter (degrees by default).

fluxfloat

Initial integrated gaussian flux or gaussian peak value if peak is True. If None, peak value is estimated.

fwhm(float,float)

Initial gaussian fwhm (fwhm_y,fwhm_x). If None, they are estimated. The unit is given by unit_fwhm (arcseconds by default).

circularbool

True: circular gaussian, False: elliptical gaussian

contfloat

continuum value, 0 by default.

fit_backbool

False: continuum value is fixed, True: continuum value is a fit parameter.

rotfloat

Initial rotation in degree. If None, rotation is fixed to 0.

peakbool

If true, flux contains a gaussian peak value.

factorint

If factor<=1, gaussian value is computed in the center of each pixel. If factor>1, for each pixel, gaussian value is the sum of the gaussian values on the factor*factor pixels divided by the pixel area.

weightbool

If weight is True, the weight is computed as the inverse of variance.

unit_centerastropy.units.Unit

type of the center and position coordinates. Degrees by default (use None for coordinates in pixels).

unit_fwhmastropy.units.Unit

FWHM unit. Arcseconds by default (use None for radius in pixels)

maxiterint

The maximum number of iterations during the sum of square minimization.

plotbool

If True, the gaussian is plotted.

verbosebool

If True, the Gaussian parameters are printed at the end of the method.

full_outputbool

True-zero to return a mpdaf.obj.Gauss2D object containing the gauss image.

Returns:
outmpdaf.obj.Gauss2D
gaussian_filter(sigma=3, interp='no', inplace=False)[source]

Return an image containing Gaussian filter applied to the current image.

Uses scipy.ndimage.gaussian_filter.

Parameters:
sigmafloat

Standard deviation for Gaussian kernel

interp‘no’ | ‘linear’ | ‘spline’

if ‘no’, data median value replaced masked values. if ‘linear’, linear interpolation of the masked values. if ‘spline’, spline interpolation of the masked values.

inplacebool

If False, return a filtered copy of the image (the default). If True, filter the original image in-place, and return that.

Returns:
outImage
get_axis_increments(unit=None)[source]

Return the displacements on the sky that result from incrementing the array indexes of the image by one along the Y and X axes, respectively.

In MPDAF, images are sampled on a regular grid of square pixels that represent a flat projection of the celestial sphere. The get_axis_increments() method returns the angular width and height of these pixels on the sky, with signs that indicate whether the angle increases or decreases as one increments the array indexes. To keep plots consistent, regardless of the rotation angle of the image on the sky, the returned height is always positive, but the returned width is negative if a plot of the image with pixel 0,0 at the bottom left would place east anticlockwise of north, and positive otherwise.

Parameters:
unitastropy.units.Unit

The angular units of the returned values.

Returns:
outnumpy.ndarray

(dy,dx). These are the angular increments of pixels along the Y and X axes of the image. The returned values are either in the unit specified by the ‘unit’ input parameter, or in the unit specified by the self.unit property.

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

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

get_end(unit=None)[source]

Return [y,x] corresponding to pixel (-1,-1).

Parameters:
unitastropy.units.Unit

type of the world coordinates

Returns:
outfloat array
get_range(unit=None)[source]

Return the minimum and maximum right-ascensions and declinations in the image array.

Specifically a list is returned with the following contents:

[dec_min, ra_min, dec_max, ra_max]

Note that if the Y axis of the image is not parallel to the declination axis, then the 4 returned values will all come from different corners of the image. In particular, note that this means that the coordinates [dec_min,ra_min] and [dec_max,ra_max] will only coincide with pixels in the image if the Y axis is aligned with the declination axis. Otherwise they will be outside the bounds of the image.

Parameters:
unitastropy.units.Unit

The units of the returned angles.

Returns:
outnumpy.ndarray

The range of right ascensions and declinations, arranged as [dec_min, ra_min, dec_max, ra_max]. The returned values are either in the units specified in the ‘unit’ input parameter, or in the units stored in the self.unit property.

get_rot(unit=Unit('deg'))[source]

Return the rotation angle of the image, 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_spatial_fmax(rot=None)[source]

Return the spatial-frequency band-limits of the image along the Y and X axes.

See the documentation of set_spatial_fmax() for an explanation of what the band-limits are used for.

If no band limits have been specified yet, this function has the side-effect of setting them to the band-limits dictated by the sampling interval of the image array. Specifically, an X axis with a sampling interval of dx can sample spatial frequencies of up to 0.5/dx cycles per unit of dx without aliasing.

Parameters:
rotfloat or None

Either None, to request band-limits that pertain to the Y and X axes of the current image without any rotation, or, if the band-limits pertain to a rotated version of the image, the rotation angle of its Y axis westward of north (degrees). This is defined such that if image.wcs.get_rot() is passed to this function, the band limits for the Y and X axes of the current image axes will be returned.

Returns:
outnumpy.ndarray

The spatial-frequency band-limits of the image along the Y and X axes of the image in cycles per self.wcs.unit.

get_start(unit=None)[source]

Return [y,x] corresponding to pixel (0,0).

Parameters:
unitastropy.units.Unit

type of the world coordinates

Returns:
outfloat array
get_stat_hdu(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(unit=None)[source]

Return the angular height and width of a pixel along the Y and X axes of the image array.

In MPDAF, images are sampled on a regular grid of square pixels that represent a flat projection of the celestial sphere. The get_step() method returns the angular width and height of these pixels on the sky.

See also get_axis_increments().

Parameters:
unitastropy.units.Unit

The angular units of the returned values.

Returns:
outnumpy.ndarray

(dy,dx). These are the angular height and width of pixels along the Y and X axes of the image. The returned values are either in the unit specified by the ‘unit’ input parameter, or in the unit specified by the self.unit property.

get_wcs_header()

Return a FITS header containing coordinate descriptions.

info()

Print information.

inside(coord, unit=Unit('deg'))[source]

Return True if coord is inside image.

Parameters:
coord(float,float)

coordinates (y,x).

unitastropy.units.Unit

Type of the coordinates (degrees by default)

Returns:
outbool
mask_ellipse(center, radius, posangle, unit_center=Unit('deg'), unit_radius=Unit('arcsec'), inside=True)[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 rotation 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.

unit_centerastropy.units.Unit

The units of the center coordinates. Degrees by default (use None for coordinates in pixels).

unit_radiusastropy.units.Unit

The units of the radius argument. Arcseconds by default. (use None for radius in pixels)

insidebool

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

mask_polygon(poly, unit=Unit('deg'), 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.

unitastropy.units.Unit

The units of the polygon coordinates (by default in degrees). Use unit=None to have polygon coordinates in pixels.

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(center, radius, unit_center=Unit('deg'), unit_radius=Unit('arcsec'), inside=True, 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.

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.

insidebool

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

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 ellipse are along the X and Y axes of the image.

mask_selection(ksel)

Mask selected pixels.

Parameters:
kseloutput of np.where

Elements depending on a condition

mask_variance(threshold)

Mask pixels with a variance above a threshold value.

Parameters:
thresholdfloat

Threshold value.

moffat_fit(pos_min=None, pos_max=None, center=None, fwhm=None, flux=None, n=2.0, circular=False, cont=0, fit_back=True, rot=0, peak=False, factor=1, weight=True, plot=False, unit_center=Unit('deg'), unit_fwhm=Unit('arcsec'), verbose=True, full_output=0, fit_n=True, maxiter=0)[source]

Perform moffat fit on image.

Parameters:
pos_min(float,float)

Minimum y and x values. Their unit is given by the unit_center parameter (degrees by default).

pos_max(float,float)

Maximum y and x values. Their unit is given by the unit_center parameter (degrees by default).

center(float,float)

Initial moffat center (y_peak,x_peak). If None it is estimated. The unit is given by the unit_center parameter (degrees by default).

fluxfloat

Initial integrated gaussian flux or gaussian peak value if peak is True. If None, peak value is estimated.

fwhm(float,float)

Initial gaussian fwhm (fwhm_y,fwhm_x). If None, they are estimated. Their unit is given by the unit_fwhm parameter (arcseconds by default).

nint

Initial atmospheric scattering coefficient.

circularbool

True: circular moffat, False: elliptical moffat

contfloat

continuum value, 0 by default.

fit_backbool

False: continuum value is fixed, True: continuum value is a fit parameter.

rotfloat

Initial angle position in degree.

peakbool

If true, flux contains a gaussian peak value.

factorint

If factor<=1, gaussian is computed in the center of each pixel. If factor>1, for each pixel, gaussian value is the sum of the gaussian values on the factor*factor pixels divided by the pixel area.

weightbool

If weight is True, the weight is computed as the inverse of variance.

plotbool

If True, the gaussian is plotted.

unit_centerastropy.units.Unit

type of the center and position coordinates. Degrees by default (use None for coordinates in pixels).

unit_fwhmastropy.units.Unit

FWHM unit. Arcseconds by default (use None for radius in pixels)

full_outputbool

True to return a mpdaf.obj.Moffat2D object containing the moffat image.

fit_nbool

False: n value is fixed, True: n value is a fit parameter.

maxiterint

The maximum number of iterations during the sum of square minimization.

Returns:
outmpdaf.obj.Moffat2D
moments(unit=Unit('arcsec'))[source]

Return [width_y, width_x] first moments of the 2D gaussian.

Parameters:
unitastropy.units.Unit

Unit of the returned moments (arcseconds by default). If None, moments will be in pixels.

Returns:
outfloat array
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).

norm(typ='flux', value=1.0)[source]

Normalize in place total flux to value (default 1).

Parameters:
type‘flux’ | ‘sum’ | ‘max’

If ‘flux’,the flux is normalized and the pixel area is taken into account.

If ‘sum’, the flux is normalized to the sum of flux independently of pixel size.

If ‘max’, the flux is normalized so that the maximum of intensity will be ‘value’.

valuefloat

Normalized value (default 1).

peak(center=None, radius=0, unit_center=Unit('deg'), unit_radius=Unit('arcsec'), dpix=2, background=None, plot=False)[source]

Find image peak location.

Used scipy.ndimage.measurements.maximum_position and scipy.ndimage.measurements.center_of_mass.

Parameters:
center(float,float)

Center (y,x) of the explored region. If center is None, the full image is explored.

radiusfloat or (float,float)

Radius defined the explored region.

unit_centerastropy.units.Unit

Type of the center coordinates. Degrees by default (use None for coordinates in pixels).

unit_radiusastropy.units.Unit

Radius unit. Arcseconds by default (use None for radius in pixels)

dpixint

Half size of the window (in pixels) to compute the center of gravity.

backgroundfloat

Background value. If None, it is computed.

plotbool

If True, the peak center is overplotted on the image.

Returns:
outdict {‘y’, ‘x’, ‘p’, ‘q’, ‘data’}

Containing the peak position and the peak intensity.

peak_detection(nstruct, niter, threshold=None)[source]

Return a list of peak locations.

Parameters:
nstructint

Size of the structuring element used for the erosion.

niterint

Number of iterations used for the erosion and the dilatation.

thresholdfloat

Threshold value. If None, it is initialized with background value.

Returns:
outnp.array
plot(title=None, scale='linear', vmin=None, vmax=None, zscale=False, colorbar=None, var=False, show_xlabel=False, show_ylabel=False, ax=None, unit=Unit('deg'), use_wcs=False, **kwargs)[source]

Plot the image with axes labeled in pixels.

If either axis has just one pixel, plot a line instead of an image.

Colors are assigned to each pixel value as follows. First each pixel value, pv, is normalized over the range vmin to vmax, to have a value nv, that goes from 0 to 1, as follows:

nv = (pv - vmin) / (vmax - vmin)

This value is then mapped to another number between 0 and 1 which determines a position along the colorbar, and thus the color to give the displayed pixel. The mapping from normalized values to colorbar position, color, can be chosen using the scale argument, from the following options:

  • ‘linear’: color = nv

  • ‘log’: color = log(1000 * nv + 1) / log(1000 + 1)

  • ‘sqrt’: color = sqrt(nv)

  • ‘arcsinh’: color = arcsinh(10*nv) / arcsinh(10.0)

A colorbar can optionally be drawn. If the colorbar argument is given the value ‘h’, then a colorbar is drawn horizontally, above the plot. If it is ‘v’, the colorbar is drawn vertically, to the right of the plot.

By default the image is displayed in its own plot. Alternatively to make it a subplot of a larger figure, a suitable matplotlib.axes.Axes object can be passed via the ax argument. Note that unless matplotlib interative mode has previously been enabled by calling matplotlib.pyplot.ion(), the plot window will not appear until the next time that matplotlib.pyplot.show() is called. So to arrange that a new window appears as soon as Image.plot() is called, do the following before the first call to Image.plot():

import matplotlib.pyplot as plt
plt.ion()
Parameters:
titlestr

An optional title for the figure (None by default).

scale‘linear’ | ‘log’ | ‘sqrt’ | ‘arcsinh’

The stretch function to use mapping pixel values to colors (The default is ‘linear’). The pixel values are first normalized to range from 0 for values <= vmin, to 1 for values >= vmax, then the stretch algorithm maps these normalized values, nv, to a position p from 0 to 1 along the colorbar, as follows: linear: p = nv log: p = log(1000 * nv + 1) / log(1000 + 1) sqrt: p = sqrt(nv) arcsinh: p = arcsinh(10*nv) / arcsinh(10.0)

vminfloat

Pixels that have values <= vmin are given the color at the dark end of the color bar. Pixel values between vmin and vmax are given colors along the colorbar according to the mapping algorithm specified by the scale argument.

vmaxfloat

Pixels that have values >= vmax are given the color at the bright end of the color bar. If None, vmax is set to the maximum pixel value in the image.

zscalebool

If True, vmin and vmax are automatically computed using the IRAF zscale algorithm.

colorbarstr

If ‘h’, a horizontal colorbar is drawn above the image. If ‘v’, a vertical colorbar is drawn to the right of the image. If None (the default), no colorbar is drawn.

varbool

If true variance array is shown in place of data array

axmatplotlib.axes.Axes

An optional Axes instance in which to draw the image, or None to have one created using matplotlib.pyplot.gca().

unitastropy.units.Unit

The units to use for displaying world coordinates (degrees by default). In the interactive plot, when the mouse pointer is over a pixel in the image the coordinates of the pixel are shown using these units, along with the pixel value.

use_wcsbool

If True, use astropy.visualization.wcsaxes to get axes with world coordinates.

kwargsmatplotlib.artist.Artist

Optional extra keyword/value arguments to be passed to the ax.imshow() function.

Returns:
outmatplotlib AxesImage
rebin(factor, margin='center', inplace=False)[source]

Combine neighboring pixels to reduce the size of an image 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)

The integer reduction factor along the y and x array axes. Note the conventional python ordering of the axes.

margin‘center’|’right’|’left’|’origin’

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

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

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

‘center’:

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

‘right’:

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

inplacebool

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

Returns:
outImage
regrid(newdim, refpos, refpix, newinc, flux=False, order=1, interp='no', unit_pos=Unit('deg'), unit_inc=Unit('arcsec'), antialias=True, inplace=False, cutoff=0.25, window='blackman')[source]

Resample an image of the sky to select its angular resolution, to specify the position of the sky in the image array, and optionally to reflect one or more of its axes.

This function can be used to decrease or increase the resolution of an image. It can also shift the contents of an image to place a specific (dec,ra) position at a specific fractional pixel position. Finally, it can be used to invert the direction of one or both of the array axes on the sky.

When this function is used to resample an image to a lower resolution, a low-pass anti-aliasing filter is applied to the image before it is resampled, to remove all spatial frequencies below half the new sampling rate. This is required to satisfy the Nyquist sampling constraint. It prevents high spatial-frequency noise and edges from being aliased to lower frequency artefacts in the resampled image. The removal of this noise improves the signal to noise ratio of the resampled image.

Parameters:
newdimint or (int,int)

The desired new dimensions. Python notation: (ny,nx)

refpos(float, float)

The sky position (dec,ra) to place at the pixel specified by the refpix argument.

If refpix and refpos are both None, the sky position at the bottom corner of the input image is placed at the bottom left corner of the output image. Note that refpix and refpos must either both be given values, or both be None.

refpix(float, float)

The [Y, X] indexes of the output pixel where the sky position, refpos, should be placed. Y and X are interpreted as floating point indexes, where integer values indicate pixel centers and integer values +/- 0.5 indicate the edges of pixels.

If refpix and refpos are both None, the sky position at the bottom corner of the input image is placed at the bottom left corner of the output image. Note that refpix and refpos must either both be given values, or both be None.

newincfloat or (float, float)

The signed increments of the angle on the sky from one pixel to the next, given as either a single increment for both image axes, or two numbers (dy,dx) for the Y and X axes respectively.

The signs of these increments are interpreted as described in the documentation of the Image.get_axis_increments() function. In particular, note that dy is typically positive and dx is usually negative, such that when the image is plotted, east appears anticlockwise of north, and east is towards the left of the plot when the image rotation angle is zero.

If either of the signs of the two newinc numbers is different from the sign of the increments of the original image (queryable with image.get_axis_increments()), then the image will be reflected about that axis. In this case the value of the refpix argument should be chosen with care, because otherwise the sampled part of the image may end up being reflected outside the limits of the image array, and the result will be a blank image.

If only one number is given for newinc then both axes are given the same resolution, but the signs of the increments are kept the same as the pixel increments of the original image.

fluxbool

This tells the function whether the pixel units of the image are flux densities (flux=True), such as erg/s/cm2/Hz, or whether they are per-steradian brightness units (flux=False), such as erg/s/cm2/Hz/steradian. It needs to know this when it changes the pixel size, because when pixel sizes change, resampled flux densities need to be corrected for the change in the area per pixel, where resampled brightnesses don’t.

orderint

The order of the spline interpolation. This can take any value from 0-5. The default is 1 (linear interpolation). When this function is used to lower the resolution of an image, the low-pass anti-aliasing filter that is applied, makes linear interpolation sufficient. Conversely, when this function is used to increase the image resolution, order=3 might be useful. Higher orders than this will tend to introduce ringing artefacts.

interp‘no’ | ‘linear’ | ‘spline’

If ‘no’, replace masked data with the median image value. If ‘linear’, replace masked values using a linear interpolation between neighboring values. if ‘spline’, replace masked values using a spline interpolation between neighboring values.

unit_posastropy.units.Unit

The units of the refpos coordinates. Degrees by default.

unit_incastropy.units.Unit

The units of newinc. Arcseconds by default.

antialiasbool

By default, when the resolution of an image axis is about to be reduced, a low pass filter is first applied to suppress high spatial frequencies that can not be represented by the reduced sampling interval. If this is not done, high-frequency noise and sharp edges get folded back to lower frequencies, where they increase the noise level of the image and introduce ringing artefacts next to sharp edges, such as CCD saturation spikes. This filtering can be disabled by passing False to the antialias argument.

inplacebool

If False, return a resampled copy of the image (the default). If True, resample the original image in-place, and return that.

cutofffloat

Mask each output pixel where at least this fraction of the pixel was interpolated from dummy values given to masked input pixels.

windowstr

The type of window function to use for antialiasing in the Fourier plane. The following windows are supported:

blackman

This window suppresses ringing better than any other window, at the expense of lowered image resolution. In the image plane, the PSF of this window is approximately gaussian, with a standard deviation of around 0.96*newstep, and a FWHM of about 2.3*newstep.

gaussian

A truncated gaussian window. This has a smaller PSF than the blackman window, however gaussians never fall to zero, so either significant ringing will be seen due to truncation of the gaussian, or low-level aliasing will occur, depending on the spatial frequency coverage of the image beyond the folding frequency. It can be a good choice for images that only contain smoothly varying features. It is equivalent to a convolution of the image with both an airy profile and a gaussian of standard deviation 0.724*newstep (FWHM 1.704*newstep).

rectangle

This window simply zeros all spatial frequencies above the highest that can be correctly sampled by the new pixel size. This gives the best resolution of any of the windows, but this is marred by the strong sidelobes of the resulting airy-profile, especially near bright point sources and CCD saturation lines.

Returns:
outImage

The resampled image is returned.

resample(newdim, newstart, newstep, flux=False, order=1, interp='no', unit_start=Unit('deg'), unit_step=Unit('arcsec'), antialias=True, inplace=False, window='blackman')[source]

Resample an image of the sky to select its angular resolution and to specify which sky position appears at the center of pixel [0,0].

This function is a simplified interface to the mpdaf.obj.Image.regrid function, which it calls with the following arguments:

regrid(newdim, newstart, [0.0, 0.0],
       [abs(newstep[0]),-abs(newstep[1])]
       flux=flux, order=order, interp=interp, unit_pos=unit_start,
       unit_inc=unit_step, inplace=inplace)

When this function is used to resample an image to a lower resolution, a low-pass anti-aliasing filter is applied to the image before it is resampled, to remove all spatial frequencies below half the new sampling rate. This is required to satisfy the Nyquist sampling constraint. It prevents high spatial-frequency noise and edges from being folded into lower frequency artefacts in the resampled image. The removal of this noise improves the signal to noise ratio of the resampled image.

Parameters:
newdimint or (int,int)

The desired new dimensions. Python notation: (ny,nx)

newstartfloat or (float, float)

The sky position (dec,ra) that should appear at the center of pixel [0,0].

If None, the value of self.get_start() is substituted, so that the sky position that appears at the center of pixel [0,0] is unchanged by the resampling operation.

newstepfloat or (float, float)

The desired angular size of the image pixels on the sky. The size is expressed as either one number to request square pixels on the sky with that width and height, or two numbers that specify the height and width of rectangular pixels on the sky. In the latter case, the two numbers are the size along the Y axis of the image array followed by the size along the X axis.

fluxbool

This tells the function whether the pixel units of the image are flux densities (flux=True), such as erg/s/cm2/Hz, or whether they are per-steradian brightness units (flux=False), such as erg/s/cm2/Hz/steradian. It needs to know this when it changes the pixel size, because when pixel sizes change, resampled flux densities need to be corrected for the change in the area per pixel, where resampled brightnesses don’t.

orderint

The order of the spline interpolation. This can take any value from 0-5. The default is 1 (linear interpolation). When this function is used to lower the resolution of an image, the low-pass anti-aliasing filter that is applied, makes linear interpolation sufficient. Conversely, when this function is used to increase the image resolution, order=3 might be useful. Higher orders than this will tend to introduce ringing artefacts.

interp‘no’ | ‘linear’ | ‘spline’

If ‘no’, replace masked data with the median image value. If ‘linear’, replace masked values using a linear interpolation between neighboring values. if ‘spline’, replace masked values using a spline interpolation between neighboring values.

unit_startastropy.units.Unit

The angular units of the newstart coordinates. Degrees by default.

unit_stepastropy.units.Unit

The angular units of the step argument. Arcseconds by default.

antialiasbool

By default, when the resolution of an image axis is about to be reduced, a low pass filter is first applied to suppress high spatial frequencies that can not be represented by the reduced sampling interval. If this is not done, high-frequency noise and sharp edges get folded back to lower frequencies, where they increase the noise level of the image and introduce ringing artefacts next to sharp edges, such as CCD saturation spikes. This filtering can be disabled by passing False to the antialias argument.

inplacebool

If False, return a rotated copy of the image (the default). If True, rotate the original image in-place, and return that.

windowstr

The type of window function to use for antialiasing in the Fourier plane. The following windows are supported:

blackman

This window suppresses ringing better than any other window, at the expense of lowered image resolution. In the image plane, the PSF of this window is approximately gaussian, with a standard deviation of around 0.96*newstep, and a FWHM of about 2.3*newstep.

gaussian

A truncated gaussian window. This has a smaller PSF than the blackman window, however gaussians never fall to zero, so either significant ringing will be seen due to truncation of the gaussian, or low-level aliasing will occur, depending on the spatial frequency coverage of the image beyond the folding frequency. It can be a good choice for images that only contain smoothly varying features. It is equivalent to a convolution of the image with both an airy profile and a gaussian of standard deviation 0.724*newstep (FWHM 1.704*newstep).

rectangle

This window simply zeros all spatial frequencies above the highest that can be correctly sampled by the new pixel size. This gives the best resolution of any of the windows, but this is marred by the strong sidelobes of the resulting airy-profile, especially near bright point sources and CCD saturation lines.

Returns:
outImage

The resampled image.

rotate(theta=0.0, interp='no', reshape=False, order=1, pivot=None, unit=Unit('deg'), regrid=None, flux=False, cutoff=0.25, inplace=False)[source]

Rotate the sky within an image in the sense of a rotation from north to east.

For example if the image rotation angle that is currently returned by image.get_rot() is zero, image.rotate(10.0) will rotate the northward direction of the image 10 degrees eastward of where it was, and self.get_rot() will thereafter return 10.0.

Uses scipy.ndimage.affine_transform.

Parameters:
thetafloat

The angle to rotate the image (degrees). Positive angles rotate features in the image in the sense of a rotation from north to east.

interp‘no’ | ‘linear’ | ‘spline’

If ‘no’, replace masked data with the median value of the image. This is the default. If ‘linear’, replace masked values using a linear interpolation between neighboring values. if ‘spline’, replace masked values using a spline interpolation between neighboring values.

reshapebool

If True, the size of the output image array is adjusted so that the input image is contained completely in the output. The default is False.

orderint

The order of the prefilter that is applied by the affine transform function. Prefiltering is not really needed for band-limited images, but this option is retained for backwards compatibility with an older version of the image.rotate method. In general orders > 1 tend to generate ringing at sharp edges, such as those of CCD saturation spikes, so this argument is best left with its default value of 1.

pivotfloat,float or None

When the reshape option is True, or the pivot argument is None, the image is rotated around its center. Alternatively, when the reshape option is False, the pivot argument can be used to indicate which pixel index [y,x] the image will be rotated around. Integer pixel indexes specify the centers of pixels. Non-integer values can be used to indicate positions between pixel centers.

On the sky, the rotation always occurs around the coordinate reference position of the observation. However the rotated sky is then mapped onto the pixel array of the image in such a way as to keep the sky position of the pivot pixel at the same place. This makes the image appear to rotate around that pixel.

unitastropy.units.Unit

The angular units of the rotation angle, theta.

regridbool

When this option is True, the pixel sizes along each axis are adjusted to avoid undersampling or oversampling any direction in the original image that would otherwise be rotated onto a lower or higher resolution axis. This is particularly important for images whose pixels have different angular dimensions along the X and Y axes, but it can also be important for images with square pixels, because the diagonal of an image with square pixels has higher resolution than the axes of that image.

If this option is left with its default value of None, then it is given the value of the reshape option.

fluxbool

This tells the function whether the pixel units of the image are flux densities (flux=True), such as erg/s/cm2/Hz, or whether they are per-steradian brightness units (flux=False), such as erg/s/cm2/Hz/steradian. It needs to know this when it changes the pixel size, because when pixel sizes change, resampled flux densities need to be corrected for the change in the area per pixel, where resampled brightnesses don’t.

cutofffloat

Mask each output pixel where at least this fraction of the pixel was interpolated from dummy values given to masked input pixels.

inplacebool

If False, return a rotated copy of the image (the default). If True, rotate the original image in-place, and return that.

Returns:
outImage
segment(shape=(2, 2), minsize=20, minpts=None, background=20, interp='no', median=None)[source]

Segment the image in a number of smaller images.

Returns a list of images. Uses scipy.ndimage.generate_binary_structure, scipy.ndimage.grey_dilation, scipy.ndimage.measurements.label, and scipy.ndimage.measurements.find_objects.

Parameters:
shape(int,int)

Shape used for connectivity.

minsizeint

Minimmum size of the images.

minptsint

Minimmum number of points in the object.

backgroundfloat

Under this value, flux is considered as background.

interp‘no’ | ‘linear’ | ‘spline’

if ‘no’, data median value replaced masked values. if ‘linear’, linear interpolation of the masked values. if ‘spline’, spline interpolation of the masked values.

median(int,int) or None

If not None (default), size of the window to apply a median filter on the image.

Returns:
outlist of Image
set_spatial_fmax(newfmax=None, rot=None)[source]

Specify the spatial-frequency band-limits of the image along the Y and X axis. This function completely replaces any existing band-limits. See also update_spatial_fmax().

The recorded limits are used to avoid redundantly performing anti-aliasing measures such as low-pass filtering an image before resampling to a lower resolution, or decreasing pixel sizes before rotating high resolution axes onto low resolution axes.

Parameters:
newfmaxnumpy.ndarray

The new frequency limits along the Y and X axes or a band-limiting ellipse, specified in cycles per the angular unit in self.wcs.unit.

rotfloat or None

Either None, to specify band-limits that pertain to the Y and X axes of the current image without any rotation, or, if the band-limits pertain to a rotated version of the image, the rotation angle of its Y axis westward of north (degrees). This is defined such that if image.wcs.get_rot() is passed to this function, the band-limit newfmax[0] will be along the Y axis of the image and newfmax[1] will be along its X axis.

set_wcs(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.

sqrt(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.

subimage(center, size, unit_center=Unit('deg'), unit_size=Unit('arcsec'), minsize=2.0)[source]

Return a view on a square or rectangular part.

This method returns a square or rectangular sub-image whose center and size are specified in world coordinates. Note that this is a view on the original map and that both will be modified at the same time. If you need to modify only the sub-image, copy() the result of the method.

Parameters:
center(float,float)

The center (dec, ra) of the square region. If this position is not within the parent image, None is returned.

sizefloat or (float,float)

The width of a square region, or the width and height of a rectangular region.

unit_centerastropy.units.Unit

The units of the center coordinates. Degrees are assumed by default. To specify the center in pixels, assign None to unit_center.

unit_sizeastropy.units.Unit

The units of the size and minsize arguments. Arcseconds are assumed by default (use None to specify sizes in pixels).

minsizefloat

The minimum width of the output image along both the Y and X axes. This function returns None if size is smaller than minsize, or if the part of the square that lies within the parent image is smaller than minsize along either axis.

Returns:
outImage
to_ds9(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(y_min, y_max, x_min, x_max, mask=True, unit=Unit('deg'), inplace=False)[source]

Return a sub-image that contains a specified area of the sky.

The ranges x_min to x_max and y_min to y_max, specify a rectangular region of the sky in world coordinates. The truncate function returns the sub-image that just encloses this region. Note that if the world coordinate axes are not parallel to the array axes, the region will appear to be a rotated rectangle within the sub-image. In such cases, the corners of the sub-image will contain pixels that are outside the region. By default these pixels are masked. However this can be disabled by changing the optional mask argument to False.

Parameters:
y_minfloat

The minimum Y-axis world-coordinate of the selected region. The Y-axis is usually Declination, which may not be parallel to the Y-axis of the image array.

y_maxfloat

The maximum Y-axis world coordinate of the selected region.

x_minfloat

The minimum X-axis world-coordinate of the selected region. The X-axis is usually Right Ascension, which may not be parallel to the X-axis of the image array.

x_maxfloat

The maximum X-axis world coordinate of the selected region.

maskbool

If True, any pixels in the sub-image that remain outside the range x_min to x_max and y_min to y_max, will be masked.

unitastropy.units.Unit

The units of the X and Y world-coordinates (degrees by default).

inplacebool

If False, return a truncated copy of the image (the default). If True, truncate the original image in-place, and return that.

Returns:
outImage
unmask()

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

update_spatial_fmax(newfmax, rot=None)[source]

Update the spatial-frequency band-limits recorded for the current image.

See the documentation of set_spatial_fmax() for an explanation of what the band-limits are used for.

If either of the new limits is less than an existing band-limit, and the rotation angle of the new limits is the same as the angle of the recorded limits, then the smaller limits replace the originals.

If either of the new limits is smaller than the existing limits, but the rotation angle for the new limits differs from the recorded limits, then both of the original limits are discarded and replaced by the new ones at the specified angle.

Parameters:
newfmaxnumpy.ndarray

The frequency limits along the Y and X axes, respectively, specified in cycles per the angular unit in self.wcs.unit.

rotfloat or None

Either None, to specify band-limits that pertain to the Y and X axes of the current image without any rotation, or, if the band-limits pertain to a rotated version of the image, the rotation angle of its Y axis westward of north (degrees). This is defined such that if image.wcs.get_rot() is passed to this function, the band-limit newfmax[0] will be along the Y axis of the image and newfmax[1] will be along its X axis.

write(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.