PSFMap#
- class gammapy.irf.PSFMap(psf_map, exposure_map=None)[source]#
Bases:
gammapy.irf.core.IRFMap
Class containing the Map of PSFs and allowing to interact with it.
- Parameters
Examples
from astropy.coordinates import SkyCoord from gammapy.maps import WcsGeom, MapAxis from gammapy.data import Observation, FixedPointingInfo from gammapy.irf import load_irf_dict_from_file from gammapy.makers import MapDatasetMaker # Define observation pointing_position = SkyCoord(0, 0, unit="deg", frame="galactic") pointing = FixedPointingInfo( fixed_icrs=pointing_position.icrs, ) filename = "$GAMMAPY_DATA/cta-1dc/caldb/data/cta/1dc/bcf/South_z20_50h/irf_file.fits" irfs = load_irf_dict_from_file(filename) obs = Observation.create(pointing=pointing, irfs=irfs, livetime="1h") # Define energy axis. Note that the name is fixed. energy_axis = MapAxis.from_energy_bounds("0.1 TeV", "10 TeV", nbin=3, name="energy_true") # Define rad axis. Again note the axis name rad_axis = MapAxis.from_bounds(0, 0.5, nbin=100, name="rad", unit="deg") # Create WcsGeom geom = WcsGeom.create( binsz=0.25, width="5 deg", skydir=pointing_position, axes=[rad_axis, energy_axis] ) maker = MapDatasetMaker() psf = maker.make_psf(geom=geom, observation=obs) # Get a PSF kernel at the center of the image upsample_geom = geom.upsample(factor=10).drop("rad") psf_kernel = psf.get_psf_kernel(geom=upsample_geom)
Attributes Summary
Mask safe for the map.
Methods Summary
containment
(rad, energy_true[, position])Containment at given coordinates.
containment_radius
(fraction, energy_true[, ...])Containment at given coordinates.
containment_radius_map
(energy_true[, fraction])Containment radius map.
copy
()Copy IRF map.
cutout
(position, width[, mode])Cutout IRF map.
downsample
(factor[, axis_name, weights])Downsample the spatial dimension by a given factor.
from_gauss
(energy_axis_true[, rad_axis, ...])Create all-sky PSF map from Gaussian width.
from_geom
(geom)Create PSF map from geometry.
from_hdulist
(hdulist[, hdu, hdu_bands, ...])Create from
HDUList
.get_psf_kernel
(geom[, position, max_radius, ...])Return a PSF kernel at the given position.
Normalize PSF map.
peek
([figsize])Quick-look summary plots.
plot_containment_radius_vs_energy
([ax, fraction])Plot containment fraction as a function of energy.
plot_psf_vs_rad
([ax, energy_true])Plot PSF vs radius.
read
(filename[, format, hdu, checksum])Read an IRF_map from file and create corresponding object.
sample_coord
(map_coord[, random_state, ...])Apply PSF corrections on the coordinates of a set of simulated events.
slice_by_idx
(slices)Slice sub dataset.
stack
(other[, weights, nan_to_num])Stack IRF map with another one in place.
to_hdulist
([format])Convert to
HDUList
.to_image
([spectrum, keepdims])Reduce to a 2D map after weighing with the associated exposure and a spectrum.
to_region_nd_map
(region)Extract IRFMap in a given region or position.
write
(filename[, overwrite, format, checksum])Write IRF map to fits.
Attributes Documentation
- energy_name#
- mask_safe_image#
Mask safe for the map.
- psf_map#
- required_axes = ['rad', 'energy_true']#
- tag = 'psf_map'#
Methods Documentation
- copy()#
Copy IRF map.
- cutout(position, width, mode='trim')#
Cutout IRF map.
- Parameters
- position
SkyCoord
Center position of the cutout region.
- widthtuple of
Angle
Angular sizes of the region in (lon, lat) in that specific order. If only one value is passed, a square region is extracted.
- mode{‘trim’, ‘partial’, ‘strict’}, optional
Mode option for Cutout2D, for details see
Cutout2D
. Default is “trim”.
- position
- Returns
- cutout
IRFMap
Cutout IRF map.
- cutout
- downsample(factor, axis_name=None, weights=None)#
Downsample the spatial dimension by a given factor.
- classmethod from_gauss(energy_axis_true, rad_axis=None, sigma=<Quantity 0.1 deg>, geom=None)[source]#
Create all-sky PSF map from Gaussian width.
This is used for testing and examples.
The width can be the same for all energies or be an array with one value per energy node. It does not depend on position.
- classmethod from_geom(geom)[source]#
Create PSF map from geometry.
- Parameters
- geom
Geom
PSF map geometry.
- geom
- Returns
- psf_map
PSFMap
Point spread function map.
- psf_map
- classmethod from_hdulist(hdulist, hdu=None, hdu_bands=None, exposure_hdu=None, exposure_hdu_bands=None, format='gadf')#
Create from
HDUList
.- Parameters
- hdulist
HDUList
HDU list.
- hdustr, optional
Name or index of the HDU with the IRF map. Default is None.
- hdu_bandsstr, optional
Name or index of the HDU with the IRF map BANDS table. Default is None.
- exposure_hdustr, optional
Name or index of the HDU with the exposure map data. Default is None.
- exposure_hdu_bandsstr, optional
Name or index of the HDU with the exposure map BANDS table. Default is None.
- format{“gadf”, “gtpsf”}, optional
File format. Default is “gadf”.
- hdulist
- Returns
- irf_map
IRFMap
IRF map.
- irf_map
- get_psf_kernel(geom, position=None, max_radius=None, containment=0.999, factor=None, precision_factor=12)[source]#
Return a PSF kernel at the given position.
The PSF is returned in the form a WcsNDMap defined by the input Geom.
- Parameters
- geom
Geom
Target geometry to use.
- position
SkyCoord
, optional Target position. Should be a single coordinate. By default, the center position is used.
- max_radius
Angle
, optional Maximum angular size of the kernel map.
- containmentfloat, optional
Containment fraction to use as size of the kernel. The max. radius across all energies is used. The radius can be overwritten using the
max_radius
argument. Default is 0.999.- factorint, optional
Oversampling factor to compute the PSF. Default is None and it will be computed automatically.
- precision_factorint, optional
Factor between the bin half-width of the geom and the median R68% containment radius. Used only for the oversampling method. Default is 10.
- geom
- Returns
- kernel
PSFKernel
The resulting kernel.
- kernel
- plot_containment_radius_vs_energy(ax=None, fraction=(0.68, 0.95), **kwargs)[source]#
Plot containment fraction as a function of energy.
The method plots the containment radius at the center of the map.
- plot_psf_vs_rad(ax=None, energy_true=<Quantity [ 0.1, 1., 10. ] TeV>, **kwargs)[source]#
Plot PSF vs radius.
The method plots the profile at the center of the map.
- classmethod read(filename, format='gadf', hdu=None, checksum=False)#
Read an IRF_map from file and create corresponding object.
- Parameters
- filenamestr or
Path
File name.
- format{“gadf”, “gtpsf”}, optional
File format. Default is “gadf”.
- hdustr or int
HDU location. Default is None.
- checksumbool
If True checks both DATASUM and CHECKSUM cards in the file headers. Default is False.
- filenamestr or
- Returns
- irf_map
PSFMap
,EDispMap
orEDispKernelMap
IRF map.
- irf_map
- sample_coord(map_coord, random_state=0, chunk_size=10000)[source]#
Apply PSF corrections on the coordinates of a set of simulated events.
- Parameters
- map_coord
MapCoord
object. Sequence of coordinates and energies of sampled events.
- random_state{int, ‘random-seed’, ‘global-rng’,
RandomState
}, optional Defines random number generator initialisation. Passed to
get_random_state
. Default is 0.- chunk_sizeint
If set, this will slice the input MapCoord into smaller chunks of chunk_size elements. Default is 10000.
- map_coord
- Returns
- corr_coord
MapCoord
Sequence of PSF-corrected coordinates of the input map_coord map.
- corr_coord
- slice_by_idx(slices)#
Slice sub dataset.
The slicing only applies to the maps that define the corresponding axes.
- Parameters
- slicesdict
Dictionary of axes names and integers or
slice
object pairs. Contains one element for each non-spatial dimension. For integer indexing the corresponding axes is dropped from the map. Axes not specified in the dictionary are kept unchanged.
- Returns
- map_out
IRFMap
Sliced IRF map object.
- map_out
- stack(other, weights=None, nan_to_num=True)#
Stack IRF map with another one in place.
- to_hdulist(format='gadf')#
Convert to
HDUList
.- Parameters
- format{“gadf”, “gtpsf”}, optional
File format. Default is “gadf”.
- Returns
- hdu_list
HDUList
HDU list.
- hdu_list
- to_image(spectrum=None, keepdims=True)[source]#
Reduce to a 2D map after weighing with the associated exposure and a spectrum.
- Parameters
- spectrum
SpectralModel
, optional Spectral model to compute the weights. Default is power-law with spectral index of 2.
- keepdimsbool, optional
If True, the energy axis is kept with one bin. If False, the axis is removed.
- spectrum
- Returns
- to_region_nd_map(region)#
Extract IRFMap in a given region or position.
If a region is given a mean IRF is computed, if a position is given the IRF is interpolated.
- write(filename, overwrite=False, format='gadf', checksum=False)#
Write IRF map to fits.
- Parameters
- filenamestr or
Path
Filename to write to.
- overwritebool, optional
Overwrite existing file. Default is False.
- format{“gadf”, “gtpsf”}, optional
File format. Default is “gadf”.
- checksumbool, optional
When True adds both DATASUM and CHECKSUM cards to the headers written to the file. Default is False.
- filenamestr or