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
psf_mapMap

the input PSF Map. Should be a Map with 2 non spatial axes. rad and true energy axes should be given in this specific order.

exposure_mapMap

Associated exposure map. Needs to have a consistent map geometry.

Examples

from astropy.coordinates import SkyCoord
from gammapy.maps import WcsGeom, MapAxis
from gammapy.data import Observation
from gammapy.irf import load_cta_irfs
from gammapy.makers import MapDatasetMaker

# Define observation
pointing = SkyCoord("0d", "0d")
irfs = load_cta_irfs("$GAMMAPY_DATA/cta-1dc/caldb/data/cta/1dc/bcf/South_z20_50h/irf_file.fits")
obs = Observation.create(pointing=pointing, irfs=irfs, livetime="1h")

# Create WcsGeom
# 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")

geom = WcsGeom.create(
    binsz=0.25, width="5 deg", skydir=pointing, 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
geom=exposure_geom.upsample(factor=10).drop("rad")
psf_kernel = psf_map.get_psf_kernel(geom=geom)

Attributes Summary

mask_safe_image

Mask safe for the map

psf_map

required_axes

tag

Methods Summary

containment(rad, energy_true[, position])

Containment at given coords

containment_radius(fraction, energy_true[, ...])

Containment at given coords

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 geom.

from_hdulist(hdulist[, hdu, hdu_bands, ...])

Create from HDUList.

get_psf_kernel(geom[, position, max_radius, ...])

Returns a PSF kernel at the given position.

normalize()

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])

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 2-D 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])

Write IRF map to fits

Attributes Documentation

mask_safe_image#

Mask safe for the map

psf_map#
required_axes = ['rad', 'energy_true']#
tag = 'psf_map'#

Methods Documentation

containment(rad, energy_true, position=None)[source]#

Containment at given coords

Parameters
radQuantity

Rad value

energy_trueQuantity

Energy true value

positionSkyCoord

Sky position. By default the center of the map is chosen

Returns
containmentQuantity

Containment values

containment_radius(fraction, energy_true, position=None)[source]#

Containment at given coords

Parameters
fractionfloat

Containment fraction

energy_trueQuantity

Energy true value

positionSkyCoord

Sky position. By default the center of the map is chosen

Returns
containmentQuantity

Containment values

containment_radius_map(energy_true, fraction=0.68)[source]#

Containment radius map.

Parameters
energy_trueQuantity

Energy at which to compute the containment radius

fractionfloat

Containment fraction (range: 0 to 1)

Returns
containment_radius_mapMap

Containment radius map

copy()#

Copy IRF map

cutout(position, width, mode='trim')#

Cutout IRF map.

Parameters
positionSkyCoord

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’}

Mode option for Cutout2D, for details see Cutout2D.

Returns
cutoutIRFMap

Cutout IRF map.

downsample(factor, axis_name=None, weights=None)#

Downsample the spatial dimension by a given factor.

Parameters
factorint

Downsampling factor.

axis_namestr

Which axis to downsample. By default spatial axes are downsampled.

weightsMap

Map with weights downsampling.

Returns
mapIRFMap

Downsampled irf map.

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.

Parameters
energy_axis_trueMapAxis

True energy axis.

rad_axisMapAxis

Offset angle wrt source position axis.

sigmaAngle

Gaussian width.

geomGeom

Image geometry. By default an allsky geometry is created.

Returns
psf_mapPSFMap

Point spread function map.

classmethod from_geom(geom)[source]#

Create psf map from geom.

Parameters
geomGeom

PSF map geometry.

Returns
psf_mapPSFMap

Point spread function map.

classmethod from_hdulist(hdulist, hdu=None, hdu_bands=None, exposure_hdu=None, exposure_hdu_bands=None, format='gadf')#

Create from HDUList.

Parameters
hdulistHDUList

HDU list.

hdustr

Name or index of the HDU with the IRF map.

hdu_bandsstr

Name or index of the HDU with the IRF map BANDS table.

exposure_hdustr

Name or index of the HDU with the exposure map data.

exposure_hdu_bandsstr

Name or index of the HDU with the exposure map BANDS table.

format{“gadf”, “gtpsf”}

File format

Returns
irf_mapIRFMap

IRF map.

get_psf_kernel(geom, position=None, max_radius=None, containment=0.999, factor=4)[source]#

Returns a PSF kernel at the given position.

The PSF is returned in the form a WcsNDMap defined by the input Geom.

Parameters
geomGeom

Target geometry to use

positionSkyCoord

Target position. Should be a single coordinate. By default the center position is used.

max_radiusAngle

maximum angular size of the kernel map

containmentfloat

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.

factorint

oversampling factor to compute the PSF

Returns
kernelPSFKernel

the resulting kernel

normalize()[source]#

Normalize PSF map

peek(figsize=(12, 10))[source]#

Quick-look summary plots.

Parameters
figsizetuple

Size of figure.

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.

Parameters
axAxes

Axes to plot on.

fractionlist of float or ndarray

Containment fraction between 0 and 1.

**kwargsdict

Keyword arguments passed to plot

Returns
axAxes

Axes to plot on.

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.

Parameters
axAxes

Axes to plot on.

energy_trueQuantity

Energies where to plot the PSF.

**kwargsdict

Keyword arguments pass to plot.

Returns
axAxes

Axes to plot on.

classmethod read(filename, format='gadf', hdu=None)#

Read an IRF_map from file and create corresponding object”

Parameters
filenamestr or Path

File name

format{“gadf”, “gtpsf”}

File format

hdustr or int

HDU location

Returns
irf_mapPSFMap, EDispMap or EDispKernelMap

IRF map

sample_coord(map_coord, random_state=0)[source]#

Apply PSF corrections on the coordinates of a set of simulated events.

Parameters
map_coordMapCoord object.

Sequence of coordinates and energies of sampled events.

random_state{int, ‘random-seed’, ‘global-rng’, RandomState}

Defines random number generator initialisation. Passed to get_random_state.

Returns
corr_coordMapCoord object.

Sequence of PSF-corrected coordinates of the input map_coord map.

slice_by_idx(slices)#

Slice sub dataset.

The slicing only applies to the maps that define the corresponding axes.

Parameters
slicesdict

Dict 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 dict are kept unchanged.

Returns
map_outIRFMap

Sliced irf map object.

stack(other, weights=None, nan_to_num=True)#

Stack IRF map with another one in place.

Parameters
otherIRFMap

IRF map to be stacked with this one.

weightsMap

Map with stacking weights.

nan_to_num: bool

Non-finite values are replaced by zero if True (default).

to_hdulist(format='gadf')#

Convert to HDUList.

Parameters
format{“gadf”, “gtpsf”}

File format

Returns
hdu_listHDUList

HDU list.

to_image(spectrum=None, keepdims=True)[source]#

Reduce to a 2-D map after weighing with the associated exposure and a spectrum

Parameters
spectrumSpectralModel, 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

Returns
psf_outPSFMap

PSFMap with the energy axis summed over

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.

Parameters
regionSkyRegion or SkyCoord

Region or position where to get the map.

Returns
irfIRFMap

IRF map with region geometry.

write(filename, overwrite=False, format='gadf')#

Write IRF map to fits

Parameters
filenamestr or Path

Filename to write to

overwritebool

Whether to overwrite

format{“gadf”, “gtpsf”}

File format