MapDatasetOnOff#
- class gammapy.datasets.MapDatasetOnOff(models=None, counts=None, counts_off=None, acceptance=None, acceptance_off=None, exposure=None, mask_fit=None, psf=None, edisp=None, name=None, mask_safe=None, gti=None, meta_table=None)[source]#
Bases:
gammapy.datasets.map.MapDataset
Map dataset for on-off likelihood fitting. Uses wstat statistics.
- Parameters
- models
Models
Source sky models.
- counts
WcsNDMap
Counts cube
- counts_off
WcsNDMap
Ring-convolved counts cube
- acceptance
WcsNDMap
Acceptance from the IRFs
- acceptance_off
WcsNDMap
Acceptance off
- exposure
WcsNDMap
Exposure cube
- mask_fit
WcsNDMap
Mask to apply to the likelihood for fitting.
- psf
PSFKernel
PSF kernel
- edisp
EDispKernel
Energy dispersion
- mask_safe
WcsNDMap
Mask defining the safe data range.
- gti
GTI
GTI of the observation or union of GTI if it is a stacked observation
- meta_table
Table
Table listing information on observations used to create the dataset. One line per observation for stacked datasets.
- namestr
Name of the dataset.
- models
See also
Attributes Summary
Exposure ratio between signal and background regions
Computed as alpha * n_off
A lazy FITS data descriptor.
Shape of the counts or background data (tuple)
A lazy FITS data descriptor.
Energy range maps defined by the mask_safe and mask_fit.
Energy range maps defined by the mask_fit only.
Energy range maps defined by the mask_safe only.
Largest energy range among all pixels, defined by mask_safe and mask_fit.
Model evaluators
Observed excess: counts-background
A lazy FITS data descriptor.
Map geometries
Combined fit and safe mask
A lazy FITS data descriptor.
Reduced mask fit
Reduced mask
A lazy FITS data descriptor.
Mask safe for edisp maps
Reduced mask safe
Mask safe for psf maps
Models set on the dataset (
Models
).A lazy FITS data descriptor.
Methods Summary
copy
([name])A deep copy.
create
(geom[, energy_axis_true, migra_axis, ...])Create a MapDataset object with zero filled maps.
cutout
(position, width[, mode, name])Cutout map dataset.
downsample
(factor[, axis_name, name])Downsample map dataset.
fake
(npred_background[, random_state])Simulate fake counts (on and off) for the current model and reduced IRFs.
from_dict
(data[, lazy, cache])Create from dicts and models list generated from YAML serialization.
from_geoms
(geom, geom_exposure[, geom_psf, ...])Create an empty
MapDatasetOnOff
object according to the specified geometriesfrom_hdulist
(hdulist[, name, format])Create map dataset from list of HDUs.
from_map_dataset
(dataset, acceptance, ...[, ...])Create on off dataset from a map dataset.
info_dict
([in_safe_data_range])Basic info dict with summary statistics
npred
()Total predicted source and background counts
Predicted background counts estimated from the marginalized likelihood estimate.
Predicted counts in the off region; mu_bkg/alpha
npred_signal
([model_name])Model predicted signal counts.
pad
()Pad the spatial dimensions of the dataset.
peek
([figsize])Quick-look summary plots.
plot_residuals
([ax_spatial, ax_spectral, ...])Plot spatial and spectral residuals in two panels.
plot_residuals_spatial
([ax, method, ...])Plot spatial residuals.
plot_residuals_spectral
([ax, method, region])Plot spectral residuals.
read
(filename[, name, lazy, cache, format])Read a dataset from file.
resample_energy_axis
(energy_axis[, name])Resample MapDatasetOnOff over reconstructed energy edges.
Reset data cache to free memory space
residuals
([method])Compute residuals map.
slice_by_energy
([energy_min, energy_max, name])Select and slice datasets in energy range
slice_by_idx
(slices[, name])Slice sub dataset.
stack
(other[, nan_to_num])Stack another dataset in place.
Likelihood per bin given the current model parameters
stat_sum
()Total likelihood given the current model parameters.
to_dict
()Convert to dict for YAML serialization.
Convert map dataset to list of HDUs.
to_image
([name])Create images by summing over the reconstructed energy axis.
to_map_dataset
([name])Convert a MapDatasetOnOff to MapDataset
to_masked
([name, nan_to_num])Return masked dataset
to_region_map_dataset
(region[, name])Integrate the map dataset in a given region.
to_spectrum_dataset
(on_region[, ...])Return a ~gammapy.datasets.SpectrumDatasetOnOff from on_region.
write
(filename[, overwrite])Write Dataset to file.
Attributes Documentation
- alpha#
Exposure ratio between signal and background regions
See WStat : Poisson data with background measurement
- Returns
- alpha
Map
Alpha map
- alpha
- background#
Computed as alpha * n_off
See WStat : Poisson data with background measurement
- Returns
- background
Map
Background map
- background
- background_model#
- counts#
A lazy FITS data descriptor.
- Parameters
- cachebool
Whether to cache the data.
- data_shape#
Shape of the counts or background data (tuple)
- edisp#
A lazy FITS data descriptor.
- Parameters
- cachebool
Whether to cache the data.
- energy_range#
Energy range maps defined by the mask_safe and mask_fit.
- energy_range_fit#
Energy range maps defined by the mask_fit only.
- energy_range_safe#
Energy range maps defined by the mask_safe only.
- energy_range_total#
Largest energy range among all pixels, defined by mask_safe and mask_fit.
- evaluators#
Model evaluators
- excess#
Observed excess: counts-background
- exposure#
A lazy FITS data descriptor.
- Parameters
- cachebool
Whether to cache the data.
- geoms#
Map geometries
- Returns
- geomsdict
Dict of map geometries involved in the dataset.
- mask#
Combined fit and safe mask
- mask_fit#
A lazy FITS data descriptor.
- Parameters
- cachebool
Whether to cache the data.
- mask_fit_image#
Reduced mask fit
- mask_image#
Reduced mask
- mask_safe#
A lazy FITS data descriptor.
- Parameters
- cachebool
Whether to cache the data.
- mask_safe_edisp#
Mask safe for edisp maps
- mask_safe_image#
Reduced mask safe
- mask_safe_psf#
Mask safe for psf maps
- name#
- psf#
A lazy FITS data descriptor.
- Parameters
- cachebool
Whether to cache the data.
- stat_type = 'wstat'#
- tag = 'MapDatasetOnOff'#
Methods Documentation
- copy(name=None)#
A deep copy.
- Parameters
- namestr
Name of the copied dataset
- Returns
- dataset
Dataset
Copied datasets.
- dataset
- classmethod create(geom, energy_axis_true=None, migra_axis=None, rad_axis=None, binsz_irf=None, reference_time='2000-01-01', name=None, meta_table=None, **kwargs)#
Create a MapDataset object with zero filled maps.
- Parameters
- geom
WcsGeom
Reference target geometry in reco energy, used for counts and background maps
- energy_axis_true
MapAxis
True energy axis used for IRF maps
- migra_axis
MapAxis
If set, this provides the migration axis for the energy dispersion map. If not set, an EDispKernelMap is produced instead. Default is None
- rad_axis
MapAxis
Rad axis for the psf map
- binsz_irffloat
IRF Map pixel size in degrees.
- reference_time
Time
the reference time to use in GTI definition
- namestr
Name of the returned dataset.
- meta_table
Table
Table listing information on observations used to create the dataset. One line per observation for stacked datasets.
- geom
- Returns
- empty_maps
MapDataset
A MapDataset containing zero filled maps
- empty_maps
Examples
>>> from gammapy.datasets import MapDataset >>> from gammapy.maps import WcsGeom, MapAxis
>>> energy_axis = MapAxis.from_energy_bounds(1.0, 10.0, 4, unit="TeV") >>> energy_axis_true = MapAxis.from_energy_bounds( 0.5, 20, 10, unit="TeV", name="energy_true" ) >>> geom = WcsGeom.create( skydir=(83.633, 22.014), binsz=0.02, width=(2, 2), frame="icrs", proj="CAR", axes=[energy_axis] ) >>> empty = MapDataset.create(geom=geom, energy_axis_true=energy_axis_true, name="empty")
- cutout(position, width, mode='trim', name=None)[source]#
Cutout map dataset.
- Parameters
- position
SkyCoord
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
.- namestr
Name of the new dataset.
- position
- Returns
- cutout
MapDatasetOnOff
Cutout map dataset.
- cutout
- downsample(factor, axis_name=None, name=None)[source]#
Downsample map dataset.
The PSFMap and EDispKernelMap are not downsampled, except if a corresponding axis is given.
- Parameters
- factorint
Downsampling factor.
- axis_namestr
Which non-spatial axis to downsample. By default only spatial axes are downsampled.
- namestr
Name of the downsampled dataset.
- Returns
- dataset
MapDatasetOnOff
Downsampled map dataset.
- dataset
- fake(npred_background, random_state='random-seed')[source]#
Simulate fake counts (on and off) for the current model and reduced IRFs.
This method overwrites the counts defined on the dataset object.
- Parameters
- npred_background
Map
Expected number of background counts in the on region
- random_state{int, ‘random-seed’, ‘global-rng’,
RandomState
} Defines random number generator initialisation. Passed to
get_random_state
.
- npred_background
- classmethod from_dict(data, lazy=False, cache=True)#
Create from dicts and models list generated from YAML serialization.
- classmethod from_geoms(geom, geom_exposure, geom_psf=None, geom_edisp=None, reference_time='2000-01-01', name=None, **kwargs)[source]#
Create an empty
MapDatasetOnOff
object according to the specified geometries- Parameters
- geom
gammapy.maps.WcsGeom
geometry for the counts, counts_off, acceptance and acceptance_off maps
- geom_exposure
gammapy.maps.WcsGeom
geometry for the exposure map
- geom_psf
gammapy.maps.WcsGeom
geometry for the psf map
- geom_edisp
gammapy.maps.WcsGeom
geometry for the energy dispersion kernel map. If geom_edisp has a migra axis, this will create an EDispMap instead.
- reference_time
Time
the reference time to use in GTI definition
- namestr
Name of the returned dataset.
- geom
- Returns
- empty_maps
MapDatasetOnOff
A MapDatasetOnOff containing zero filled maps
- empty_maps
- classmethod from_hdulist(hdulist, name=None, format='gadf')[source]#
Create map dataset from list of HDUs.
- Parameters
- hdulist
HDUList
List of HDUs.
- namestr
Name of the new dataset.
- format{“gadf”}
Format the hdulist is given in.
- hdulist
- Returns
- dataset
MapDatasetOnOff
Map dataset.
- dataset
- classmethod from_map_dataset(dataset, acceptance, acceptance_off, counts_off=None, name=None)[source]#
Create on off dataset from a map dataset.
- Parameters
- dataset
MapDataset
Spectrum dataset defining counts, edisp, aeff, livetime etc.
- acceptance
Map
Relative background efficiency in the on region.
- acceptance_off
Map
Relative background efficiency in the off region.
- counts_off
Map
Off counts map . If the dataset provides a background model, and no off counts are defined. The off counts are deferred from counts_off / alpha.
- namestr
Name of the returned dataset.
- dataset
- Returns
- dataset
MapDatasetOnOff
Map dataset on off.
- dataset
- info_dict(in_safe_data_range=True)[source]#
Basic info dict with summary statistics
If a region is passed, then a spectrum dataset is extracted, and the corresponding info returned.
- Parameters
- in_safe_data_rangebool
Whether to sum only in the safe energy range
- Returns
- info_dictdict
Dictionary with summary info.
- npred()#
Total predicted source and background counts
- Returns
- npred
Map
Total predicted counts
- npred
- npred_background()[source]#
Predicted background counts estimated from the marginalized likelihood estimate.
See WStat : Poisson data with background measurement
- Returns
- npred_background
Map
Predicted background counts
- npred_background
- npred_off()[source]#
Predicted counts in the off region; mu_bkg/alpha
See WStat : Poisson data with background measurement
- Returns
- npred_off
Map
Predicted off counts
- npred_off
- npred_signal(model_name=None)#
Model predicted signal counts.
If a model name is passed, predicted counts from that component are returned. Else, the total signal counts are returned.
- Parameters
- model_name: str
Name of SkyModel for which to compute the npred for. If none, the sum of all components (minus the background model) is returned
- Returns
- npred_sig:
gammapy.maps.Map
Map of the predicted signal counts
- npred_sig:
- pad()[source]#
Pad the spatial dimensions of the dataset.
The padding only applies to counts, masks, background and exposure.
Counts, background and masks are padded with zeros, exposure is padded with edge value.
- Parameters
- pad_width{sequence, array_like, int}
Number of pixels padded to the edges of each axis.
- namestr
Name of the padded dataset.
- Returns
- dataset
MapDataset
Padded map dataset.
- dataset
- peek(figsize=(12, 10))#
Quick-look summary plots.
- Parameters
- figsizetuple
Size of the figure.
- plot_residuals(ax_spatial=None, ax_spectral=None, kwargs_spatial=None, kwargs_spectral=None)#
Plot spatial and spectral residuals in two panels.
Calls
plot_residuals_spatial
andplot_residuals_spectral
. The spectral residuals are extracted from the provided region, and the normalization used for its computation can be controlled using the method parameter. The region outline is overlaid on the residuals map. If no region is passed, the residuals are computed for the entire map- Parameters
- ax_spatial
WCSAxes
Axes to plot spatial residuals on.
- ax_spectral
Axes
Axes to plot spectral residuals on.
- kwargs_spatialdict
Keyword arguments passed to
plot_residuals_spatial
.- kwargs_spectraldict
Keyword arguments passed to
plot_residuals_spectral
. The region should be passed as a dictionary key
- ax_spatial
- Returns
Examples
>>> from regions import CircleSkyRegion >>> from astropy.coordinates import SkyCoord >>> import astropy.units as u >>> from gammapy.datasets import MapDataset >>> dataset = MapDataset.read("$GAMMAPY_DATA/cta-1dc-gc/cta-1dc-gc.fits.gz") >>> reg = CircleSkyRegion(SkyCoord(0,0, unit="deg", frame="galactic"), radius=1.0 * u.deg) >>> kwargs_spatial = {"cmap": "RdBu_r", "vmin":-5, "vmax":5, "add_cbar": True} >>> kwargs_spectral = {"region":reg, "markerfacecolor": "blue", "markersize": 8, "marker": "s"} # noqa: E501 >>> dataset.plot_residuals(kwargs_spatial=kwargs_spatial, kwargs_spectral=kwargs_spectral)
- plot_residuals_spatial(ax=None, method='diff', smooth_kernel='gauss', smooth_radius='0.1 deg', **kwargs)#
Plot spatial residuals.
The normalization used for the residuals computation can be controlled using the method parameter.
- Parameters
- ax
WCSAxes
Axes to plot on.
- method{“diff”, “diff/model”, “diff/sqrt(model)”}
Normalization used to compute the residuals, see
MapDataset.residuals
.- smooth_kernel{“gauss”, “box”}
Kernel shape.
- smooth_radius: `~astropy.units.Quantity`, str or float
Smoothing width given as quantity or float. If a float is given, it is interpreted as smoothing width in pixels.
- **kwargsdict
Keyword arguments passed to
imshow
.
- ax
- Returns
- ax
WCSAxes
WCSAxes object.
- ax
Examples
>>> from gammapy.datasets import MapDataset >>> dataset = MapDataset.read("$GAMMAPY_DATA/cta-1dc-gc/cta-1dc-gc.fits.gz") >>> kwargs = {"cmap": "RdBu_r", "vmin":-5, "vmax":5, "add_cbar": True} >>> dataset.plot_residuals_spatial(method="diff/sqrt(model)", **kwargs)
- plot_residuals_spectral(ax=None, method='diff', region=None, **kwargs)#
Plot spectral residuals.
The residuals are extracted from the provided region, and the normalization used for its computation can be controlled using the method parameter.
The error bars are computed using the uncertainty on the excess with a symmetric assumption.
- Parameters
- ax
Axes
Axes to plot on.
- method{“diff”, “diff/sqrt(model)”}
Normalization used to compute the residuals, see
SpectrumDataset.residuals
.- region: `~regions.SkyRegion` (required)
Target sky region.
- **kwargsdict
Keyword arguments passed to
errorbar
.
- ax
- Returns
- ax
Axes
Axes object.
- ax
Examples
>>> from gammapy.datasets import MapDataset >>> dataset = MapDataset.read("$GAMMAPY_DATA/cta-1dc-gc/cta-1dc-gc.fits.gz") >>> kwargs = {"markerfacecolor": "blue", "markersize":8, "marker":'s'} >>> dataset.plot_residuals_spectral(method="diff/sqrt(model)", **kwargs)
- classmethod read(filename, name=None, lazy=False, cache=True, format='gadf')#
Read a dataset from file.
- Parameters
- filenamestr
Filename to read from.
- namestr
Name of the new dataset.
- lazybool
Whether to lazy load data into memory
- cachebool
Whether to cache the data after loading.
- format{“gadf”}
Format of the dataset file.
- Returns
- dataset
MapDataset
Map dataset.
- dataset
- resample_energy_axis(energy_axis, name=None)[source]#
Resample MapDatasetOnOff over reconstructed energy edges.
Counts are summed taking into account safe mask.
- Parameters
- energy_axis
MapAxis
New reco energy axis.
- name: str
Name of the new dataset.
- energy_axis
- Returns
- dataset:
SpectrumDataset
Resampled spectrum dataset .
- dataset:
- reset_data_cache()#
Reset data cache to free memory space
- residuals(method='diff', **kwargs)#
Compute residuals map.
- Parameters
- method: {“diff”, “diff/model”, “diff/sqrt(model)”}
- Method used to compute the residuals. Available options are:
“diff” (default): data - model
“diff/model”: (data - model) / model
“diff/sqrt(model)”: (data - model) / sqrt(model)
- **kwargsdict
Keyword arguments forwarded to
Map.smooth()
- Returns
- residuals
gammapy.maps.Map
Residual map.
- residuals
- slice_by_energy(energy_min=None, energy_max=None, name=None)#
Select and slice datasets in energy range
- Parameters
- energy_min, energy_max
Quantity
Energy bounds to compute the flux point for.
- namestr
Name of the sliced dataset.
- energy_min, energy_max
- Returns
- dataset
MapDataset
Sliced Dataset
- dataset
Examples
>>> from gammapy.datasets import MapDataset >>> dataset = MapDataset.read("$GAMMAPY_DATA/cta-1dc-gc/cta-1dc-gc.fits.gz") >>> sliced = dataset.slice_by_energy(energy_min="1 TeV", energy_max="5 TeV") >>> sliced.data_shape (3, 240, 320)
- slice_by_idx(slices, name=None)[source]#
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.- namestr
Name of the sliced dataset.
- Returns
- map_out
Map
Sliced map object.
- map_out
- stack(other, nan_to_num=True)[source]#
Stack another dataset in place.
The
acceptance
of the stacked dataset is normalized to 1, and the stackedacceptance_off
is scaled so that:\[\alpha_\text{stacked} = \frac{1}{a_\text{off}} = \frac{\alpha_1\text{OFF}_1 + \alpha_2\text{OFF}_2}{\text{OFF}_1 + OFF_2}\]- Parameters
- other
MapDatasetOnOff
Other dataset
- nan_to_num: bool
Non-finite values are replaced by zero if True (default).
- other
- to_dict()#
Convert to dict for YAML serialization.
- to_hdulist()[source]#
Convert map dataset to list of HDUs.
- Returns
- hdulist
HDUList
Map dataset list of HDUs.
- hdulist
- to_image(name=None)#
Create images by summing over the reconstructed energy axis.
- Parameters
- namestr
Name of the new dataset.
- Returns
- dataset
MapDataset
orSpectrumDataset
Dataset integrated over non-spatial axes.
- dataset
- to_map_dataset(name=None)[source]#
Convert a MapDatasetOnOff to MapDataset
The background model template is taken as alpha * counts_off
- Parameters
- name: str
Name of the new dataset
- Returns
- dataset:
MapDataset
Map dataset with cash statistics
- dataset:
- to_masked(name=None, nan_to_num=True)#
Return masked dataset
- Parameters
- namestr
Name of the masked dataset.
- nan_to_num: bool
Non-finite values are replaced by zero if True (default).
- Returns
- dataset
MapDataset
orSpectrumDataset
Masked dataset
- dataset
- to_region_map_dataset(region, name=None)#
Integrate the map dataset in a given region.
Counts and background of the dataset are integrated in the given region, taking the safe mask into accounts. The exposure is averaged in the region again taking the safe mask into account. The PSF and energy dispersion kernel are taken at the center of the region.
- Parameters
- region
SkyRegion
Region from which to extract the spectrum
- namestr
Name of the new dataset.
- region
- Returns
- dataset
MapDataset
the resulting reduced dataset
- dataset
- to_spectrum_dataset(on_region, containment_correction=False, name=None)[source]#
Return a ~gammapy.datasets.SpectrumDatasetOnOff from on_region.
Counts and OFF counts are summed in the on_region.
Acceptance is the average of all acceptances while acceptance OFF is taken such that number of excess is preserved in the on_region.
Effective area is taken from the average exposure.
The energy dispersion kernel is obtained at the on_region center. Only regions with centers are supported.
The models are not exported to the ~gammapy.dataset.SpectrumDatasetOnOff. It must be set after the dataset extraction.
- Parameters
- on_region
SkyRegion
the input ON region on which to extract the spectrum
- containment_correctionbool
Apply containment correction for point sources and circular on regions
- namestr
Name of the new dataset.
- on_region
- Returns
- dataset
SpectrumDatasetOnOff
the resulting reduced dataset
- dataset
- write(filename, overwrite=False)#
Write Dataset to file.
A MapDataset is serialised using the GADF format with a WCS geometry. A SpectrumDataset uses the same format, with a RegionGeom.
- Parameters
- filenamestr
Filename to write to.
- overwritebool
Overwrite file if it exists.