RadialProfile#
- class photutils.profiles.RadialProfile(data, xycen, radii, *, error=None, mask=None, method='exact', subpixels=5)[source]#
Bases:
ProfileBaseClass to create a radial profile using concentric circular annulus apertures.
The radial profile represents the azimuthally-averaged flux in circular annuli apertures as a function of radius.
For this class, the input radii represent the edges of the radial bins. This differs from the
CurveOfGrowthclass, where the input radii are the radii of the circular apertures used to compute the cumulative flux.The output
radiusattribute represents the bin centers.- Parameters:
- data2D
ndarray The 2D data array. The data should be background-subtracted. Non-finite values (e.g., NaN or inf) in the
dataorerrorarray are automatically masked.- xycentuple of 2 floats
The
(x, y)pixel coordinate of the source center.- radii1D float
ndarray An array of radii defining the edges of the radial bins.
radiimust be strictly increasing with a minimum value greater than or equal to zero, and contain at least 2 values. The radial spacing does not need to be constant. The outputradiusattribute will be defined at the bin centers.- error2D
ndarray, optional The 1-sigma errors of the input
data.erroris assumed to include all sources of error, including the Poisson error of the sources (seecalc_total_error).errormust have the same shape as the inputdata. Non-finite values (e.g., NaN or inf) in thedataorerrorarray are automatically masked.- mask2D bool
ndarray, optional A boolean mask with the same shape as
datawhere aTruevalue indicates the corresponding element ofdatais masked. Masked data are excluded from all calculations.- method{‘exact’, ‘center’, ‘subpixel’}, optional
The method used to determine the overlap of the aperture on the pixel grid:
'exact'(default): The exact fractional overlap of the aperture and each pixel is calculated. The aperture weights will contain values between 0 and 1.'center': A pixel is considered to be entirely in or out of the aperture depending on whether its center is in or out of the aperture. The aperture weights will contain values only of 0 (out) and 1 (in).'subpixel': A pixel is divided into subpixels (see thesubpixelskeyword), each of which are considered to be entirely in or out of the aperture depending on whether its center is in or out of the aperture. Ifsubpixels=1, this method is equivalent to'center'. The aperture weights will contain values between 0 and 1.
- subpixelsint, optional
For the
'subpixel'method, resample pixels by this factor in each dimension. That is, each pixel is divided intosubpixels**2subpixels. This keyword is ignored unlessmethod='subpixel'.
- data2D
See also
Notes
If the minimum of
radiiis zero, then a circular aperture with radius equal toradii[1]will be used for the innermost aperture.Examples
>>> import numpy as np >>> from astropy.modeling.models import Gaussian2D >>> from photutils.centroids import centroid_2dg >>> from photutils.datasets import make_noise_image >>> from photutils.profiles import RadialProfile
Create an artificial single source. Note that this image does not have any background.
>>> gmodel = Gaussian2D(42.1, 47.8, 52.4, 4.7, 4.7, 0) >>> yy, xx = np.mgrid[0:100, 0:100] >>> data = gmodel(xx, yy) >>> bkg_sig = 2.1 >>> noise = make_noise_image(data.shape, mean=0., stddev=bkg_sig, seed=0) >>> data += noise >>> error = np.zeros_like(data) + bkg_sig
Create the radial profile.
>>> xycen = centroid_2dg(data) >>> edge_radii = np.arange(25) >>> rp = RadialProfile(data, xycen, edge_radii, error=error)
>>> print(rp.radius) [ 0.5 1.5 2.5 3.5 4.5 5.5 6.5 7.5 8.5 9.5 10.5 11.5 12.5 13.5 14.5 15.5 16.5 17.5 18.5 19.5 20.5 21.5 22.5 23.5]
>>> print(rp.profile) [ 4.30187860e+01 4.02502046e+01 3.57758011e+01 3.16071235e+01 2.61511082e+01 2.10539746e+01 1.63701300e+01 1.16674718e+01 8.12828014e+00 5.78962699e+00 3.59342666e+00 2.35353336e+00 1.20355937e+00 7.67093923e-01 4.24650784e-01 8.67989701e-02 5.11484374e-02 -9.82041768e-02 2.37482124e-02 -3.66602855e-02 6.84802299e-02 1.72239596e-01 -3.86056497e-02 7.30423743e-02]
>>> print(rp.profile_error) [1.18479813 0.68404352 0.52985783 0.4478116 0.39493271 0.35723008 0.32860388 0.30591356 0.28735575 0.27181133 0.25854415 0.24704749 0.23695963 0.22801451 0.22001149 0.21279603 0.20624688 0.20026744 0.19477961 0.18971954 0.18503438 0.18068002 0.17661928 0.17282057]
Plot the radial profile, including the raw data profile.
(
Source code,png,hires.png,pdf,svg)
Normalize the profile and plot the normalized radial profile.
(
Source code,png,hires.png,pdf,svg)
Plot three of the annulus apertures on the data.
(
Source code,png,hires.png,pdf,svg)
Fit 1D Gaussian and Moffat models to the normalized radial profile.
>>> rp.normalize() >>> rp.gaussian_fit <Gaussian1D(amplitude=0.98231088, mean=0., stddev=4.67512776)> >>> rp.moffat_fit <Moffat1D(amplitude=..., x_0=0., gamma=..., alpha=...)> >>> print(rp.gaussian_fwhm) 11.009084813327846 >>> print(rp.moffat_fwhm) 10.868426507928344
Plot the fitted 1D Gaussian and Moffat models on the radial profile.
(
Source code,png,hires.png,pdf,svg)
Attributes Summary
A list of the circular annulus apertures used to measure the radial profile, as
CircularAnnulusobjects.The unmasked area in each circular annulus (or aperture) as a function of radius as a 1D
ndarray.The raw data profile as a 1D
ndarray.The radii of the raw data profile as a 1D
ndarray.The fitted 1D Gaussian to the radial profile as a
Gaussian1Dmodel.The full-width at half-maximum (FWHM) in pixels of the 1D Gaussian fitted to the radial profile.
The fitted 1D Gaussian profile to the radial profile as a 1D
ndarray.The fitted 1D Moffat to the radial profile as a
Moffat1Dmodel.The full-width at half-maximum (FWHM) in pixels of the 1D Moffat fitted to the radial profile.
The fitted 1D Moffat profile to the radial profile as a 1D
ndarray.The radial profile as a 1D
ndarray.The radial profile errors as a 1D
ndarray.The profile radius (bin centers) in pixels as a 1D
ndarray.Methods Summary
normalize([method])Normalize the profile.
plot([ax])Plot the profile.
plot_error([ax])Plot the profile errors.
Unnormalize the profile back to the original state before any calls to
normalize.Attributes Documentation
- apertures#
A list of the circular annulus apertures used to measure the radial profile, as
CircularAnnulusobjects.If the minimum of
radiiis zero, then the innermost element will be aCircularAperturewith radius equal toradii[1].
- area#
The unmasked area in each circular annulus (or aperture) as a function of radius as a 1D
ndarray.
- gaussian_fit#
The fitted 1D Gaussian to the radial profile as a
Gaussian1Dmodel.The cached fit is automatically invalidated when the profile normalization is changed, so the fit is always consistent with the current profile.
- gaussian_fwhm#
The full-width at half-maximum (FWHM) in pixels of the 1D Gaussian fitted to the radial profile.
The cached value is automatically invalidated when the profile normalization is changed.
Returns
Noneif the fit failed (e.g., the profile is entirely non-finite or masked).
- gaussian_profile#
The fitted 1D Gaussian profile to the radial profile as a 1D
ndarray.The cached profile is automatically invalidated when the profile normalization is changed.
Returns
Noneif the fit failed (e.g., the profile is entirely non-finite or masked).
- moffat_fit#
The fitted 1D Moffat to the radial profile as a
Moffat1Dmodel.The cached fit is automatically invalidated when the profile normalization is changed, so the fit is always consistent with the current profile.
- moffat_fwhm#
The full-width at half-maximum (FWHM) in pixels of the 1D Moffat fitted to the radial profile.
The cached value is automatically invalidated when the profile normalization is changed.
Returns
Noneif the fit failed (e.g., the profile is entirely non-finite or masked).
- moffat_profile#
The fitted 1D Moffat profile to the radial profile as a 1D
ndarray.The cached profile is automatically invalidated when the profile normalization is changed.
Returns
Noneif the fit failed (e.g., the profile is entirely non-finite or masked).
- profile_error#
The radial profile errors as a 1D
ndarray.If no
errorarray was provided, an empty array with shape(0,)is returned.
- radius#
The profile radius (bin centers) in pixels as a 1D
ndarray.The returned radius values are defined as the arithmetic means of the input radial-bins edges (
radii).For logarithmically-spaced input
radii, one could instead use a radius array defined using the geometric mean of the bin edges, i.e.np.sqrt(radii[:-1] * radii[1:]).
Methods Documentation
- normalize(method='max')#
Normalize the profile.
- Parameters:
- method{‘max’, ‘sum’}, optional
The method used to normalize the profile:
'max'(default): The profile is normalized such that its maximum value is 1.'sum': The profile is normalized such that its sum (integral) is 1.
- plot(ax=None, **kwargs)#
Plot the profile.
- Parameters:
- ax
matplotlib.axes.AxesorNone, optional The matplotlib axes on which to plot. If
None, then the currentAxesinstance is used.- **kwargsdict, optional
Any keyword arguments accepted by
matplotlib.pyplot.plot.
- ax
- Returns:
- lineslist of
Line2D A list of lines representing the plotted data.
- lineslist of
- plot_error(ax=None, **kwargs)#
Plot the profile errors.
- Parameters:
- ax
matplotlib.axes.AxesorNone, optional The matplotlib axes on which to plot. If
None, then the currentAxesinstance is used.- **kwargsdict, optional
Any keyword arguments accepted by
matplotlib.pyplot.fill_between.
- ax
- Returns:
- poly
matplotlib.collections.PolyCollectionorNone A
PolyCollectioncontaining the plotted polygons, orNoneif no errors were input.
- poly