PSFMap#

class gammapy.irf.PSFMap[source]#

Bases: 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, 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

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.

from_gauss(energy_axis_true[, rad_axis, ...])

Create all-sky PSF map from Gaussian width.

from_geom(geom)

Create PSF map from geometry.

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

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

sample_coord(map_coord[, random_state, ...])

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

to_image([spectral_model, keepdims])

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

Attributes Documentation

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

Methods Documentation

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

Containment at given coordinates.

Parameters:
radQuantity

Rad value.

energy_trueQuantity

Energy true value.

positionSkyCoord, optional

Sky position. If None, the center of the map is chosen. Default is None.

Returns:
containmentQuantity

Containment values.

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

Containment at given coordinates.

Parameters:
fractionfloat

Containment fraction.

energy_trueQuantity

Energy true value.

positionSkyCoord, optional

Sky position. If None, the center of the map is chosen. Default is None.

Returns:
containmentQuantity

Containment values.

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

Containment radius map.

Parameters:
energy_trueQuantity

Energy true at which to compute the containment radius.

fractionfloat, optional

Containment fraction (range: 0 to 1). Default is 0.68.

Returns:
containment_radius_mapMap

Containment radius 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 all-sky geometry is created.

Returns:
psf_mapPSFMap

Point spread function map.

classmethod from_geom(geom)[source]#

Create PSF map from geometry.

Parameters:
geomGeom

PSF map geometry.

Returns:
psf_mapPSFMap

Point spread function map.

get_psf_kernel(geom, position=None, max_radius=None, containment=0.999, 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:
geomGeom

Target geometry to use.

positionSkyCoord, optional

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

max_radiusAngle, optional

Maximum angular size of the kernel map. Default is None and it will be computed for the containment fraction set.

containmentfloat, optional

Containment fraction to use as size of the kernel. The radius can be overwritten using the max_radius argument. Default is 0.999.

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

Returns:
kernelPSFKernel or list of PSFKernel

The resulting kernel.

normalize()[source]#

Normalize PSF map.

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

Quick-look summary plots.

This method creates a figure with four subplots:

  • Containment radius at center of map plot : Containment radius as a function of energy for containment fractions of 65% and 95%.

  • PSF at center of map plot : PSF vs radius.

  • Exposure 2D map : exposure summed over true energy

  • Containment radius map : 2D sky map of the 68% containment radius at 1 TeV

Parameters:
figsizetuple, optional

Size of figure. Default is (12, 10).

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, optional

Matplotlib axes. Default is None.

fractionlist of float or ndarray

Containment fraction between 0 and 1.

**kwargsdict

Keyword arguments passed to plot

Returns:
axAxes

Matplotlib axes.

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, optional

Matplotlib axes. Default is None.

energy_trueQuantity, optional

Energies where to plot the PSF. Default is [0.1, 1, 10] TeV.

**kwargsdict

Keyword arguments pass to plot.

Returns:
axAxes

Matplotlib axes.

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

Returns:
corr_coordMapCoord

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

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

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

Parameters:
spectral_modelSpectralModel, 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.

__init__(psf_map, exposure_map=None)[source]#
classmethod __new__(*args, **kwargs)#