FluxPointsDataset#
- class gammapy.datasets.FluxPointsDataset(models=None, data=None, mask_fit=None, mask_safe=None, name=None, meta_table=None)[source]#
Bases:
gammapy.datasets.core.Dataset
Bundle a set of flux points with a parametric model, to compute fit statistic function using chi2 statistics.
For more information see Datasets (DL4).
- Parameters
- models
Models
Models (only spectral part needs to be set).
- data
FluxPoints
Flux points. Must be sorted along the energy axis.
- mask_fit
numpy.ndarray
Mask to apply for fitting.
- mask_safe
numpy.ndarray
Mask defining the safe data range. By default, upper limit values are excluded.
- meta_table
Table
Table listing information on observations used to create the dataset. One line per observation for stacked datasets.
- models
Examples
Load flux points from file and fit with a power-law model:
>>> from gammapy.modeling import Fit >>> from gammapy.modeling.models import PowerLawSpectralModel, SkyModel >>> from gammapy.estimators import FluxPoints >>> from gammapy.datasets import FluxPointsDataset >>> filename = "$GAMMAPY_DATA/tests/spectrum/flux_points/diff_flux_points.fits" >>> dataset = FluxPointsDataset.read(filename) >>> model = SkyModel(spectral_model=PowerLawSpectralModel()) >>> dataset.models = model
Make the fit:
>>> fit = Fit() >>> result = fit.run([dataset]) >>> print(result) OptimizeResult backend : minuit method : migrad success : True message : Optimization terminated successfully. nfev : 135 total stat : 25.21 CovarianceResult backend : minuit method : hesse success : True message : Hesse terminated successfully.
>>> print(result.parameters.to_table()) type name value unit ... frozen is_norm link prior ---- --------- ---------- -------------- ... ------ ------- ---- ----- index 2.2159e+00 ... False False amplitude 2.1619e-13 TeV-1 s-1 cm-2 ... False True reference 1.0000e+00 TeV ... True False
Note: In order to reproduce the example, you need the tests datasets folder. You may download it with the command:
gammapy download datasets --tests --out $GAMMAPY_DATA
Attributes Summary
Good time interval info (
GTI
).Combined fit and safe mask.
Methods Summary
copy
([name])A deep copy.
Shape of the flux points data (tuple).
Compute predicted flux.
from_dict
(data, **kwargs)Create flux point dataset from dict.
plot_fit
([ax_spectrum, ax_residuals, ...])Plot flux points, best fit model and residuals in two panels.
plot_residuals
([ax, method])Plot flux point residuals.
plot_spectrum
([ax, kwargs_fp, kwargs_model, ...])Plot flux points and model.
read
(filename[, name])Read pre-computed flux points and create a dataset.
residuals
([method])Compute flux point residuals.
Fit statistic array.
stat_sum
()Total statistic given the current model parameters and priors.
to_dict
()Convert to dict for YAML serialization.
write
(filename[, overwrite, checksum])Write flux point dataset to file.
Attributes Documentation
- gti#
Good time interval info (
GTI
).
- mask#
Combined fit and safe mask.
- models#
- name#
- stat_type = 'chi2'#
- tag = 'FluxPointsDataset'#
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 from_dict(data, **kwargs)[source]#
Create flux point dataset from dict.
- Parameters
- datadict
Dictionary containing data to create dataset from.
- Returns
- ——-
- dataset
FluxPointsDataset
Flux point datasets.
- plot_fit(ax_spectrum=None, ax_residuals=None, kwargs_spectrum=None, kwargs_residuals=None, axis_name='energy')[source]#
Plot flux points, best fit model and residuals in two panels.
Calls
plot_spectrum
andplot_residuals
.- Parameters
- ax_spectrum
Axes
, optional Axes to plot flux points and best fit model on. Default is None.
- ax_residuals
Axes
, optional Axes to plot residuals on. Default is None.
- kwargs_spectrumdict, optional
Keyword arguments passed to
plot_spectrum
. Default is None.- kwargs_residualsdict, optional
Keyword arguments passed to
plot_residuals
. Default is None.- Returns
- ——-
- ax_spectrum, ax_residuals
Axes
Flux points, best fit model and residuals plots.
- ax_spectrum
Examples
>>> from gammapy.modeling.models import PowerLawSpectralModel, SkyModel >>> from gammapy.estimators import FluxPoints >>> from gammapy.datasets import FluxPointsDataset
>>> #load precomputed flux points >>> filename = "$GAMMAPY_DATA/tests/spectrum/flux_points/diff_flux_points.fits" >>> flux_points = FluxPoints.read(filename) >>> model = SkyModel(spectral_model=PowerLawSpectralModel()) >>> dataset = FluxPointsDataset(model, flux_points) >>> #configuring optional parameters >>> kwargs_spectrum = {"kwargs_model": {"color":"red", "ls":"--"}, "kwargs_fp":{"color":"green", "marker":"o"}} >>> kwargs_residuals = {"color": "blue", "markersize":4, "marker":'s', } >>> dataset.plot_fit(kwargs_residuals=kwargs_residuals, kwargs_spectrum=kwargs_spectrum)
- plot_residuals(ax=None, method='diff', **kwargs)[source]#
Plot flux point residuals.
- Parameters
- ax
Axes
, optional Axes to plot on. Default is None.
- method{“diff”, “diff/model”}
Normalization used to compute the residuals, see
FluxPointsDataset.residuals
. Default is “diff”.- **kwargsdict
Keyword arguments passed to
errorbar
.
- ax
- Returns
- ax
Axes
Axes object.
- ax
- plot_spectrum(ax=None, kwargs_fp=None, kwargs_model=None, axis_name='energy')[source]#
Plot flux points and model.
- Parameters
- ax
Axes
, optional Axes to plot on. Default is None.
- kwargs_fpdict, optional
Keyword arguments passed to
gammapy.estimators.FluxPoints.plot
to configure the plot style. Default is None.- kwargs_modeldict, optional
Keyword arguments passed to
gammapy.modeling.models.SpectralModel.plot
andgammapy.modeling.models.SpectralModel.plot_error
to configure the plot style. Default is None.- axis_namestr
Axis along which to plot the flux points for multiple axes. Default is energy.
- ax
- Returns
- ax
Axes
Axes object.
- ax
Examples
>>> from gammapy.modeling.models import PowerLawSpectralModel, SkyModel >>> from gammapy.estimators import FluxPoints >>> from gammapy.datasets import FluxPointsDataset
>>> #load precomputed flux points >>> filename = "$GAMMAPY_DATA/tests/spectrum/flux_points/diff_flux_points.fits" >>> flux_points = FluxPoints.read(filename) >>> model = SkyModel(spectral_model=PowerLawSpectralModel()) >>> dataset = FluxPointsDataset(model, flux_points) >>> #configuring optional parameters >>> kwargs_model = {"color":"red", "ls":"--"} >>> kwargs_fp = {"color":"green", "marker":"o"} >>> dataset.plot_spectrum(kwargs_fp=kwargs_fp, kwargs_model=kwargs_model)
- classmethod read(filename, name=None)[source]#
Read pre-computed flux points and create a dataset.
- Parameters
- filenamestr
Filename to read from.
- namestr, optional
Name of the new dataset. Default is None.
- Returns
- ——-
- dataset
FluxPointsDataset
FluxPointsDataset.
- residuals(method='diff')[source]#
Compute flux point residuals.
- Parameters
- method: {“diff”, “diff/model”}
- Method used to compute the residuals. Available options are:
diff
(default): data - model.diff/model
: (data - model) / model.
Default is “diff”.
- Returns
- residuals
ndarray
Residuals array.
- residuals
- stat_sum()#
Total statistic given the current model parameters and priors.
- to_dict()#
Convert to dict for YAML serialization.
- write(filename, overwrite=False, checksum=False, **kwargs)[source]#
Write flux point dataset to file.
- 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 FITS file. Applies only if filename has .fits suffix. Default is False.
- **kwargsdict, optional
Keyword arguments passed to
write
.