radial_profile¶
-
gammapy.image.radial_profile(image, center, radius)[source]¶ Image radial profile.
TODO: show example and explain handling of “overflow” and “underflow” bins (see
radial_profile_label_imagedocstring).Calls
numpy.digitizeto compute a label image and thenscipy.ndimage.sumto do measurements.Parameters: image :
SkyImageImage
center :
SkyCoordCenter position
radius :
AngleOffset bin edge array.
Returns: table :
TableTable with the following fields that define the binning:
RADIUS_BIN_ID: Integer bin ID (starts at1).RADIUS_MEAN: Radial bin centerRADIUS_MIN: Radial bin minimum edgeRADIUS_MAX: Radial bin maximum edge
And the following measurements from the pixels in each bin:
N_PIX: Number of pixelsSUM: Sum of pixel valuesMEAN: Mean of pixel values, computed asSUM / N_PIX
Examples
Make some example data:
from astropy.coordinates import Angle from gammapy.image import SkyImage image = SkyImage.empty() image.fill(value=1) center = image.center radius = Angle([0.1, 0.2, 0.4, 0.5, 1.0], 'deg')
Compute and print a radial profile:
from gammapy.image import radial_profile table = radial_profile(image, center, radius) table.pprint()
If your measurement represents counts, you could e.g. use this method to compute errors:
import numpy as np table['SUM_ERR'] = np.sqrt(table['SUM']) table['MEAN_ERR'] = table['SUM_ERR'] / table['N_PIX']
If you need to do special measurements or error computation in each bin with access to the pixel values, you could get the label image and then do the measurements yourself:
labels = radial_profile_label_image(image, center, radius) labels.show()