Skip to content

Commit

Permalink
Merge pull request #1741 from larrybradley/psfphot-model
Browse files Browse the repository at this point in the history
Check dimensionality of PSF model
  • Loading branch information
larrybradley committed May 17, 2024
2 parents 473c6e0 + 49837f2 commit 655481f
Show file tree
Hide file tree
Showing 3 changed files with 71 additions and 17 deletions.
6 changes: 6 additions & 0 deletions CHANGES.rst
Original file line number Diff line number Diff line change
Expand Up @@ -37,6 +37,12 @@ API Changes
functions from ``photutils.datasets.make``, import functions from
``photutils.datasets``. [#1726]

- ``photutils.psf``

- ``PSFPhotometry`` and ``IterativePSFPhotometry`` now raise a
``ValueError`` if the input ``psf_model`` is not two-dimensional
with ``n_inputs=2`` and ``n_outputs=1``. [#1741]


1.12.0 (2024-04-12)
-------------------
Expand Down
47 changes: 31 additions & 16 deletions photutils/psf/photometry.py
Original file line number Diff line number Diff line change
Expand Up @@ -43,12 +43,14 @@ class PSFPhotometry:
Parameters
----------
psf_model : `astropy.modeling.Fittable2DModel`
psf_model : 2D `astropy.modeling.Model`
The PSF model to fit to the data. The model must have parameters
named ``x_0``, ``y_0``, and ``flux``, corresponding to the
center (x, y) position and flux, or it must have 'x_name',
'y_name', and 'flux_name' attributes that map to the x, y, and
flux parameters (i.e., a model output from `prepare_psf_model`).
The model must be two-dimensional such that it accepts 2 inputs
(e.g., x and y) and provides 1 output.
fit_shape : int or length-2 array_like
The rectangular shape around the center of a star that will
Expand Down Expand Up @@ -219,13 +221,22 @@ def _validate_psf_model(self):
"""
Validate the input PSF model.
The PSF model must be a subclass of `astropy.modeling.Model`. It
must also be two-dimensional and have a single output.
The PSF model must have parameters called 'x_0', 'y_0', and
'flux' or it must have 'x_name', 'y_name', and 'flux_name'
attributes (i.e., output from `prepare_psf_model`). Otherwise, a
`ValueError` is raised.
"""
if not isinstance(self.psf_model, Model):
raise TypeError('psf_model must be an Astropy Model subclass.')

if self.psf_model.n_inputs != 2 or self.psf_model.n_outputs != 1:
raise ValueError('psf_model must be two-dimensional with '
'n_inputs=2 and n_outputs=1.')

# check for required PSF model parameters
_ = self._psf_param_names

@staticmethod
Expand Down Expand Up @@ -641,7 +652,7 @@ def _make_fit_results(self, models, infos):
fit_infos.extend([info] * npsf_models) # views
fit_param_errs.extend(self._split_param_errs(param_err, nfitparam))

if len(fit_models) != len(fit_infos):
if len(fit_models) != len(fit_infos): # pragma: no cover
raise ValueError('fit_models and fit_infos have different lengths')

# change the sorting from group_id to source id order
Expand All @@ -666,9 +677,9 @@ def _fit_sources(self, data, init_params, *, error=None, mask=None):
ungroup_idx = np.argsort(sources['id'].value)
self._group_results['ungroup_indices'] = ungroup_idx
sources = sources.groups
if self.progress_bar:
if self.progress_bar: # pragma: no cover
desc = 'Fit source/group'
sources = add_progress_bar(sources, desc=desc) # pragma: no cover
sources = add_progress_bar(sources, desc=desc)

# Save the fit_info results for these keys if they are present.
# Some of these keys are returned only by some fitters. These
Expand Down Expand Up @@ -795,10 +806,10 @@ def _calc_fit_metrics(self, data, source_tbl):

for npixfit, residuals in zip(self.fit_results['npixfit'],
fit_residuals):
if len(residuals) != npixfit:
if len(residuals) != npixfit: # pragma: no cover
raise ValueError('size of residuals does not match npixfit')

if len(fit_residuals) != len(source_tbl):
if len(fit_residuals) != len(source_tbl): # pragma: no cover
raise ValueError('fit_residuals does not match the source '
'table length')

Expand Down Expand Up @@ -1027,14 +1038,14 @@ def __call__(self, data, *, mask=None, error=None, init_params=None):

# create output table
fit_sources = self._model_params_to_table(fit_models) # ungrouped
if len(init_params) != len(fit_sources):
if len(init_params) != len(fit_sources): # pragma: no cover
raise ValueError('init_params and fit_sources tables have '
'different lengths')
source_tbl = hstack((init_params, fit_sources))

param_errors = self._param_errors_to_table()
if len(param_errors) > 0:
if len(param_errors) != len(source_tbl):
if len(param_errors) != len(source_tbl): # pragma: no cover
raise ValueError('param errors and fit sources tables have '
'different lengths')
source_tbl = hstack((source_tbl, param_errors))
Expand Down Expand Up @@ -1104,9 +1115,9 @@ def make_model_image(self, shape, psf_shape, *, include_localbkg=False):
data = np.zeros(shape)
xname, yname = self._psf_param_names[0:2]

if self.progress_bar:
if self.progress_bar: # pragma: no cover
desc = 'Model image'
fit_models = add_progress_bar(fit_models, desc=desc) # pragma: no cover
fit_models = add_progress_bar(fit_models, desc=desc)

# fit_models must be a list of individual, not grouped, PSF
# models, i.e., there should be one PSF model (which may be
Expand Down Expand Up @@ -1189,10 +1200,14 @@ class IterativePSFPhotometry:
Parameters
----------
psf_model : `astropy.modeling.Fittable2DModel`
The PSF model to fit to the data. The model needs to have
three parameters named ``x_0``, ``y_0``, and ``flux``,
corresponding to the center (x, y) position and flux.
psf_model : 2D `astropy.modeling.Model`
The PSF model to fit to the data. The model must have parameters
named ``x_0``, ``y_0``, and ``flux``, corresponding to the
center (x, y) position and flux, or it must have 'x_name',
'y_name', and 'flux_name' attributes that map to the x, y, and
flux parameters (i.e., a model output from `prepare_psf_model`).
The model must be two-dimensional such that it accepts 2 inputs
(e.g., x and y) and provides 1 output.
fit_shape : int or length-2 array_like
The rectangular shape around the center of a star that will
Expand Down Expand Up @@ -1563,9 +1578,9 @@ def make_model_image(self, shape, psf_shape):
data = np.zeros(shape)
xname, yname = self.fit_results[0]._psf_param_names[0:2]

if self.psfphot.progress_bar:
if self.psfphot.progress_bar: # pragma: no cover
desc = 'Model image'
fit_models = add_progress_bar(fit_models, desc=desc) # pragma: no cover
fit_models = add_progress_bar(fit_models, desc=desc)

# fit_models must be a list of individual, not grouped, PSF
# models, i.e., there should be one PSF model (which may be
Expand Down
35 changes: 34 additions & 1 deletion photutils/psf/tests/test_photometry.py
Original file line number Diff line number Diff line change
Expand Up @@ -7,7 +7,7 @@
import numpy as np
import pytest
from astropy.modeling.fitting import LMLSQFitter, SimplexLSQFitter
from astropy.modeling.models import Gaussian2D
from astropy.modeling.models import Gaussian1D, Gaussian2D, custom_model
from astropy.nddata import NDData, StdDevUncertainty, overlap_slices
from astropy.table import QTable, Table
from astropy.utils.exceptions import (AstropyDeprecationWarning,
Expand All @@ -33,6 +33,27 @@ def test_inputs():
with pytest.raises(TypeError, match=match):
_ = PSFPhotometry(1, 3)

match = 'psf_model must be two-dimensional'
with pytest.raises(ValueError, match=match):
psf_model = Gaussian1D()
_ = PSFPhotometry(psf_model, 3)

@custom_model
def my_model(x, y, flux=1, x_0=0, y_0=0, sigma=1):
return flux, flux * 2
m = my_model()
m.n_outputs = 2
psf_model = my_model()
match = 'psf_model must be two-dimensional'
with pytest.raises(ValueError, match=match):
psf_model = Gaussian1D()
_ = PSFPhotometry(psf_model, 3)

match = 'Invalid PSF model - could not find PSF parameter names'
with pytest.raises(ValueError, match=match):
psf_model = Gaussian2D()
_ = PSFPhotometry(psf_model, 3)

match = 'fit_shape must have an odd value for both axes'
for shape in ((0, 0), (4, 3)):
with pytest.raises(ValueError, match=match):
Expand Down Expand Up @@ -264,6 +285,18 @@ def test_psf_photometry_compound(test_data):
assert isinstance(phot, QTable)
assert len(phot) == len(sources)

psf_model2 = psf_model.copy()
psf_model2.xname = psf_model2.x_name
psf_model2.yname = psf_model2.y_name
psf_model2.fluxname = psf_model2.flux_name
del psf_model2.x_name, psf_model2.y_name, psf_model2.flux_name

psfphot2 = PSFPhotometry(psf_model2, fit_shape, finder=finder,
aperture_radius=4)
phot2 = psfphot2(data, error=error)
assert isinstance(phot2, QTable)
assert len(phot2) == len(sources)

# test results when fit does not converge (fitter_maxiters=10)
match = r'One or more fit\(s\) may not have converged.'
with pytest.warns(AstropyUserWarning, match=match):
Expand Down

0 comments on commit 655481f

Please sign in to comment.