PSFKernel

class gammapy.cube.PSFKernel(psf_kernel_map)[source]

Bases: object

PSF kernel for Map.

This is a container class to store a PSF kernel that can be used to convolve WcsNDMap objects. It is usually computed from an EnergyDependentTablePSF.

Parameters:

psf_kernel_map : Map

PSF kernel stored in a Map

Examples

import numpy as np
from gammapy.maps import Map, WcsGeom, MapAxis
from gammapy.irf import EnergyDependentMultiGaussPSF
from gammapy.cube import PSFKernel
from astropy import units as u

# Define energy axis
energy_axis = MapAxis.from_edges(np.logspace(-1., 1., 4), unit='TeV', name='energy')

# Create WcsGeom and map
geom = WcsGeom.create(binsz=0.02*u.deg, width=2.0*u.deg, axes=[energy_axis])
some_map = Map.from_geom(geom)
# Fill map at two positions
some_map.fill_by_coord([[0.2,0.4],[-0.1,0.6],[0.5,3.6]])

# Extract EnergyDependentTablePSF from CTA 1DC IRF
filename = '$GAMMAPY_DATA/cta-1dc/caldb/data/cta/1dc/bcf/South_z20_50h/irf_file.fits'
psf = EnergyDependentMultiGaussPSF.read(filename, hdu='POINT SPREAD FUNCTION')
table_psf = psf.to_energy_dependent_table_psf(theta=0.5*u.deg)

psf_kernel = PSFKernel.from_table_psf(table_psf,geom, max_radius=1*u.deg)

# Do the convolution
some_map_convolved = some_map.convolve(psf_kernel)

some_map_convolved.get_image_by_coord(dict(energy=0.6*u.TeV)).plot()

Attributes Summary

data Access the PSFKernel numpy array
psf_kernel_map The map object holding the kernel (Map)

Methods Summary

from_gauss(geom, sigma[, max_radius, …]) Create Gaussian PSF.
from_table_psf(table_psf, geom[, …]) Create a PSF kernel from a TablePSF or an EnergyDependentTablePSF on a given MapGeom.
read(*args, **kwargs) Read kernel Map from file.
write(*args, **kwargs) Write the Map object which contains the PSF kernel to file.

Attributes Documentation

data

Access the PSFKernel numpy array

psf_kernel_map

The map object holding the kernel (Map)

Methods Documentation

classmethod from_gauss(geom, sigma, max_radius=None, containment_fraction=0.99, factor=4)[source]

Create Gaussian PSF.

This is used for testing and examples. The map geometry parameters (pixel size, energy bins) are taken from geom. The Gaussian width sigma is a scalar,

TODO : support array input if it should vary along the energy axis.

Parameters:

geom : WcsGeom

Map geometry

sigma : Angle

Gaussian width.

max_radius : Angle

Desired kernel map size.

factor : int

Oversample factor to compute the PSF

Returns:

kernel : PSFKernel

the kernel Map with reduced geometry according to the max_radius

classmethod from_table_psf(table_psf, geom, max_radius=None, factor=4)[source]

Create a PSF kernel from a TablePSF or an EnergyDependentTablePSF on a given MapGeom.

If the MapGeom is not an image, the same kernel will be present on all axes.

The PSF is estimated by oversampling defined by a given factor.

Parameters:

table_psf : TablePSF or EnergyDependentTablePSF

the input table PSF

geom : WcsGeom

the target geometry. The PSF kernel will be centered on the central pixel. The geometry axes should contain an axis with name “energy”

max_radius : Angle

the maximum radius of the PSF kernel.

factor : int

the oversample factor to compute the PSF

Returns:

kernel : PSFKernel

the kernel Map with reduced geometry according to the max_radius

classmethod read(*args, **kwargs)[source]

Read kernel Map from file.

write(*args, **kwargs)[source]

Write the Map object which contains the PSF kernel to file.