Skip to content

Commit

Permalink
Fix axis 18 (#20)
Browse files Browse the repository at this point in the history
* Axis order changed

* Camera.py docstring improved

* clipping below zero pixel values

* Correcting test_psf.py

* camera.py docstring improved

* below zero assert added

* below zero and infinit asserts added

* combining pairs

* `@pytest.mark.parameterize` implemented

* Some problems with ...
  • Loading branch information
bazkiaei authored Sep 27, 2018
1 parent 7b6fcbf commit 81f4855
Show file tree
Hide file tree
Showing 4 changed files with 111 additions and 75 deletions.
3 changes: 2 additions & 1 deletion gunagala/camera.py
Original file line number Diff line number Diff line change
Expand Up @@ -38,7 +38,8 @@ class Camera:
Pixel pitch. Square pixels are assumed.
Resolution : astropy.units.Quantity
Two element Quantity containing the number of pixels across
the image sensor in both horizontal & vertical directions.
the image sensor in both vertical & horizontal directions.
(y, x)
read_noise astropy.units.Quantity
Intrinsic noise of image sensor and readout electronics, in
electrons/pixel units.
Expand Down
6 changes: 3 additions & 3 deletions gunagala/imager.py
Original file line number Diff line number Diff line change
Expand Up @@ -248,9 +248,9 @@ def __init__(self, optic, camera, filters, psf, sky, num_imagers=1, num_per_comp

# Construct a simple template WCS to store the focal plane configuration parameters
self.wcs = WCS(naxis=2)
self.wcs._naxis1, self.wcs._naxis2 = self.camera.resolution.value.astype(int)
self.wcs.wcs.crpix = [(self.camera.resolution[0].value - 1)/2,
(self.camera.resolution[1].value - 1)/2]
self.wcs._naxis2, self.wcs._naxis1 = self.camera.resolution.value.astype(int)
self.wcs.wcs.crpix = [(self.camera.resolution[1].value - 1)/2,
(self.camera.resolution[0].value - 1)/2]
self.wcs.wcs.cdelt = [self.pixel_scale.to(u.degree / u.pixel).value,
self.pixel_scale.to(u.degree / u.pixel).value]
self.wcs.wcs.ctype = ['RA---TAN', 'DEC--TAN']
Expand Down
22 changes: 13 additions & 9 deletions gunagala/psf.py
Original file line number Diff line number Diff line change
Expand Up @@ -12,7 +12,7 @@
from astropy.modeling.functional_models import Moffat2D

from gunagala import utils

from warnings import warn

class PSF():
"""
Expand Down Expand Up @@ -177,9 +177,9 @@ def pixellated(self, size=(21, 21), offsets=(0.0, 0.0)):
Parameters
----------
size : (int, int) optional
x, y size of the pixellated PSF to calculate. Default value (21, 21).
y, x size of the pixellated PSF to calculate. Default value (21, 21).
offset : tuple of floats, optional
x and y axis offsets of the centre of the PSF from the centre
y and x axis offsets of the centre of the PSF from the centre
of the returned image, in pixels.
Returns
Expand All @@ -195,11 +195,11 @@ def pixellated(self, size=(21, 21), offsets=(0.0, 0.0)):
raise ValueError("`size` must be > 0, got {}!".format(size))

# Update PSF centre coordinates
self.x_0 = offsets[0]
self.y_0 = offsets[1]
self.x_0 = offsets[1]
self.y_0 = offsets[0]

xrange = (-(size[0] - 1) / 2, (size[0] + 1) / 2)
yrange = (-(size[1] - 1) / 2, (size[1] + 1) / 2)
xrange = (-(size[1] - 1) / 2, (size[1] + 1) / 2)
yrange = (-(size[0] - 1) / 2, (size[0] + 1) / 2)

return discretize_model(self, xrange, yrange, mode='oversample', factor=10)

Expand Down Expand Up @@ -360,7 +360,6 @@ def pixellated(self, size=(21, 21), offsets=(0.0, 0.0)):
"""
size = np.array(size, dtype=np.int)
offsets = np.array(offsets)

# Only want to caclulate resampled PSF for positions that fall within the PSF data,
# otherwise end up filling the RAM with lots of double precision zeros.

Expand Down Expand Up @@ -413,9 +412,14 @@ def pixellated(self, size=(21, 21), offsets=(0.0, 0.0)):
resampled_psf = utils.bin_array(resampled_psf, self._oversampling)
# Renormalise to correct for the effect of resampling
resampled_psf = resampled_psf / self._resampling_factor**2
# Check and if there are some below zero pixel values, clip them to zero and renormalize resampled_psf.
if (resampled_psf < 0).any():
warn("Warning: below zero values in resampled PSF. Clipping to zero.")
previous_total = resampled_psf.sum()
resampled_psf[resampled_psf < 0] = 0
resampled_psf *= previous_total / resampled_psf.sum()
# Insert into output array in the correct place.
pixellated[y0:y1,x0:x1] = resampled_psf

return pixellated

def _get_n_pix(self):
Expand Down
155 changes: 93 additions & 62 deletions gunagala/tests/test_psf.py
Original file line number Diff line number Diff line change
Expand Up @@ -6,13 +6,13 @@


@pytest.fixture(scope='module')
def psf():
def psf_moffat():
psf = MoffatPSF(FWHM=1 / 30 * u.arcminute, shape=4.7)
return psf


@pytest.fixture(scope='module')
def pix_psf():
def psf_pixellated():
psf_data = np.array([[0.0, 0.0, 0.1, 0.0, 0.0],
[0.0, 0.3, 0.7, 0.4, 0.0],
[0.1, 0.8, 1.0, 0.6, 0.1],
Expand All @@ -32,84 +32,115 @@ def test_base():
psf_base = PSF(FWHM=1 / 30 * u.arcminute)


def test_moffat(psf):
assert isinstance(psf, MoffatPSF)
assert isinstance(psf, PSF)
def test_moffat(psf_moffat):
assert isinstance(psf_moffat, MoffatPSF)
assert isinstance(psf_moffat, PSF)


def test_pix(pix_psf):
assert isinstance(pix_psf, PixellatedPSF)
assert isinstance(pix_psf, PSF)
def test_pix(psf_pixellated):
assert isinstance(psf_pixellated, PixellatedPSF)
assert isinstance(psf_pixellated, PSF)


def test_FWHM(psf):
assert psf.FWHM == 2 * u.arcsecond
psf.FWHM = 4 * u.arcsecond
assert psf.FWHM == 1 / 15 * u.arcminute
def test_FWHM(psf_moffat):
assert psf_moffat.FWHM == 2 * u.arcsecond
psf_moffat.FWHM = 4 * u.arcsecond
assert psf_moffat.FWHM == 1 / 15 * u.arcminute
with pytest.raises(ValueError):
psf.FWHM = -1 * u.degree
psf.FWHM = 2 * u.arcsecond
psf_moffat.FWHM = -1 * u.degree
psf_moffat.FWHM = 2 * u.arcsecond


def test_pixel_scale(psf):
psf.pixel_scale = 2.85 * u.arcsecond / u.pixel
assert psf.pixel_scale == 2.85 * u.arcsecond / u.pixel
def test_pixel_scale(psf_moffat):
psf_moffat.pixel_scale = 2.85 * u.arcsecond / u.pixel
assert psf_moffat.pixel_scale == 2.85 * u.arcsecond / u.pixel


def test_pixel_scale_pix(pix_psf):
pix_psf.pixel_scale = (1 / 3) * u.arcsecond / u.pixel
assert pix_psf.pixel_scale == (1 / 3) * u.arcsecond / u.pixel
pix_psf.pixel_scale = (2 / 3) * u.arcsecond / u.pixel
def test_pixel_scale_pix(psf_pixellated):
psf_pixellated.pixel_scale = (1 / 3) * u.arcsecond / u.pixel
assert psf_pixellated.pixel_scale == (1 / 3) * u.arcsecond / u.pixel
psf_pixellated.pixel_scale = (2 / 3) * u.arcsecond / u.pixel


def test_n_pix(psf):
assert psf.n_pix == 4.25754067000986 * u.pixel
moffat = psf_moffat()
pixellated = psf_pixellated()


def test_n_pix_pix(pix_psf):
assert pix_psf.n_pix / u.pixel == pytest.approx(21.01351017)
# @pytest.mark.parametrize("psf, expected_n_pix", [
# (moffat, 4.25754067000986),
# (pixellated, pytest.approx(21.06994544))],
# ids=["moffat", "pixellated"]
# )
# def test_n_pix(psf, expected_n_pix):
# assert psf.n_pix == expected_n_pix


def test_peak(psf):
assert psf.peak == 0.7134084656751443 / u.pixel
def test_n_pix(psf_moffat):
assert psf_moffat.n_pix == 4.25754067000986 * u.pixel


def test_peak_pix(pix_psf):
assert pix_psf.peak * u.pixel == pytest.approx(0.08073066)
def test_n_pix_pix(psf_pixellated):
assert psf_pixellated.n_pix / u.pixel == pytest.approx(21.069945447)


def test_shape(psf):
assert psf.shape == 4.7
psf.shape = 2.5
assert psf.shape == 2.5
with pytest.raises(ValueError):
psf.shape = 0.5
psf.shape = 4.7


def test_pixellated(psf):
pixellated = psf.pixellated()
assert isinstance(pixellated, np.ndarray)
assert pixellated.shape == (21, 21)
pixellated = psf.pixellated(size=(7.2, 7.2))
assert isinstance(pixellated, np.ndarray)
assert pixellated.shape == (7, 7)
pixellated = psf.pixellated(offsets=(0.3, -0.7))
assert isinstance(pixellated, np.ndarray)
assert pixellated.shape == (21, 21)
def test_peak(psf_moffat):
assert psf_moffat.peak == 0.7134084656751443 / u.pixel


def test_peak_pix(psf_pixellated):
assert psf_pixellated.peak * u.pixel == pytest.approx(0.08073066)


def test_shape(psf_moffat):
assert psf_moffat.shape == 4.7
psf_moffat.shape = 2.5
assert psf_moffat.shape == 2.5
with pytest.raises(ValueError):
psf.pixellated(size=(1.3, -1.3))


def test_pixellated(pix_psf):
pixellated = pix_psf.pixellated()
assert isinstance(pixellated, np.ndarray)
assert pixellated.shape == (21, 21)
pixellated = pix_psf.pixellated(size=(7.2, 7.2))
assert isinstance(pixellated, np.ndarray)
assert pixellated.shape == (7, 7)
pixellated = pix_psf.pixellated(offsets=(0.3, -0.7))
assert isinstance(pixellated, np.ndarray)
assert pixellated.shape == (21, 21)
psf_moffat.shape = 0.5
psf_moffat.shape = 4.7


@pytest.mark.parametrize("psf, image_size", [
(moffat, (21, 21)),
(pixellated, (21, 21))],
ids=["moffat", "pixellated"]
)
def test_pixellated_square(psf, image_size):
assert isinstance(psf.pixellated(), np.ndarray)
assert psf.pixellated().shape == image_size
assert (psf.pixellated() >= 0).all()
assert np.isfinite(psf.pixellated()).all()


@pytest.mark.parametrize("psf, image_size", [
(moffat, (7, 9)),
(pixellated, (7, 9))],
ids=["moffat", "pixellated"]
)
def test_pixellated_rectangle(psf, image_size):
assert isinstance(psf.pixellated(size=(7.2, 9.2)), np.ndarray)
assert psf.pixellated(size=(7.2, 9.2)).shape == image_size
assert (psf.pixellated(size=(7.2, 9.2)) >= 0).all()
assert np.isfinite(psf.pixellated(size=(7.2, 9.2))).all()


@pytest.mark.parametrize("psf, image_size", [
(moffat, (21, 21)),
(pixellated, (21, 21))],
ids=["moffat", "pixellated"]
)
def test_pixellated_offsets(psf, image_size):
assert isinstance(psf.pixellated(), np.ndarray)
assert psf.pixellated().shape == image_size
assert (psf.pixellated() >= 0).all()
assert np.isfinite(psf.pixellated()).all()


@pytest.mark.parametrize("psf, test_size", [
(moffat, (1.3, -1.3)),
(pixellated, (-1.3, 1.3))],
ids=["moffat", "pixellated"]
)
def test_pixellated(psf, test_size):
with pytest.raises(ValueError):
pix_psf.pixellated(size=(-1.3, 1.3))
psf.pixellated(size=test_size)

0 comments on commit 81f4855

Please sign in to comment.