SpectrumDatasetOnOff#
- class gammapy.datasets.SpectrumDatasetOnOff(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:
gammapy.datasets.spectrum.PlotMixin
,gammapy.datasets.map.MapDatasetOnOff
Spectrum dataset for 1D 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 fit mask can be added to exclude bins during the analysis.
It uses the Wstat statistic (see
wstat
).For more information see Datasets (DL4).
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 fit mask.
Reduced mask.
A lazy FITS data descriptor.
Safe mask for edisp maps.
Reduced safe mask.
Safe mask 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
(*args, **kwargs)Not supported for
SpectrumDatasetOnOff
.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, **kwargs)Create spectrum dataset from dict.
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.
from_spectrum_dataset
(**kwargs)Create a SpectrumDatasetOnOff from a
SpectrumDataset
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_names, stack])Model predicted signal counts.
pad
()Not implemented for MapDatasetOnOff.
peek
([figsize])Quick-look summary plots.
plot_counts
([ax, kwargs_counts, ...])Plot counts and background.
plot_excess
([ax, kwargs_excess, ...])Plot excess and predicted signal.
plot_fit
([ax_spectrum, ax_residuals, ...])Plot spectrum and residuals in two panels.
plot_masks
([ax, kwargs_fit, kwargs_safe])Plot safe mask and fit mask.
plot_residuals
([ax_spatial, ax_spectral, ...])Plot spatial and spectral residuals in two panels.
plot_residuals_spatial
(*args, **kwargs)Not supported for
SpectrumDatasetOnOff
.plot_residuals_spectral
([ax, method, region])Plot spectral residuals.
read
(filename[, format, checksum])Read 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.
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.
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
([name])Convert a SpectrumDatasetOnOff to a SpectrumDataset.
write
(filename[, overwrite, format, checksum])Write spectrum dataset on off 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
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#
- name#
- psf#
A lazy FITS data descriptor.
- Parameters
- cachebool
Whether to cache the data.
- stat_type = 'wstat'#
- tag = 'SpectrumDatasetOnOff'#
Methods Documentation
- copy(name=None)#
A deep copy.
- Parameters
- namestr, optional
Name of the copied dataset. Default is None.
- Returns
- dataset
Dataset
Copied datasets.
- dataset
- 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
- geom
WcsGeom
Reference target geometry in reco energy, used for counts and background maps.
- energy_axis_true
MapAxis
, optional True energy axis used for IRF maps. Default is None.
- migra_axis
MapAxis
, 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_axis
MapAxis
, optional Rad axis for the PSF map. Default is None.
- binsz_irffloat
IRF Map pixel size in degrees. Default is BINSZ_IRF_DEFAULT.
- reference_time
Time
The reference time to use in GTI definition. Default is “2000-01-01”.
- namestr, optional
Name of the returned dataset. Default is None.
- meta_table
Table
, 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.
- 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(*args, **kwargs)[source]#
Not supported for
SpectrumDatasetOnOff
.
- downsample(factor, axis_name=None, name=None)#
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
- dataset
MapDatasetOnOff
Downsampled map dataset.
- dataset
- fake(npred_background, random_state='random-seed')#
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
. Default is “random-seed”.
- npred_background
- classmethod from_dict(data, **kwargs)[source]#
Create spectrum dataset from dict. Reads file from the disk as specified in the dict.
- Parameters
- datadict
Dictionary containing data to create dataset from.
- Returns
- dataset
SpectrumDatasetOnOff
Spectrum dataset on off.
- dataset
- classmethod from_geoms(geom, geom_exposure=None, geom_psf=None, geom_edisp=None, reference_time='2000-01-01', name=None, **kwargs)#
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
, optional Geometry for the exposure map. Default is None.
- geom_psf
gammapy.maps.WcsGeom
, optional Geometry for the PSF map. Default is None.
- geom_edisp
gammapy.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_time
Time
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.
- geom
- Returns
- empty_maps
MapDatasetOnOff
A MapDatasetOnOff containing zero filled maps.
- empty_maps
- classmethod from_hdulist(hdulist, name=None, format='gadf')#
Create map dataset from list of HDUs.
- Parameters
- hdulist
HDUList
List of HDUs.
- namestr, optional
Name of the new dataset. Default is None.
- format{“gadf”}
Format the hdulist is given in. Default is “gadf”.
- hdulist
- Returns
- dataset
MapDatasetOnOff
Map dataset.
- dataset
- classmethod from_map_dataset(dataset, acceptance, acceptance_off, counts_off=None, name=None)#
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
, 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.
- dataset
- Returns
- dataset
MapDatasetOnOff
Map dataset on off.
- dataset
- classmethod from_spectrum_dataset(**kwargs)[source]#
Create a SpectrumDatasetOnOff from a
SpectrumDataset
dataset.- Parameters
- dataset
SpectrumDataset
Spectrum dataset defining counts, edisp, exposure etc.
- acceptance
array
or float Relative background efficiency in the on region.
- acceptance_off
array
or float Relative background efficiency in the off region.
- counts_off
RegionNDMap
Off counts spectrum. If the dataset provides a background model, and no off counts are defined. The off counts are deferred from counts_off / alpha.
- dataset
- Returns
- dataset
SpectrumDatasetOnOff
Spectrum dataset on off.
- dataset
- info_dict(in_safe_data_range=True)#
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
- npred
Map
Total predicted counts.
- npred
- npred_background()#
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()#
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_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_sig
gammapy.maps.Map
Map of the predicted signal counts.
- npred_sig
- pad()#
Not implemented for MapDatasetOnOff.
- peek(figsize=(16, 4))#
Quick-look summary plots.
- Parameters
- figsizetuple
Size of the figure. Default is (16, 4).
- plot_counts(ax=None, kwargs_counts=None, kwargs_background=None, **kwargs)#
Plot counts and background.
- Parameters
- ax
Axes
, optional Axes to plot on. Default is None.
- kwargs_countsdict, optional
Keyword arguments passed to
hist
for the counts. Default is None.- kwargs_backgrounddict, optional
Keyword arguments passed to
hist
for the background. Default is None.- **kwargsdict, optional
Keyword arguments passed to both
hist
.
- ax
- Returns
- ax
Axes
Axes object.
- ax
- plot_excess(ax=None, kwargs_excess=None, kwargs_npred_signal=None, **kwargs)#
Plot excess and predicted signal.
The error bars are computed with a symmetric assumption on the excess.
- Parameters
- ax
Axes
, optional Axes to plot on. Default is None.
- kwargs_excessdict, optional
Keyword arguments passed to
errorbar
for the excess. Default is None.- kwargs_npred_signaldict, optional
Keyword arguments passed to
hist
for the predicted signal. Default is None.- **kwargsdict, optional
Keyword arguments passed to both plot methods.
- ax
- Returns
- ax
Axes
Axes object.
- ax
Examples
>>> #Creating a spectral dataset >>> from gammapy.datasets import SpectrumDatasetOnOff >>> from gammapy.modeling.models import PowerLawSpectralModel, SkyModel >>> filename = "$GAMMAPY_DATA/joint-crab/spectra/hess/pha_obs23523.fits" >>> dataset = SpectrumDatasetOnOff.read(filename) >>> p = PowerLawSpectralModel() >>> dataset.models = SkyModel(spectral_model=p) >>> #Plot the excess in blue and the npred in black dotted lines >>> kwargs_excess = {"color": "blue", "markersize":8, "marker":'s', } >>> kwargs_npred_signal = {"color": "black", "ls":"--"} >>> dataset.plot_excess(kwargs_excess=kwargs_excess, kwargs_npred_signal=kwargs_npred_signal)
- plot_fit(ax_spectrum=None, ax_residuals=None, kwargs_spectrum=None, kwargs_residuals=None)#
Plot spectrum and residuals in two panels.
Calls
plot_excess
andplot_residuals_spectral
.- Parameters
- ax_spectrum
Axes
, optional Axes to plot spectrum on. Default is None.
- ax_residuals
Axes
, optional Axes to plot residuals on. Default is None.
- kwargs_spectrumdict, optional
Keyword arguments passed to
plot_excess
. Default is None.- kwargs_residualsdict, optional
Keyword arguments passed to
plot_residuals_spectral
. Default is None.
- ax_spectrum
- Returns
- ax_spectrum, ax_residuals
Axes
Spectrum and residuals plots.
- ax_spectrum, ax_residuals
Examples
>>> #Creating a spectral dataset >>> from gammapy.datasets import SpectrumDatasetOnOff >>> from gammapy.modeling.models import PowerLawSpectralModel, SkyModel >>> filename = "$GAMMAPY_DATA/joint-crab/spectra/hess/pha_obs23523.fits" >>> dataset = SpectrumDatasetOnOff.read(filename) >>> p = PowerLawSpectralModel() >>> dataset.models = SkyModel(spectral_model=p) >>> # optional configurations >>> kwargs_excess = {"color": "blue", "markersize":8, "marker":'s', } >>> kwargs_npred_signal = {"color": "black", "ls":"--"} >>> kwargs_spectrum = {"kwargs_excess":kwargs_excess, "kwargs_npred_signal":kwargs_npred_signal} >>> kwargs_residuals = {"color": "black", "markersize":4, "marker":'s', } >>> dataset.plot_fit(kwargs_residuals=kwargs_residuals, kwargs_spectrum=kwargs_spectrum)
- plot_masks(ax=None, kwargs_fit=None, kwargs_safe=None)#
Plot safe mask and fit mask.
- Parameters
- ax
Axes
, optional Axes to plot on. Default is None.
- kwargs_fitdict, optional
Keyword arguments passed to
plot_mask()
for mask fit. Default is None.- kwargs_safedict, optional
Keyword arguments passed to
plot_mask()
for mask safe. Default is None.
- ax
- Returns
- ax
Axes
Axes object.
- ax
Examples
>>> # Reading a spectral dataset >>> from gammapy.datasets import SpectrumDatasetOnOff >>> filename = "$GAMMAPY_DATA/joint-crab/spectra/hess/pha_obs23523.fits" >>> dataset = SpectrumDatasetOnOff.read(filename) >>> dataset.mask_fit = dataset.mask_safe.copy() >>> dataset.mask_fit.data[40:46] = False # setting dummy mask_fit for illustration >>> # Plot the masks on top of the counts histogram >>> kwargs_safe = {"color":"green", "alpha":0.2} #optinonal arguments to configure >>> kwargs_fit = {"color":"pink", "alpha":0.2} >>> ax=dataset.plot_counts() >>> dataset.plot_masks(ax=ax, kwargs_fit=kwargs_fit, kwargs_safe=kwargs_safe)
- 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
, optional Axes to plot spatial residuals on. Default is None.
- ax_spectral
Axes
, 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.
- 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"} >>> dataset.plot_residuals(kwargs_spatial=kwargs_spatial, kwargs_spectral=kwargs_spectral)
- plot_residuals_spatial(*args, **kwargs)[source]#
Not supported for
SpectrumDatasetOnOff
.
- 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
, 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”.- region
SkyRegion
(required) Target sky region. Default is None.
- **kwargsdict, optional
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, format='ogip', checksum=False, **kwargs)[source]#
Read from file.
For OGIP formats, filename is the name of a PHA file. The BKG, ARF, and RMF file names must be set in the PHA header and the files must be present in the same folder. For details, see
OGIPDatasetReader.read
.For the GADF format, a MapDataset serialisation is used.
- Parameters
- filename
Path
or str OGIP PHA file to read.
- format{“ogip”, “ogip-sherpa”, “gadf”}
Format to use. Default is “ogip”.
- checksumbool
If True checks both DATASUM and CHECKSUM cards in the file headers. Default is False.
- kwargsdict, optional
Keyword arguments passed to
MapDataset.read
.
- filename
- resample_energy_axis(energy_axis, name=None)#
Resample MapDatasetOnOff over reconstructed energy edges.
Counts are summed taking into account safe mask.
- Parameters
- energy_axis
MapAxis
New reco energy axis.
- namestr, optional
Name of the new dataset. Default is None.
- 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).
Default is “diff”.
- **kwargsdict, optional
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
, optional Energy bounds to compute the flux point for. Default is None.
- namestr, optional
Name of the sliced dataset. Default is None.
- 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)#
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_out
Map
Sliced map object.
- map_out
- stack(other, nan_to_num=True)#
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_numbool
Non-finite values are replaced by zero if True. Default is True.
- other
- stat_array()#
Statistic function value per bin given the current model parameters.
- stat_sum()#
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_image(name=None)#
Create images by summing over the reconstructed energy axis.
- Parameters
- namestr, optional
Name of the new dataset. Default is None.
- Returns
- dataset
MapDataset
orSpectrumDataset
Dataset integrated over non-spatial axes.
- dataset
- to_map_dataset(name=None)#
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
- dataset
MapDataset
Map dataset with cash statistics.
- dataset
- 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
- 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 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
- region
SkyRegion
Region from which to extract the spectrum.
- namestr, optional
Name of the new dataset. Default is None.
- region
- Returns
- dataset
MapDataset
The resulting reduced dataset.
- dataset
- to_spectrum_dataset(name=None)[source]#
Convert a SpectrumDatasetOnOff to a SpectrumDataset. The background model template is taken as alpha*counts_off.
- Parameters
- namestr, optional
Name of the new dataset. Default is None.
- Returns
- dataset
SpectrumDataset
SpectrumDataset with Cash statistic.
- dataset
- write(filename, overwrite=False, format='ogip', checksum=False)[source]#
Write spectrum dataset on off to file.
Can be serialised either as a
MapDataset
with aRegionGeom
following the GADF specifications, or as per the OGIP format. For OGIP formats specs, seeOGIPDatasetWriter
.- Parameters
- filename
Path
or str Filename to write to.
- overwritebool, optional
Overwrite existing file. Default is False.
- format{“ogip”, “ogip-sherpa”, “gadf”}
Format to use. Default is “ogip”.
- checksumbool
When True adds both DATASUM and CHECKSUM cards to the headers written to the file. Default is False.
- filename