EDispMap¶
- 
class gammapy.irf.EDispMap(edisp_map, exposure_map=None)[source]¶
- Bases: - gammapy.irf.core.IRFMap- Energy dispersion map. - Parameters
 - Examples - import numpy as np from astropy import units as u from astropy.coordinates import SkyCoord from gammapy.maps import WcsGeom, MapAxis from gammapy.irf import EnergyDispersion2D, EffectiveAreaTable2D from gammapy.makers.utils import make_edisp_map, make_map_exposure_true_energy # Define energy dispersion map geometry energy_axis = MapAxis.from_edges(np.logspace(-1, 1, 4), unit="TeV", name="energy") migra_axis = MapAxis.from_edges(np.linspace(0, 3, 100), name="migra") pointing = SkyCoord(0, 0, unit="deg") max_offset = 4 * u.deg geom = WcsGeom.create( binsz=0.25 * u.deg, width=10 * u.deg, skydir=pointing, axes=[migra_axis, energy_axis], ) # Extract EnergyDispersion2D from CTA 1DC IRF filename = "$GAMMAPY_DATA/cta-1dc/caldb/data/cta/1dc/bcf/South_z20_50h/irf_file.fits" edisp2D = EnergyDispersion2D.read(filename, hdu="ENERGY DISPERSION") aeff2d = EffectiveAreaTable2D.read(filename, hdu="EFFECTIVE AREA") # Create the exposure map exposure_geom = geom.to_image().to_cube([energy_axis]) exposure_map = make_map_exposure_true_energy(pointing, "1 h", aeff2d, exposure_geom) # create the EDispMap for the specified pointing edisp_map = make_edisp_map(edisp2D, pointing, geom, max_offset, exposure_map) # Get an Energy Dispersion (1D) at any position in the image pos = SkyCoord(2.0, 2.5, unit="deg") energy = np.logspace(-1.0, 1.0, 10) * u.TeV edisp = edisp_map.get_edisp_kernel(pos=pos, energy=energy) # Write map to disk edisp_map.write("edisp_map.fits") - Attributes Summary - Mask safe for the map - Methods Summary - 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_diagonal_response(energy_axis_true[, …])- Create an allsky EDisp map with diagonal response. - from_geom(geom)- Create edisp map from geom. - from_hdulist(hdulist[, hdu, hdu_bands, …])- Create from - HDUList.- get_edisp_kernel([position, energy_axis])- Get energy dispersion at a given position. - Normalize PSF map - read(filename[, format, hdu])- Read an IRF_map from file and create corresponding object” - sample_coord(map_coord[, random_state])- Apply the energy dispersion 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_edisp_kernel_map(energy_axis)- Convert to map with edisp kernels - to_hdulist([format])- Convert to - HDUList.- to_region_nd_map(region)- Extract IRFMap in a given region or position - write(filename[, overwrite, format])- Write IRF map to fits - Attributes Documentation - 
edisp_map¶
 - 
mask_safe_image¶
- Mask safe for the map 
 - 
required_axes= ['migra', 'energy_true']¶
 - 
tag= 'edisp_map'¶
 - Methods Documentation - 
copy()¶
- Copy IRF map 
 - 
cutout(position, width, mode='trim')¶
- Cutout IRF map. - Parameters
- Returns
- cutoutIRFMap
- Cutout IRF map. 
 
- cutout
 
 - 
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. 
 
- map
 
 - 
classmethod from_diagonal_response(energy_axis_true, migra_axis=None)[source]¶
- Create an allsky EDisp map with diagonal response. 
 - 
classmethod from_geom(geom)[source]¶
- Create edisp map from geom. - By default a diagonal edisp matrix is created. - Parameters
- geomGeom
- Edisp map geometry. 
 
- geom
- Returns
- edisp_mapEDispMap
- Energy dispersion map. 
 
- edisp_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 
 
- hdulist
- Returns
- irf_mapIRFMap
- IRF map. 
 
- irf_map
 
 - 
get_edisp_kernel(position=None, energy_axis=None)[source]¶
- Get energy dispersion at a given position. - Parameters
- positionSkyCoord
- the target position. Should be a single coordinates 
- energy_axisMapAxis
- Reconstructed energy axis 
 
- position
- Returns
- edispEnergyDispersion
- the energy dispersion (i.e. rmf object) 
 
- edisp
 
 - 
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 
 
- filenamestr or 
- Returns
- irf_mapPSFMap,EDispMaporEDispKernelMap
- IRF map 
 
- irf_map
 
 - 
sample_coord(map_coord, random_state=0)[source]¶
- Apply the energy dispersion corrections on the coordinates of a set of simulated events. - Parameters
- map_coordMapCoordobject.
- 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.
 
- map_coord
- Returns
- MapCoord.
- Sequence of Edisp-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 - sliceobject 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. 
 
- map_out
 
 - 
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). 
 
- other
 
 - 
to_edisp_kernel_map(energy_axis)[source]¶
- Convert to map with edisp kernels - Parameters
- energy_axisMapAxis
- Reconstructed energy axis. 
 
- energy_axis
- Returns
- edispEDispKernelMap
- Energy dispersion kernel map. 
 
- edisp
 
 - 
to_hdulist(format='gadf')¶
- Convert to - HDUList.- Parameters
- format{“gadf”, “gtpsf”}
- File format 
 
- Returns
- hdu_listHDUList
- HDU list. 
 
- hdu_list
 
 - 
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
- regionSkyRegionorSkyCoord
- Region or position where to get the map. 
 
- region
- Returns
- irfIRFMap
- IRF map with region geometry. 
 
- irf
 
 - 
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 
 
- filenamestr or