make_psf_model

photutils.psf.make_psf_model(model, *, x_name=None, y_name=None, flux_name=None, normalize=True, dx=50, dy=50, subsample=100, use_dblquad=False)[source]

Make a PSF model that can be used with the PSF photometry classes (PSFPhotometry or IterativePSFPhotometry) from an Astropy fittable 2D model.

If the x_name, y_name, or flux_name keywords are input, this function will map the those model parameter names to x_0, y_0, or flux, respectively.

If any of these keywords are None, then a new parameter will be added to the model corresponding to the missing parameter. The new position parameters will be set to 0, and the new flux parameter will be set to 1.

The output PSF model will have x_name, y_name, and flux_name attributes that contain the name of the corresponding model parameter.

Note

This function is needed only in cases where the 2D PSF model does not have x_0, y_0, and flux parameters.

It is not needed for any of the PSF models provided by Photutils (e.g., GriddedPSFModel, IntegratedGaussianPRF).

Parameters:
modelFittable2DModel

An Astropy fittable 2D model to use as a PSF.

x_namestr or None, optional

The name of the model parameter that corresponds to the x center of the PSF. If None, the model will be assumed to be centered at x=0, and a new model parameter called xpos_0 will be added for the x position.

y_namestr or None, optional

The name of the model parameter that corresponds to the y center of the PSF. If None, the model will be assumed to be centered at y=0, and a new parameter called ypos_1 will be added for the y position.

flux_namestr or None, optional

The name of the model parameter that corresponds to the total flux of a source. If None, a new model parameter called flux_3 will be added for model flux.

normalizebool, optional

If True, the input model will be integrated and rescaled so that its sum integrates to 1. This normalization occurs only once for the input model. If the total flux of model somehow depends on (x, y) position, then one will need to correct the fitted model fluxes for this effect.

dx, dyodd int, optional

The size of the integration grid in x and y for normalization. Must be odd. These keywords are ignored if normalize is False or use_dblquad is True.

subsampleint, optional

The subsampling factor for the integration grid along each axis for normalization. Each pixel will be sampled subsample x subsample times. This keyword is ignored if normalize is False or use_dblquad is True.

use_dblquadbool, optional

If True, then use scipy.integrate.dblquad to integrate the model for normalization. This is much slower than the default integration of the evaluated model, but it is more accurate. This keyword is ignored if normalize is False.

Returns:
resultCompoundModel

A PSF model that can be used with the PSF photometry classes. The returned model will always be an Astropy compound model.

Notes

To normalize the model, by default it is discretized on a grid of size dx x dy from the model center with a subsampling factor of subsample. The model is then integrated over the grid using trapezoidal integration.

If the use_dblquad keyword is set to True, then the model is integrated using scipy.integrate.dblquad. This is much slower than the default integration of the evaluated model, but it is more accurate. Also, note that the dblquad integration can sometimes fail, e.g., return zero for a non-zero model. This can happen when the model function is sharply localized relative to the size of the integration interval.