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, meta=None)[source]#

Bases: MapDataset

Map dataset for on-off likelihood fitting.

It bundles together the binned on and off counts, the binned IRFs as well as the on and off acceptances.

A safe mask and a fit mask can be added to exclude bins during the analysis.

It uses the Wstat statistic (see wstat), therefore no background model is needed.

For more information see Datasets (DL4).

Parameters:
modelsModels

Source sky models.

countsWcsNDMap

Counts cube.

counts_offWcsNDMap

Ring-convolved counts cube.

acceptanceWcsNDMap

Acceptance from the IRFs.

acceptance_offWcsNDMap

Acceptance off.

exposureWcsNDMap

Exposure cube.

mask_fitWcsNDMap

Mask to apply to the likelihood for fitting.

psfPSFKernel

PSF kernel.

edispEDispKernel

Energy dispersion.

mask_safeWcsNDMap

Mask defining the safe data range.

gtiGTI

GTI of the observation or union of GTI if it is a stacked observation.

meta_tableTable

Table listing information on observations used to create the dataset. One line per observation for stacked datasets.

namestr

Name of the dataset.

metaMapDatasetMetaData

Associated meta data container

See also

MapDataset, SpectrumDataset, FluxPointsDataset.

Attributes Summary

alpha

Exposure ratio between signal and background regions.

background

Computed as alpha * n_off.

background_model

counts

A lazy FITS data descriptor.

data_shape

Shape of the counts or background data (tuple).

edisp

A lazy FITS data descriptor.

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.

geoms

Map geometries.

gti

mask

Combined fit and safe mask.

mask_fit

A lazy FITS data descriptor.

mask_fit_image

Reduced fit mask.

mask_image

Reduced mask.

mask_safe

A lazy FITS data descriptor.

mask_safe_edisp

Safe mask for edisp maps.

mask_safe_image

Reduced safe mask.

mask_safe_psf

Safe mask for PSF maps.

meta

meta_table

models

Models set on the dataset (Models).

name

psf

A lazy FITS data descriptor.

stat_type

tag

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

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

npred_background()

Predicted background counts estimated from the marginalized likelihood estimate.

npred_off()

Predicted counts in the off region; mu_bkg/alpha.

npred_signal([model_names, stack])

Model predicted signal counts.

pad()

Not implemented for MapDatasetOnOff.

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()

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.

stat_array()

Statistic function value per bin given the current model parameters.

stat_sum()

Total statistic function value given the current model parameters.

to_dict()

Convert to dict for YAML serialization.

to_hdulist()

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 a 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, checksum])

Write Dataset to file.

Attributes Documentation

alpha#

Exposure ratio between signal and background regions.

See WStat : Poisson data with background measurement.

Returns:
alphaMap

Alpha map.

background#

Computed as alpha * n_off.

See WStat : Poisson data with background measurement.

Returns:
backgroundMap

Background map.

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

Dictionary of map geometries involved in the dataset.

gti = None#
mask#

Combined fit and safe mask.

mask_fit#

A lazy FITS data descriptor.

Parameters:
cachebool

Whether to cache the data.

mask_fit_image#

Reduced fit mask.

mask_image#

Reduced mask.

mask_safe#

A lazy FITS data descriptor.

Parameters:
cachebool

Whether to cache the data.

mask_safe_edisp#

Safe mask for edisp maps.

mask_safe_image#

Reduced safe mask.

mask_safe_psf#

Safe mask for PSF maps.

meta#
meta_table = None#
models#

Models set on the dataset (Models).

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

Name of the copied dataset. Default is None.

Returns:
datasetDataset

Copied datasets.

classmethod create(geom, energy_axis_true=None, migra_axis=None, rad_axis=None, binsz_irf=<Quantity 0.2 deg>, reference_time='2000-01-01', name=None, meta_table=None, reco_psf=False, **kwargs)#

Create a MapDataset object with zero filled maps.

Parameters:
geomWcsGeom

Reference target geometry in reco energy, used for counts and background maps.

energy_axis_trueMapAxis, optional

True energy axis used for IRF maps. Default is None.

migra_axisMapAxis, optional

If set, this provides the migration axis for the energy dispersion map. If not set, an EDispKernelMap is produced instead. Default is None.

rad_axisMapAxis, optional

Rad axis for the PSF map. Default is None.

binsz_irffloat

IRF Map pixel size in degrees. Default is BINSZ_IRF_DEFAULT.

reference_timeTime

The reference time to use in GTI definition. Default is “2000-01-01”.

namestr, optional

Name of the returned dataset. Default is None.

meta_tableTable, optional

Table listing information on observations used to create the dataset. One line per observation for stacked datasets. Default is None.

reco_psfbool

Use reconstructed energy for the PSF geometry. Default is False.

Returns:
empty_mapsMapDataset

A MapDataset containing zero filled 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:
positionSkyCoord

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. Default is “trim”.

namestr, optional

Name of the new dataset. Default is None.

Returns:
cutoutMapDatasetOnOff

Cutout map dataset.

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

Which non-spatial axis to downsample. By default, only spatial axes are downsampled. Default is None.

namestr, optional

Name of the downsampled dataset. Default is None.

Returns:
datasetMapDatasetOnOff

Downsampled map 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_backgroundMap

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. Default is “random-seed”.

classmethod from_dict(data, lazy=False, cache=True)#

Create from dicts and models list generated from YAML serialization.

classmethod from_geoms(geom, geom_exposure=None, 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:
geomgammapy.maps.WcsGeom

Geometry for the counts, counts_off, acceptance and acceptance_off maps.

geom_exposuregammapy.maps.WcsGeom, optional

Geometry for the exposure map. Default is None.

geom_psfgammapy.maps.WcsGeom, optional

Geometry for the PSF map. Default is None.

geom_edispgammapy.maps.WcsGeom, optional

Geometry for the energy dispersion kernel map. If geom_edisp has a migra axis, this will create an EDispMap instead. Default is None.

reference_timeTime

The reference time to use in GTI definition. Default is “2000-01-01”.

namestr, optional

Name of the returned dataset. Default is None.

**kwargsdict, optional

Keyword arguments to be passed.

Returns:
empty_mapsMapDatasetOnOff

A MapDatasetOnOff containing zero filled maps.

classmethod from_hdulist(hdulist, name=None, format='gadf')[source]#

Create map dataset from list of HDUs.

Parameters:
hdulistHDUList

List of HDUs.

namestr, optional

Name of the new dataset. Default is None.

format{“gadf”}

Format the hdulist is given in. Default is “gadf”.

Returns:
datasetMapDatasetOnOff

Map dataset.

classmethod from_map_dataset(dataset, acceptance, acceptance_off, counts_off=None, name=None)[source]#

Create on off dataset from a map dataset.

Parameters:
datasetMapDataset

Spectrum dataset defining counts, edisp, aeff, livetime etc.

acceptanceMap

Relative background efficiency in the on region.

acceptance_offMap

Relative background efficiency in the off region.

counts_offMap, optional

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. Default is None.

namestr, optional

Name of the returned dataset. Default is None.

Returns:
datasetMapDatasetOnOff

Map dataset on off.

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. Default is True.

Returns:
info_dictdict

Dictionary with summary info.

npred()#

Total predicted source and background counts.

Returns:
npredMap

Total predicted counts.

npred_background()[source]#

Predicted background counts estimated from the marginalized likelihood estimate.

See WStat : Poisson data with background measurement.

Returns:
npred_backgroundMap

Predicted background counts.

npred_off()[source]#

Predicted counts in the off region; mu_bkg/alpha.

See WStat : Poisson data with background measurement.

Returns:
npred_offMap

Predicted off counts.

npred_signal(model_names=None, stack=True)#

Model predicted signal counts.

If a list of model name is passed, predicted counts from these components are returned. If stack is set to True, a map of the sum of all the predicted counts is returned. If stack is set to False, a map with an additional axis representing the models is returned.

Parameters:
model_nameslist of str

List of name of SkyModel for which to compute the npred. If none, all the SkyModel predicted counts are computed.

stackbool

Whether to stack the npred maps upon each other.

Returns:
npred_siggammapy.maps.Map

Map of the predicted signal counts.

pad()[source]#

Not implemented for MapDatasetOnOff.

peek(figsize=(12, 8))#

Quick-look summary plots.

Parameters:
figsizetuple

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

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 and plot_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_spatialWCSAxes, optional

Axes to plot spatial residuals on. Default is None.

ax_spectralAxes, optional

Axes to plot spectral residuals on. Default is None.

kwargs_spatialdict, optional

Keyword arguments passed to plot_residuals_spatial. Default is None.

kwargs_spectraldict, optional

Keyword arguments passed to plot_residuals_spectral. The region should be passed as a dictionary key. Default is None.

Returns:
ax_spatial, ax_spectralWCSAxes, Axes

Spatial and spectral residuals plots.

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"}
>>> 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:
axWCSAxes, optional

Axes to plot on. Default is None.

method{“diff”, “diff/model”, “diff/sqrt(model)”}

Normalization used to compute the residuals, see MapDataset.residuals. Default is “diff”.

smooth_kernel{“gauss”, “box”}

Kernel shape. Default is “gauss”.

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. Default is “0.1 deg”.

**kwargsdict, optional

Keyword arguments passed to imshow.

Returns:
axWCSAxes

WCSAxes object.

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

Axes to plot on. Default is None.

method{“diff”, “diff/sqrt(model)”}

Normalization used to compute the residuals, see SpectrumDataset.residuals. Default is “diff”.

regionSkyRegion (required)

Target sky region. Default is None.

**kwargsdict, optional

Keyword arguments passed to errorbar.

Returns:
axAxes

Axes object.

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', checksum=False)#

Read a dataset from file.

Parameters:
filenamestr

Filename to read from.

namestr, optional

Name of the new dataset. Default is None.

lazybool

Whether to lazy load data into memory. Default is False.

cachebool

Whether to cache the data after loading. Default is True.

format{“gadf”}

Format of the dataset file. Default is “gadf”.

checksumbool

If True checks both DATASUM and CHECKSUM cards in the file headers. Default is False.

Returns:
datasetMapDataset

Map 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_axisMapAxis

New reco energy axis.

namestr, optional

Name of the new dataset. Default is None.

Returns:
datasetSpectrumDataset

Resampled spectrum 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).

Default is “diff”.

**kwargsdict, optional

Keyword arguments forwarded to Map.smooth().

Returns:
residualsgammapy.maps.Map

Residual map.

slice_by_energy(energy_min=None, energy_max=None, name=None)#

Select and slice datasets in energy range.

Parameters:
energy_min, energy_maxQuantity, optional

Energy bounds to compute the flux point for. Default is None.

namestr, optional

Name of the sliced dataset. Default is None.

Returns:
datasetMapDataset

Sliced 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

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

Name of the sliced dataset. Default is None.

Returns:
map_outMap

Sliced map object.

stack(other, nan_to_num=True)[source]#

Stack another dataset in place.

The acceptance of the stacked dataset is normalized to 1, and the stacked acceptance_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:
otherMapDatasetOnOff

Other dataset.

nan_to_numbool

Non-finite values are replaced by zero if True. Default is True.

stat_array()[source]#

Statistic function value per bin given the current model parameters.

stat_sum()[source]#

Total statistic function value given the current model parameters.

If the off counts are passed as None and the elements of the safe mask are False, zero will be returned. Otherwise, the stat sum will be calculated and returned.

to_dict()#

Convert to dict for YAML serialization.

to_hdulist()[source]#

Convert map dataset to list of HDUs.

Returns:
hdulistHDUList

Map dataset list of HDUs.

to_image(name=None)#

Create images by summing over the reconstructed energy axis.

Parameters:
namestr, optional

Name of the new dataset. Default is None.

Returns:
datasetMapDataset or SpectrumDataset

Dataset integrated over non-spatial axes.

to_map_dataset(name=None)[source]#

Convert a MapDatasetOnOff to a MapDataset.

The background model template is taken as alpha * counts_off.

Parameters:
namestr, optional

Name of the new dataset. Default is None.

Returns:
datasetMapDataset

Map dataset with cash statistics.

to_masked(name=None, nan_to_num=True)#

Return masked dataset.

Parameters:
namestr, optional

Name of the masked dataset. Default is None.

nan_to_numbool

Non-finite values are replaced by zero if True. Default is True.

Returns:
datasetMapDataset or SpectrumDataset

Masked 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 account. 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:
regionSkyRegion

Region from which to extract the spectrum.

namestr, optional

Name of the new dataset. Default is None.

Returns:
datasetMapDataset

The resulting reduced 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_regionSkyRegion

The input ON region on which to extract the spectrum.

containment_correctionbool

Apply containment correction for point sources and circular on regions. Default is False.

namestr, optional

Name of the new dataset. Default is None.

Returns:
datasetSpectrumDatasetOnOff

The resulting reduced dataset.

write(filename, overwrite=False, checksum=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, optional

Overwrite existing file. Default is False.

checksumbool

When True adds both DATASUM and CHECKSUM cards to the headers written to the file. Default is False.