PSFPhotometry

class photutils.psf.PSFPhotometry(psf_model, fit_shape, *, finder=None, grouper=None, fitter=<astropy.modeling.fitting.LevMarLSQFitter object>, fitter_maxiters=100, localbkg_estimator=None, aperture_radius=None, progress_bar=False)[source]

Bases: object

Class to perform PSF photometry.

This class implements a flexible PSF photometry algorithm that can find sources in an image, group overlapping sources, fit the PSF model to the sources, and subtract the fit PSF models from the image.

Parameters:
psf_modelastropy.modeling.Fittable2DModel

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).

fit_shapeint or length-2 array_like

The rectangular shape around the center of a star that will be used to define the PSF-fitting data. If fit_shape is a scalar then a square shape of size fit_shape will be used. If fit_shape has two elements, they must be in (ny, nx) order. Each element of fit_shape must be an odd number. In general, fit_shape should be set to a small size (e.g., (5, 5)) that covers the region with the highest flux signal-to-noise.

findercallable or StarFinderBase or None, optional

A callable used to identify stars in an image. The finder must accept a 2D image as input and return a Table containing the x and y centroid positions. These positions are used as the starting points for the PSF fitting. The allowed x column names are (same suffix for y): 'x_init', 'xinit', 'xcentroid', 'x_centroid', 'x_peak', 'x', 'xcen', 'x_cen', 'xpos', 'x_pos', 'x_0', and 'x0'. If None, then the initial (x, y) model positions must be input using the init_params keyword when calling the class. The (x, y) values in init_params override this keyword.

grouperSourceGrouper or callable or None, optional

A callable used to group stars. Typically, grouped stars are those that overlap with their neighbors. Stars that are grouped are fit simultaneously. The grouper must accept the x and y coordinates of the sources and return an integer array of the group id numbers (starting from 1) indicating the group in which a given source belongs. If None, then no grouping is performed, i.e. each source is fit independently. The group_id values in init_params override this keyword only for the first iteration. A warning is raised if any group size is larger than 25 sources.

fitterFitter, optional

The fitter object used to perform the fit of the model to the data.

fitter_maxitersint, optional

The maximum number of iterations in which the fitter is called for each source.

localbkg_estimatorLocalBackground or None, optional

The object used to estimate the local background around each source. If None, then no local background is subtracted. The local_bkg values in init_params override this keyword.

aperture_radiusfloat, optional

The radius of the circular aperture used to estimate the initial flux of each source. The flux_init values in init_params override this keyword.

progress_barbool, optional

Whether to display a progress bar when fitting the sources (or groups). The progress bar requires that the tqdm optional dependency be installed. Note that the progress bar does not currently work in the Jupyter console due to limitations in tqdm.

Methods Summary

__call__(data, *[, mask, error, init_params])

Perform PSF photometry.

make_model_image(shape, psf_shape, *[, ...])

Create a 2D image from the fit PSF models and optional local background.

make_residual_image(data, psf_shape, *[, ...])

Create a 2D residual image from the fit PSF models and local background.

Methods Documentation

__call__(data, *, mask=None, error=None, init_params=None)[source]

Perform PSF photometry.

Parameters:
data2D ndarray

The 2D array on which to perform photometry. Invalid data values (i.e., NaN or inf) are automatically masked.

mask2D bool ndarray, optional

A boolean mask with the same shape as data, where a True value indicates the corresponding element of data is masked.

error2D ndarray, optional

The pixel-wise 1-sigma errors of the input data. error is assumed to include all sources of error, including the Poisson error of the sources (see calc_total_error) . error must have the same shape as the input data. If a Quantity array, then data must also be a Quantity array with the same units.

init_paramsTable or None, optional

A table containing the initial guesses of the (x, y, flux) model parameters for each source. If the x and y values are not input, then the finder keyword must be defined. If the flux values are not input, then the aperture_radius keyword must be defined. Note that the initial flux values refer to the model flux parameters and are not corrected for local background values (computed using localbkg_estimator or input in a local_bkg column) The allowed column names are:

  • x_init, xinit, xcentroid, x_centroid, x_peak, x, xcen, x_cen, xpos, x_pos, x_0, and x0.

  • y_init, yinit, ycentroid, y_centroid, y_peak, y, ycen, y_cen, ypos, y_pos, y_0, and y0.

  • flux_init, flux, source_sum, segment_flux, and kron_flux.

The parameter names are searched in the input table in the above order, stopping at the first match.

The table can also have group_id and local_bkg columns. If group_id is input, the values will be used and grouper keyword will be ignored. If local_bkg is input, they will be used and the localbkg_estimator will be ignored.

Returns:
tableQTable

An astropy table with the PSF-fitting results. The table will contain the following columns:

  • id : unique identification number for the source

  • group_id : unique identification number for the source group

  • x_init, x_fit, x_err : the initial, fit, and error of the source x center

  • y_init, y_fit, y_err : the initial, fit, and error of the source y center

  • flux_init, flux_fit, flux_err : the initial, fit, and error of the source flux

  • npixfit : the number of unmasked pixels used to fit the source

  • group_size : the total number of sources that were simultaneously fit along with the given source

  • qfit : a quality-of-fit metric defined as the absolute value of the sum of the fit residuals divided by the fit flux

  • cfit : a quality-of-fit metric defined as the fit residual in the central pixel divided by the fit flux

  • flagsbitwise flag values:
    • 1 : one or more pixels in the fit_shape region were masked

    • 2 : the fit x and/or y position lies outside of the input data

    • 4 : the fit flux is less than or equal to zero

    • 8 : the fitter may not have converged

    • 16 : the fitter parameter covariance matrix was not returned

make_model_image(shape, psf_shape, *, include_localbkg=False)[source]

Create a 2D image from the fit PSF models and optional local background.

Parameters:
shape2 tuple of int

The shape of the output array.

psf_shape2 tuple of int

The shape of region around the center of the fit model to render in the output image.

include_localbkgbool, optional

Whether to include the local background in the rendered output image. Note that the local background level is included around each source over the region defined by psf_shape. Thus, regions where the psf_shape of sources overlap will have the local background added multiple times.

Returns:
array2D ndarray

The rendered image from the fit PSF models.

make_residual_image(data, psf_shape, *, include_localbkg=False)[source]

Create a 2D residual image from the fit PSF models and local background.

Parameters:
data2D ndarray

The 2D array on which photometry was performed. This should be the same array input when calling the PSF-photometry class.

psf_shape2 tuple of int

The shape of region around the center of the fit model to subtract.

include_localbkgbool, optional

Whether to include the local background in the subtracted model. Note that the local background level is subtracted around each source over the region defined by psf_shape. Thus, regions where the psf_shape of sources overlap will have the local background subtracted multiple times.

Returns:
array2D ndarray

The residual image of the data minus the local_bkg minus the fit PSF models.