fit_fwhm#
- photutils.psf.fit_fwhm(data, *, xypos=None, fwhm=None, fit_shape=None, mask=None, error=None)[source]#
Fit the FWHM of one or more sources in an image.
This convenience function uses a
CircularGaussianPRF
model to fit the sources using thePSFPhotometry
class.Non-finite values (e.g., NaN or inf) in the
data
array are automatically masked.- Parameters:
- data2D array
The 2D array of the image. The input array must be background subtracted.
- xyposarray-like, optional
The initial (x, y) pixel coordinates of the sources. If
None
, then one source will be fit with an initial position using the center-of-mass centroid of thedata
array.- fwhmfloat, optional
The initial guess for the FWHM of the Gaussian PSF model. If
None
, then the initial guess is half the mean of the x and y sizes of thefit_shape
values.- fit_shapeint or tuple of two ints, optional
The shape of the fitting region. If a scalar, then it is assumed to be a square. If
None
, then the shape of the inputdata
will be used.- maskarray-like (bool), optional
A boolean mask with the same shape as the input
data
, where aTrue
value indicates the corresponding element ofdata
is masked.- error2D array, optional
The pixel-wise Gaussian 1-sigma errors of the input
data
.error
is assumed to include all sources of error, including the Poisson error of the sources (seecalc_total_error
) .error
must have the same shape as the inputdata
. If aQuantity
array, thendata
must also be aQuantity
array with the same units.
- Returns:
- fwhm
ndarray
The FWHM of the sources. Note that the returned FWHM values are always positive.
- fwhm
See also
fit_2dgaussian
Fit a 2D Gaussian model to one or more sources in an image.
Notes
The source(s) are fit using the
fit_2dgaussian()
function, which uses aCircularGaussianPRF
model with thePSFPhotometry
class. The initial guess for the flux is the sum of the pixel values within the fitting region. Iffwhm
isNone
, then the initial guess for the FWHM is half the mean of the x and y sizes of thefit_shape
values.Examples
Fit the FWHM of a single source (e.g., a cutout image):
>>> import numpy as np >>> from photutils.psf import CircularGaussianPRF, fit_fwhm >>> yy, xx = np.mgrid[:51, :51] >>> model = CircularGaussianPRF(x_0=22.17, y_0=28.87, fwhm=3.123, flux=9.7) >>> data = model(xx, yy) >>> fwhm = fit_fwhm(data) >>> fwhm array([3.123])
Fit the FWHMs of multiple sources in an image:
>>> import numpy as np >>> from photutils.detection import DAOStarFinder >>> from photutils.psf import (CircularGaussianPRF, fit_fwhm, ... make_psf_model_image) >>> model = CircularGaussianPRF() >>> data, sources = make_psf_model_image((100, 100), model, 5, ... min_separation=25, ... model_shape=(15, 15), ... flux=(100, 200), fwhm=[3, 8]) >>> finder = DAOStarFinder(0.1, 5) >>> finder_tbl = finder(data) >>> xypos = list(zip(sources['x_0'], sources['y_0'])) >>> fwhms = fit_fwhm(data, xypos=xypos, fit_shape=7) >>> fwhms array([5.69467204, 5.21376414, 7.65508658, 3.20255356, 6.66003098])