fit_2dgaussian#
- photutils.psf.fit_2dgaussian(data, *, xypos=None, fwhm=None, fix_fwhm=True, fit_shape=None, mask=None, error=None)[source]#
Fit a 2D Gaussian model to 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.- fix_fwhmbool, optional
Whether to fix the FWHM of the Gaussian PSF model during the fitting process.
- 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:
- result
PSFPhotometry
The PSF-fitting photometry results.
- result
See also
fit_fwhm
Fit the FWHM of one or more sources in an image.
Notes
The source(s) are fit with a
CircularGaussianPRF
model using 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 a 2D Gaussian model to a image containing only one source (e.g., a cutout image):
>>> import numpy as np >>> from photutils.psf import CircularGaussianPRF, fit_2dgaussian >>> 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) >>> fit = fit_2dgaussian(data, fix_fwhm=False) >>> phot_tbl = fit.results >>> cols = ['x_fit', 'y_fit', 'fwhm_fit', 'flux_fit'] >>> for col in cols: ... phot_tbl[col].info.format = '.4f' # optional format >>> print(phot_tbl[['id'] + cols]) id x_fit y_fit fwhm_fit flux_fit --- ------- ------- -------- -------- 1 22.1700 28.8700 3.1230 9.7000
Fit a 2D Gaussian model to multiple sources in an image:
>>> import numpy as np >>> from photutils.detection import DAOStarFinder >>> from photutils.psf import (CircularGaussianPRF, fit_2dgaussian, ... 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'])) >>> psfphot = fit_2dgaussian(data, xypos=xypos, fit_shape=7, ... fix_fwhm=False) >>> phot_tbl = psfphot.results >>> len(phot_tbl) 5
Here we show only a few columns of the photometry table:
>>> cols = ['x_fit', 'y_fit', 'fwhm_fit', 'flux_fit'] >>> for col in cols: ... phot_tbl[col].info.format = '.4f' # optional format >>> print(phot_tbl[['id'] + cols]) id x_fit y_fit fwhm_fit flux_fit --- ------- ------- -------- -------- 1 61.7787 74.6905 5.6947 147.9988 2 30.2017 27.5858 5.2138 123.2373 3 10.5237 82.3776 7.6551 180.1881 4 8.4214 12.0369 3.2026 192.3530 5 76.9412 35.9061 6.6600 126.6130