Skip to content

Commit

Permalink
Merge d004813 into 51f02cd
Browse files Browse the repository at this point in the history
  • Loading branch information
mwcraig committed Jun 9, 2014
2 parents 51f02cd + d004813 commit 7944396
Show file tree
Hide file tree
Showing 3 changed files with 67 additions and 12 deletions.
35 changes: 28 additions & 7 deletions ccdproc/ccddata.py
Original file line number Diff line number Diff line change
Expand Up @@ -11,7 +11,7 @@
from astropy.io import fits, registry
from astropy.utils.compat.odict import OrderedDict
from astropy import units as u
import astropy
from astropy import log

from .utils.collections import CaseInsensitiveOrderedDict

Expand Down Expand Up @@ -39,7 +39,7 @@ class CCDData(NDData):
make copy the `data` before passing it in if that's the desired
behavior.
uncertainty : `~astropy.nddata.StdDevUncertainty`, optional
uncertainty : `~astropy.nddata.StdDevUncertainty` or `~numpy.ndarray`, optional
Uncertainties on the data.
mask : `~numpy.ndarray`, optional
Expand Down Expand Up @@ -143,10 +143,17 @@ def uncertainty(self, value):
if value is not None:
if isinstance(value, NDUncertainty):
self._uncertainty = value
self._uncertainty._parent_nddata = self
elif isinstance(value, np.ndarray):
if value.shape != self.shape:
raise ValueError("Uncertainty must have same shape as "
"data")
self._uncertainty = StdDevUncertainty(value)
log.info("Array provided for uncertainty; assuming it is a "
"StdDevUncertainty.")
else:
raise TypeError("Uncertainty must be an instance of a "
"NDUncertainty object")
"NDUncertainty object or a numpy array.")
self._uncertainty._parent_nddata = self
else:
self._uncertainty = value

Expand Down Expand Up @@ -262,8 +269,10 @@ def fits_ccddata_reader(filename, hdu=0, unit=None, **kwd):
hdu : int, optional
FITS extension from which CCDData should be initialized.
unit : astropy.units.Unit
Units of the image data
unit : astropy.units.Unit, optional
Units of the image data. If this argument is provided and there is a
unit for the image in the FITS header (the keyword ``BUNIT`` is used
as the unit, if present), this argument is used for the unit.
kwd :
Any additional keyword parameters are passed through to the FITS reader
Expand All @@ -286,7 +295,19 @@ def fits_ccddata_reader(filename, hdu=0, unit=None, **kwd):
prefix = 'Unsupported keyword: {0}.'.format(key)
raise TypeError(' '.join([prefix, msg]))
hdus = fits.open(filename, **kwd)
ccd_data = CCDData(hdus[hdu].data, meta=hdus[hdu].header, unit=unit)
hdr = hdus[hdu].header

try:
fits_unit_string = hdr['bunit']
except KeyError:
fits_unit_string = None

if unit is not None and fits_unit_string:
log.info("Using the unit {0} passed to the FITS reader instead of "
"the unit {1} in the FITS file.", unit, fits_unit_string)

use_unit = unit or fits_unit_string
ccd_data = CCDData(hdus[hdu].data, meta=hdus[hdu].header, unit=use_unit)
hdus.close()
return ccd_data

Expand Down
27 changes: 27 additions & 0 deletions ccdproc/tests/test_ccddata.py
Original file line number Diff line number Diff line change
Expand Up @@ -48,6 +48,21 @@ def test_initialize_from_FITS(ccd_data, tmpdir):
assert cd.meta[k] == v


def test_initialize_from_fits_with_unit_in_header(tmpdir):
fake_img = np.random.random(size=(100, 100))
hdu = fits.PrimaryHDU(fake_img)
hdu.header['bunit'] = u.adu.to_string()
filename = tmpdir.join('afile.fits').strpath
hdu.writeto(filename)
ccd = CCDData.read(filename)
# ccd should pick up the unit adu from the fits header...did it?
assert ccd.unit is u.adu

# An explicit unit in the read overrides any unit in the FITS file
ccd2 = CCDData.read(filename, unit="photon")
assert ccd2.unit is u.photon


def test_initialize_from_FITS_bad_keyword_raises_error(ccd_data, tmpdir):
# There are two fits.open keywords that are not permitted in ccdpro:
# do_not_scale_image_data and scale_back
Expand Down Expand Up @@ -135,6 +150,18 @@ def test_setting_bad_uncertainty_raises_error(ccd_data):
ccd_data.uncertainty = 10


def test_setting_uncertainty_with_array(ccd_data):
ccd_data.uncertainty = None
fake_uncertainty = np.sqrt(np.abs(ccd_data.data))
ccd_data.uncertainty = fake_uncertainty.copy()
np.testing.assert_array_equal(ccd_data.uncertainty.array, fake_uncertainty)


def test_setting_uncertainty_wrong_shape_raises_error(ccd_data):
with pytest.raises(ValueError):
ccd_data.uncertainty = np.random.random(size=2 * ccd_data.shape)


def test_to_hdu(ccd_data):
ccd_data.meta = {'observer': 'Edwin Hubble'}
fits_hdulist = ccd_data.to_hdu()
Expand Down
17 changes: 12 additions & 5 deletions docs/ccdproc/ccddata.rst
Original file line number Diff line number Diff line change
Expand Up @@ -45,8 +45,9 @@ A `~ccdproc.ccddata.CCDData` object can also be initialized from a FITS file:

>>> ccd = ccdproc.CCDData.read('my_file.fits', unit="adu") # doctest: +SKIP

but for the moment you need to set the unit explicitly, even if it is in the
FITS header.
If there is a unit in the FITS file (in the ``BUNIT`` keyword), that will be
used, but a unit explicitly provided in ``read`` will override any unit in the
FITS file.

There is no restriction at all on what the unit can be -- any unit in
`astropy.units` or that you create yourself will work.
Expand Down Expand Up @@ -131,16 +132,22 @@ Pixel-by-pixel uncertainty can be calculated for you:

See :ref:`create_variance` for more details.

You can also set the uncertainty directly but need to create a
You can also set the uncertainty directly, either by creating a
`~astropy.nddata.StdDevUncertainty` object first:

>>> from astropy.nddata.nduncertainty import StdDevUncertainty
>>> uncertainty = 0.1 * ccd.data # can be any array whose shape matches the data
>>> my_uncertainty = StdDevUncertainty(uncertainty)
>>> ccd.uncertainty = my_uncertainty

Using `~astropy.nddata.StdDevUncertainty` is required to enable error
propagation in `~ccdproc.ccddata.CCDData`
or by providing a `~numpy.ndarray` with the same shape as the data:

>>> ccd.uncertainty = 0.1 * ccd.data
INFO: Array provided for uncertainty; assuming it is a StdDevUncertainty. [ccdproc.ccddata]

In this case the uncertainty is assumed to be
`~astropy.nddata.StdDevUncertainty`. Using `~astropy.nddata.StdDevUncertainty`
is required to enable error propagation in `~ccdproc.ccddata.CCDData`

If you want access to the underlying uncertainty use its ``.array`` attribute:

Expand Down

0 comments on commit 7944396

Please sign in to comment.