CCDData#

class astropy.nddata.CCDData(*args, **kwd)[source]#

Bases: NDDataArray

A class describing basic CCD data.

The CCDData class is based on the NDData object and includes a data array, uncertainty frame, mask frame, flag frame, meta data, units, and WCS information for a single CCD image.

Parameters:
dataCCDDataastropy:-like or array_like

The actual data contained in this CCDData object. Note that the data will always be saved by reference, so you should make a copy of the data before passing it in if that’s the desired behavior.

uncertaintyStdDevUncertainty, VarianceUncertainty, InverseVariance, numpy.ndarray or None, optional

Uncertainties on the data. If the uncertainty is a numpy.ndarray, it it assumed to be, and stored as, a StdDevUncertainty. Default is None.

masknumpy.ndarray or None, optional

Mask for the data, given as a boolean Numpy array with a shape matching that of the data. The values must be False where the data is valid and True when it is not (like Numpy masked arrays). If data is a numpy masked array, providing mask here will causes the mask from the masked array to be ignored. Default is None.

flagsnumpy.ndarray or FlagCollection or None, optional

Flags giving information about each pixel. These can be specified either as a Numpy array of any type with a shape matching that of the data, or as a FlagCollection instance which has a shape matching that of the data. Default is None.

wcsWCS or None, optional

WCS-object containing the world coordinate system for the data. Default is None.

metadict-like object or None, optional

Metadata for this object. “Metadata” here means all information that is included with this object but not part of any other attribute of this particular object, e.g. creation date, unique identifier, simulation parameters, exposure time, telescope name, etc.

unitUnit or str, optional

The units of the data. Default is None.

Warning

If the unit is None or not otherwise specified it will raise a ValueError

psfnumpy.ndarray or None, optional

Image representation of the PSF at the center of this image. In order for convolution to be flux-preserving, this should generally be normalized to sum to unity.

Raises:
ValueError

If the uncertainty or mask inputs cannot be broadcast (e.g., match shape) onto data.

Notes

CCDData objects can be easily converted to a regular

Numpy array using numpy.asarray.

For example:

>>> from astropy.nddata import CCDData
>>> import numpy as np
>>> x = CCDData([1,2,3], unit='adu')
>>> np.asarray(x)
array([1, 2, 3])

This is useful, for example, when plotting a 2D image using matplotlib.

>>> from astropy.nddata import CCDData
>>> from matplotlib import pyplot as plt   
>>> x = CCDData([[1,2,3], [4,5,6]], unit='adu')
>>> plt.imshow(x)   
Attributes:
known_invalid_fits_unit_strings

A dictionary that maps commonly-used fits unit name strings that are technically invalid to the correct valid unit type (or unit string). This is primarily for variant names like “ELECTRONS/S” which are not formally valid, but are unambiguous and frequently enough encountered that it is convenient to map them to the correct unit.

Methods

read(*args, **kwargs)

Classmethod to create an CCDData instance based on a FITS file. This method uses fits_ccddata_reader() with the provided parameters.

write(*args, **kwargs)

Writes the contents of the CCDData instance into a new FITS file. This method uses fits_ccddata_writer() with the provided parameters.

Attributes Summary

data

ndarray-like : The stored dataset.

header

known_invalid_fits_unit_strings

psf

Image representation of the PSF for the dataset.

uncertainty

any type : Uncertainty in the dataset, if any.

unit

Unit : Unit for the dataset, if any.

wcs

any type : A world coordinate system (WCS) for the dataset, if any.

Methods Summary

add(operand[, operand2])

See astropy.nddata.NDArithmeticMixin.add.

copy()

Return a copy of the CCDData object.

divide(operand[, operand2])

See astropy.nddata.NDArithmeticMixin.divide.

multiply(operand[, operand2])

See astropy.nddata.NDArithmeticMixin.multiply.

subtract(operand[, operand2])

See astropy.nddata.NDArithmeticMixin.subtract.

to_hdu([hdu_mask, hdu_uncertainty, ...])

Creates an HDUList object from a CCDData object.

Attributes Documentation

data#
header#
known_invalid_fits_unit_strings = {'ELECTRONS': Unit("electron"), 'ELECTRONS/S': Unit("electron / s"), 'electrons': Unit("electron")}#
psf#
uncertainty#
unit#
wcs#

Methods Documentation

classmethod add(operand, operand2=None, **kwargs)#

See astropy.nddata.NDArithmeticMixin.add.

copy()[source]#

Return a copy of the CCDData object.

classmethod divide(operand, operand2=None, **kwargs)#

See astropy.nddata.NDArithmeticMixin.divide.

classmethod multiply(operand, operand2=None, **kwargs)#

See astropy.nddata.NDArithmeticMixin.multiply.

classmethod subtract(operand, operand2=None, **kwargs)#

See astropy.nddata.NDArithmeticMixin.subtract.

to_hdu(hdu_mask='MASK', hdu_uncertainty='UNCERT', hdu_flags=None, wcs_relax=True, key_uncertainty_type='UTYPE', as_image_hdu=False, hdu_psf='PSFIMAGE')[source]#

Creates an HDUList object from a CCDData object.

Parameters:
hdu_mask, hdu_uncertainty, hdu_flags, hdu_psfstr or None, optional

If it is a string append this attribute to the HDUList as ImageHDU with the string as extension name. Flags are not supported at this time. If None this attribute is not appended. Default is 'MASK' for mask, 'UNCERT' for uncertainty, 'PSFIMAGE' for psf, and None for flags.

wcs_relaxbool

Value of the relax parameter to use in converting the WCS to a FITS header using to_header. The common CTYPE RA---TAN-SIP and DEC--TAN-SIP requires relax=True for the -SIP part of the CTYPE to be preserved.

key_uncertainty_typestr, optional

The header key name for the class name of the uncertainty (if any) that is used to store the uncertainty type in the uncertainty hdu. Default is UTYPE.

New in version 3.1.

as_image_hdubool

If this option is True, the first item of the returned HDUList is a ImageHDU, instead of the default PrimaryHDU.

Returns:
hdulistHDUList
Raises:
ValueError
  • If self.mask is set but not a numpy.ndarray.

  • If self.uncertainty is set but not a astropy uncertainty type.

  • If self.uncertainty is set but has another unit then self.data.

NotImplementedError

Saving flags is not supported.