BasicPSFPhotometry

class photutils.BasicPSFPhotometry(group_maker, bkg_estimator, psf_model, fitshape, finder=None, fitter=<astropy.modeling.fitting.LevMarLSQFitter object>, aperture_radius=None)[source]

Bases: object

This class implements a PSF photometry algorithm that can find sources in an image, group overlapping sources into a single model, fit the model to the sources, and subtracting the models from the image. This is roughly equivalent to the DAOPHOT routines FIND, GROUP, NSTAR, and SUBTRACT. This implementation allows a flexible and customizable interface to perform photometry. For instance, one is able to use different implementations for grouping and finding sources by using group_maker and finder respectivelly. In addition, sky background estimation is performed by bkg_estimator.

Parameters:
group_maker : callable or GroupStarsBase

group_maker should be able to decide whether a given star overlaps with any other and label them as beloging to the same group. group_maker receives as input an Table object with columns named as id, x_0, y_0, in which x_0 and y_0 have the same meaning of xcentroid and ycentroid. This callable must return an Table with columns id, x_0, y_0, and group_id. The column group_id should cotain integers starting from 1 that indicate which group a given source belongs to. See, e.g., DAOGroup.

bkg_estimator : callable, instance of any BackgroundBase subclass, or None

bkg_estimator should be able to compute either a scalar background or a 2D background of a given 2D image. See, e.g., MedianBackground. If None, no background subtraction is performed.

psf_model : astropy.modeling.Fittable2DModel instance

PSF or PRF model to fit the data. Could be one of the models in this package like DiscretePRF, IntegratedGaussianPRF, or any other suitable 2D model. This object needs to identify three parameters (position of center in x and y coordinates and the flux) in order to set them to suitable starting values for each fit. The names of these parameters should be given as x_0, y_0 and flux. prepare_psf_model can be used to prepare any 2D model to match this assumption.

fitshape : int or length-2 array-like

Rectangular shape around the center of a star which will be used to collect the data to do the fitting. Can be an integer to be the same along both axes. E.g., 5 is the same as (5, 5), which means to fit only at the following relative pixel positions: [-2, -1, 0, 1, 2]. Each element of fitshape must be an odd number.

finder : callable or instance of any StarFinderBase subclasses or None

finder should be able to identify stars, i.e. compute a rough estimate of the centroids, in a given 2D image. finder receives as input a 2D image and returns an Table object which contains columns with names: id, xcentroid, ycentroid, and flux. In which id is an integer-valued column starting from 1, xcentroid and ycentroid are center position estimates of the sources and flux contains flux estimates of the sources. See, e.g., DAOStarFinder. If finder is None, initial guesses for positions of objects must be provided.

fitter : Fitter instance

Fitter object used to compute the optimized centroid positions and/or flux of the identified sources. See fitting for more details on fitters.

aperture_radius : float or None

The radius (in units of pixels) used to compute initial estimates for the fluxes of sources. If None, one FWHM will be used if it can be determined from the psf_model.

Notes

Note that an ambiguity arises whenever finder and init_guesses (keyword argument for do_photometry) are both not None. In this case, finder is ignored and initial guesses are taken from init_guesses. In addition, an warning is raised to remaind the user about this behavior.

If there are problems with fitting large groups, change the parameters of the grouping algorithm to reduce the number of sources in each group or input a star_groups table that only includes the groups that are relevant (e.g. manually remove all entries that coincide with artifacts).

References

[1] Stetson, Astronomical Society of the Pacific, Publications,
(ISSN 0004-6280), vol. 99, March 1987, p. 191-222. Available at: http://adsabs.harvard.edu/abs/1987PASP…99..191S

Attributes Summary

aperture_radius
fitshape

Methods Summary

__call__(image[, init_guesses]) Performs PSF photometry.
do_photometry(image[, init_guesses]) Perform PSF photometry in image.
get_residual_image() Returns an image that is the result of the subtraction between the original image and the fitted sources.
nstar(image, star_groups) Fit, as appropriate, a compound or single model to the given star_groups.

Attributes Documentation

aperture_radius
fitshape

Methods Documentation

__call__(image, init_guesses=None)[source]

Performs PSF photometry. See do_photometry for more details including the __call__ signature.

do_photometry(image, init_guesses=None)[source]

Perform PSF photometry in image.

This method assumes that psf_model has centroids and flux parameters which will be fitted to the data provided in image. A compound model, in fact a sum of psf_model, will be fitted to groups of stars automatically identified by group_maker. Also, image is not assumed to be background subtracted. If init_guesses are not None then this method uses init_guesses as initial guesses for the centroids. If the centroid positions are set as fixed in the PSF model psf_model, then the optimizer will only consider the flux as a variable.

Parameters:
image : 2D array-like, ImageHDU, HDUList

Image to perform photometry.

init_guesses: `~astropy.table.Table`

Table which contains the initial guesses (estimates) for the set of parameters. Columns ‘x_0’ and ‘y_0’ which represent the positions (in pixel coordinates) for each object must be present. ‘flux_0’ can also be provided to set initial fluxes. If ‘flux_0’ is not provided, aperture photometry is used to estimate initial values for the fluxes. Additional columns of the form ‘<parametername>_0’ will be used to set the initial guess for any parameters of the psf_model model that are not fixed.

Returns:
output_tab : Table or None

Table with the photometry results, i.e., centroids and fluxes estimations and the initial estimates used to start the fitting process. Uncertainties on the fitted parameters are reported as columns called <paramname>_unc provided that the fitter object contains a dictionary called fit_info with the key param_cov, which contains the covariance matrix. If param_cov is not present, uncertanties are not reported.

get_residual_image()[source]

Returns an image that is the result of the subtraction between the original image and the fitted sources.

Returns:
residual_image : 2D array-like, ImageHDU, HDUList
nstar(image, star_groups)[source]

Fit, as appropriate, a compound or single model to the given star_groups. Groups are fitted sequentially from the smallest to the biggest. In each iteration, image is subtracted by the previous fitted group.

Parameters:
image : numpy.ndarray

Background-subtracted image.

star_groups : Table

This table must contain the following columns: id, group_id, x_0, y_0, flux_0. x_0 and y_0 are initial estimates of the centroids and flux_0 is an initial estimate of the flux. Additionally, columns named as <param_name>_0 are required if any other parameter in the psf model is free (i.e., the fixed attribute of that parameter is False).

Returns:
result_tab : Table

Astropy table that contains photometry results.

image : numpy.ndarray

Residual image.