PSFMap¶
-
class
gammapy.cube.
PSFMap
(psf_map)[source]¶ Bases:
object
Class containing the Map of PSFs and allowing to interact with it.
Parameters: psf_map :
Map
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.
Examples
import numpy as np from gammapy.maps import Map, WcsGeom, MapAxis from gammapy.irf import EnergyDependentMultiGaussPSF from gammapy.cube import make_psf_map, PSFMap from astropy import units as u from astropy.coordinates import SkyCoord # Define energy axis. Note that the name is fixed. energy_axis = MapAxis.from_edges(np.logspace(-1., 1., 4), unit='TeV', name='energy') # Define rad axis. Again note the axis name rads = np.linspace(0., 0.5, 100) * u.deg rad_axis = MapAxis.from_edges(rads, unit='deg', name='theta') # Define parameters pointing = SkyCoord(0., 0., unit='deg') max_offset = 4 * u.deg # Create WcsGeom geom = WcsGeom.create(binsz=0.25*u.deg, width=10*u.deg, skydir=pointing, axes=[rad_axis, energy_axis]) # 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') psf3d = psf.to_psf3d(rads) # create the PSFMap for the specified pointing psf_map = make_psf_map(psf3d, pointing, geom, max_offset) # Get an EnergyDependentTablePSF at any position in the image psf_table = psf_map.get_energy_dependent_table_psf(SkyCoord(2., 2.5, unit='deg')) # Write map to disk psf_map.write('psf_map.fits')
Attributes Summary
data
the PSFMap data geom
The PSFMap MapGeom object psf_map
the PSFMap itself ( Map
)quantity
the PSFMap data as a quantity Methods Summary
containment_radius_map
(energy[, fraction])Containment radius map. get_energy_dependent_table_psf
(position)Returns EnergyDependentTable PSF at a given position get_psf_kernel
(position, geom[, max_radius, …])Returns a PSF kernel at the given position. read
(filename, **kwargs)Read a psf_map from file and create a PSFMap object write
(*args, **kwargs)Write the Map object containing the PSF Library map. Attributes Documentation
-
data
¶ the PSFMap data
-
geom
¶ The PSFMap MapGeom object
-
quantity
¶ the PSFMap data as a quantity
Methods Documentation
-
containment_radius_map
(energy, fraction=0.68)[source]¶ Containment radius map.
Parameters: energy :
Quantity
Scalar energy at which to compute the containment radius
fraction : float
the containment fraction (range: 0 to 1)
Returns: containment_radius_map :
Map
Containment radius map
-
get_energy_dependent_table_psf
(position)[source]¶ Returns EnergyDependentTable PSF at a given position
Parameters: position :
SkyCoord
the target position. Should be a single coordinates
Returns: psf_table :
EnergyDependentTablePSF
the table PSF
-
get_psf_kernel
(position, geom, max_radius=None, factor=4)[source]¶ Returns a PSF kernel at the given position.
The PSF is returned in the form a WcsNDMap defined by the input MapGeom.
Parameters: position :
SkyCoord
the target position. Should be a single coordinate
geom :
MapGeom
the target geometry to use
max_radius :
Angle
maximum angular size of the kernel map
factor : int
oversampling factor to compute the PSF
Returns: kernel :
PSFKernel
the resulting kernel
-