FluxPoints#

class gammapy.estimators.FluxPoints(data, reference_model, meta=None, gti=None, filter_success_nan=True)[source]#

Bases: gammapy.estimators.map.core.FluxMaps

Flux points container.

The supported formats are described here: SED.

In summary, the following formats and minimum required columns are:

  • Format dnde: columns e_ref and dnde

  • Format e2dnde: columns e_ref, e2dnde

  • Format flux: columns e_min, e_max, flux

  • Format eflux: columns e_min, e_max, eflux

Parameters
tableTable

Table with flux point data.

Examples

The FluxPoints object is most easily created by reading a file with flux points given in one of the formats documented above:

>>> from gammapy.estimators import FluxPoints
>>> filename = '$GAMMAPY_DATA/hawc_crab/HAWC19_flux_points.fits'
>>> flux_points = FluxPoints.read(filename)
>>> flux_points.plot() 

An instance of FluxPoints can also be created by passing an instance of astropy.table.Table, which contains the required columns, such as 'e_ref' and 'dnde'. The corresponding sed_type has to be defined in the meta data of the table:

>>> import numpy as np
>>> from astropy import units as u
>>> from astropy.table import Table
>>> from gammapy.estimators import FluxPoints
>>> from gammapy.modeling.models import PowerLawSpectralModel
>>> table = Table()
>>> pwl = PowerLawSpectralModel()
>>> e_ref = np.geomspace(1, 100, 7) * u.TeV
>>> table["e_ref"] = e_ref
>>> table["dnde"] = pwl(e_ref)
>>> table["dnde_err"] = pwl.evaluate_error(e_ref)[0]
>>> table.meta["SED_TYPE"] = "dnde"
>>> flux_points = FluxPoints.from_table(table)
>>> flux_points.plot(sed_type="flux") 

If you have flux points in a different data format, the format can be changed by renaming the table columns and adding meta data:

>>> from astropy import units as u
>>> from astropy.table import Table
>>> from gammapy.estimators import FluxPoints
>>> from gammapy.utils.scripts import make_path
>>> filename = make_path('$GAMMAPY_DATA/tests/spectrum/flux_points/flux_points_ctb_37b.txt')
>>> table = Table.read(filename ,format='ascii.csv', delimiter=' ', comment='#')
>>> table.rename_column('Differential_Flux', 'dnde')
>>> table['dnde'].unit = 'cm-2 s-1 TeV-1'
>>> table.rename_column('lower_error', 'dnde_errn')
>>> table['dnde_errn'].unit = 'cm-2 s-1 TeV-1'
>>> table.rename_column('upper_error', 'dnde_errp')
>>> table['dnde_errp'].unit = 'cm-2 s-1 TeV-1'
>>> table.rename_column('E', 'e_ref')
>>> table['e_ref'].unit = 'TeV'
>>> flux_points = FluxPoints.from_table(table, sed_type="dnde")
>>> flux_points.plot(sed_type="e2dnde") 

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
tableTable

Table with flux point data.

Attributes Summary

available_quantities

Available quantities

counts

Predicted counts null hypothesis.

dnde

Return differential flux (dnde) SED values.

dnde_err

Return differential flux (dnde) SED errors.

dnde_errn

Return differential flux (dnde) SED negative errors.

dnde_errp

Return differential flux (dnde) SED positive errors.

dnde_ref

Reference differential flux.

dnde_ul

Return differential flux (dnde) SED upper limit.

e2dnde

Return differential energy flux (e2dnde) SED values.

e2dnde_err

Return differential energy flux (e2dnde) SED errors.

e2dnde_errn

Return differential energy flux (e2dnde) SED negative errors.

e2dnde_errp

Return differential energy flux (e2dnde) SED positive errors.

e2dnde_ref

Reference differential flux * energy ** 2.

e2dnde_ul

Return differential energy flux (e2dnde) SED upper limit.

eflux

Return energy flux (eflux) SED values.

eflux_err

Return energy flux (eflux) SED errors.

eflux_errn

Return energy flux (eflux) SED negative errors.

eflux_errp

Return energy flux (eflux) SED positive errors.

eflux_ref

Reference energy flux.

eflux_ul

Return energy flux (eflux) SED upper limits.

energy_axis

Energy axis as a MapAxis.

energy_max

Energy maximum.

energy_min

Energy minimum.

energy_ref

Reference energy.

filter_success_nan

flux

Return integral flux (flux) SED values.

flux_err

Return integral flux (flux) SED values.

flux_errn

Return integral flux (flux) SED negative errors.

flux_errp

Return integral flux (flux) SED positive errors.

flux_ref

Reference integral flux.

flux_sensitivity

Sensitivity given as the flux for which the significance is self.meta["n_sigma_sensitivity].

flux_ul

Return integral flux (flux) SED upper limits.

geom

Reference map geometry as a Geom.

has_any_ts

Whether the flux estimate has either sqrt(TS) or test statistic defined.

has_stat_profiles

Whether the flux estimate has test statistic profiles.

has_success

Whether the flux estimate has the fit status.

has_ul

Whether the flux estimate has norm_ul defined.

is_convertible_to_flux_sed_type

Check whether differential SED type is convertible to integral SED type.

is_ul

Whether data is an upper limit.

n_dof

Number of degrees of freedom of the fit per energy bin.

n_sigma

n sigma

n_sigma_ul

n sigma UL.

niter

Number of iterations of fit.

norm

Norm values.

norm_err

Norm error.

norm_errn

Negative norm error.

norm_errp

Positive norm error.

norm_sensitivity

Norm sensitivity.

norm_ul

Norm upper limit.

npred

Predicted counts from best fit hypothesis.

npred_background

Predicted background counts from best fit hypothesis.

npred_excess

Predicted excess count from best fit hypothesis.

npred_excess_err

Predicted excess counts error.

npred_excess_errn

Predicted excess counts negative error.

npred_excess_errp

Predicted excess counts positive error.

npred_excess_ref

Predicted excess reference counts.

npred_excess_ul

Predicted excess counts upper limits.

reference_model

Reference model as a SkyModel.

reference_model_default

Reference model default as a SkyModel.

reference_spectral_model

Reference spectral model as a SpectralModel

sed_type_init

Initial SED type.

sed_type_plot_default

Initial SED type.

sqrt_ts

sqrt(TS) as defined by:

sqrt_ts_threshold_ul

sqrt(TS) threshold for upper limits.

stat

Fit statistic value.

stat_null

Fit statistic value for the null hypothesis.

stat_scan

Fit statistic scan value.

success

Fit success flag.

ts

Test statistic map as a Map object.

ts_scan

Test statistic scan as a Map object

Methods Summary

all_quantities(sed_type)

All quantities allowed for a given SED type.

from_hdulist(hdulist[, hdu_bands, sed_type, ...])

Create flux map dataset from list of HDUs.

from_maps(maps[, sed_type, reference_model, ...])

Create FluxMaps from a dictionary of maps.

from_stack(maps, axis[, meta])

Create flux points by stacking list of flux points.

from_table(table[, sed_type, format, ...])

Create flux points from a table.

get_flux_points([position])

Extract flux point at a given position.

iter_by_axis(axis_name)

Create a set of FluxMaps by splitting along an axis.

plot([ax, sed_type, energy_power, time_format])

Plot flux points.

plot_ts_profiles([ax, sed_type, add_cbar])

Plot fit statistic SED profiles as a density plot.

read(filename[, sed_type, format, ...])

Read precomputed flux points.

recompute_ul([n_sigma_ul])

Recompute upper limits corresponding to the given value.

resample_axis(axis_new)

Rebin the flux point object along the new axis.

slice_by_coord(slices)

Slice flux maps by coordinate values

slice_by_energy(energy_min, energy_max)

Slice flux maps by coordinate values along the energy axis.

slice_by_idx(slices)

Slice flux maps by index.

slice_by_time(time_min, time_max)

Slice flux maps by coordinate values along the time axis.

to_hdulist([sed_type, hdu_bands])

Convert flux map to list of HDUs.

to_maps([sed_type])

Return maps in a given SED type.

to_table([sed_type, format, formatted])

Create table for a given SED type.

write(filename[, sed_type, format, ...])

Write flux points.

Attributes Documentation

available_quantities#

Available quantities

counts#

Predicted counts null hypothesis.

dnde#

Return differential flux (dnde) SED values.

dnde_err#

Return differential flux (dnde) SED errors.

dnde_errn#

Return differential flux (dnde) SED negative errors.

dnde_errp#

Return differential flux (dnde) SED positive errors.

dnde_ref#

Reference differential flux.

dnde_ul#

Return differential flux (dnde) SED upper limit.

e2dnde#

Return differential energy flux (e2dnde) SED values.

e2dnde_err#

Return differential energy flux (e2dnde) SED errors.

e2dnde_errn#

Return differential energy flux (e2dnde) SED negative errors.

e2dnde_errp#

Return differential energy flux (e2dnde) SED positive errors.

e2dnde_ref#

Reference differential flux * energy ** 2.

e2dnde_ul#

Return differential energy flux (e2dnde) SED upper limit.

eflux#

Return energy flux (eflux) SED values.

eflux_err#

Return energy flux (eflux) SED errors.

eflux_errn#

Return energy flux (eflux) SED negative errors.

eflux_errp#

Return energy flux (eflux) SED positive errors.

eflux_ref#

Reference energy flux.

eflux_ul#

Return energy flux (eflux) SED upper limits.

energy_axis#

Energy axis as a MapAxis.

energy_max#

Energy maximum.

Returns
energy_maxQuantity

Upper bound of energy bin.

energy_min#

Energy minimum.

Returns
energy_minQuantity

Lower bound of energy bin.

energy_ref#

Reference energy.

Defined by energy_ref column in FluxPoints.table or computed as log center, if energy_min and energy_max columns are present in FluxEstimate.data.

Returns
energy_refQuantity

Reference energy.

filter_success_nan#
flux#

Return integral flux (flux) SED values.

flux_err#

Return integral flux (flux) SED values.

flux_errn#

Return integral flux (flux) SED negative errors.

flux_errp#

Return integral flux (flux) SED positive errors.

flux_ref#

Reference integral flux.

flux_sensitivity#

Sensitivity given as the flux for which the significance is self.meta["n_sigma_sensitivity].

flux_ul#

Return integral flux (flux) SED upper limits.

geom#

Reference map geometry as a Geom.

has_any_ts#

Whether the flux estimate has either sqrt(TS) or test statistic defined.

has_stat_profiles#

Whether the flux estimate has test statistic profiles.

has_success#

Whether the flux estimate has the fit status.

has_ul#

Whether the flux estimate has norm_ul defined.

is_convertible_to_flux_sed_type#

Check whether differential SED type is convertible to integral SED type.

is_ul#

Whether data is an upper limit.

n_dof#

Number of degrees of freedom of the fit per energy bin.

n_sigma#

n sigma

n_sigma_ul#

n sigma UL.

niter#

Number of iterations of fit.

norm#

Norm values.

norm_err#

Norm error.

norm_errn#

Negative norm error.

norm_errp#

Positive norm error.

norm_sensitivity#

Norm sensitivity.

norm_ul#

Norm upper limit.

npred#

Predicted counts from best fit hypothesis.

npred_background#

Predicted background counts from best fit hypothesis.

npred_excess#

Predicted excess count from best fit hypothesis.

npred_excess_err#

Predicted excess counts error.

npred_excess_errn#

Predicted excess counts negative error.

npred_excess_errp#

Predicted excess counts positive error.

npred_excess_ref#

Predicted excess reference counts.

npred_excess_ul#

Predicted excess counts upper limits.

reference_model#

Reference model as a SkyModel.

reference_model_default#

Reference model default as a SkyModel.

reference_spectral_model#

Reference spectral model as a SpectralModel

sed_type_init#

Initial SED type.

sed_type_plot_default#

Initial SED type.

sqrt_ts#

sqrt(TS) as defined by:

\[\begin{split}\sqrt{TS} = \left \{ \begin{array}{ll} -\sqrt{TS} & : \text{if} \ norm < 0 \\ \sqrt{TS} & : \text{else} \end{array} \right.\end{split}\]
Returns
sqrt_tsMap

sqrt(TS) map.

sqrt_ts_threshold_ul#

sqrt(TS) threshold for upper limits.

stat#

Fit statistic value.

stat_null#

Fit statistic value for the null hypothesis.

stat_scan#

Fit statistic scan value.

success#

Fit success flag.

ts#

Test statistic map as a Map object.

ts_scan#

Test statistic scan as a Map object

Methods Documentation

static all_quantities(sed_type)#

All quantities allowed for a given SED type.

Parameters
sed_type{“likelihood”, “dnde”, “e2dnde”, “flux”, “eflux”}

SED type.

Returns
listlist of str

All allowed quantities for a given SED type.

classmethod from_hdulist(hdulist, hdu_bands=None, sed_type=None, checksum=False)#

Create flux map dataset from list of HDUs.

Parameters
hdulistHDUList

List of HDUs.

hdu_bandsstr, optional

Name of the HDU with the BANDS table. Default is ‘BANDS’ If set to None, each map should have its own hdu_band. Default is None.

sed_type{“dnde”, “flux”, “e2dnde”, “eflux”, “likelihood”}, optional

Sed type. Default is None.

Returns
flux_mapsFluxMaps

Flux maps object.

classmethod from_maps(maps, sed_type=None, reference_model=None, gti=None, meta=None)#

Create FluxMaps from a dictionary of maps.

Parameters
mapsMaps

Maps object containing the input maps.

sed_typestr, optional

SED type of the input maps. If None, set to “likelihood”. Default is None.

reference_modelSkyModel, optional

Reference model to use for conversions. If None, a model consisting of a point source with a power law spectrum of index 2 is assumed. Default is None.

gtiGTI, optional

Maps GTI information. Default is None.

metadict

Meta dictionary.

Returns
flux_mapsFluxMaps

Flux maps object.

classmethod from_stack(maps, axis, meta=None)#

Create flux points by stacking list of flux points.

The first FluxPoints object in the list is taken as a reference to infer column names and units for the stacked object.

Parameters
mapslist of FluxMaps

List of maps to stack.

axisMapAxis

New axis to create.

metadict, optional

Metadata of the resulting flux points. Default is None.

Returns
flux_mapsFluxMaps

Stacked flux maps along axis.

classmethod from_table(table, sed_type=None, format=None, reference_model=None, gti=None)[source]#

Create flux points from a table. The table column names must be consistent with the sed_type.

Parameters
tableTable

Table.

sed_type{“dnde”, “flux”, “eflux”, “e2dnde”, “likelihood”}, optional

SED type. Default is None.

format{“gadf-sed”, “lightcurve”, “profile”}, optional

Table format. If None, it is extracted from the table column content. Default is None.

reference_modelSpectralModel, optional

Reference spectral model. Default is None.

gtiGTI, optional

Good time intervals. Default is None.

Returns
flux_pointsFluxPoints

Flux points.

get_flux_points(position=None)#

Extract flux point at a given position.

Parameters
positionSkyCoord

Position where the flux points are extracted.

Returns
flux_pointsFluxPoints

Flux points object.

iter_by_axis(axis_name)#

Create a set of FluxMaps by splitting along an axis.

Parameters
axis_namestr

Name of the axis to split on.

Returns
flux_mapsFluxMap

FluxMap iteration.

plot(ax=None, sed_type=None, energy_power=0, time_format='iso', **kwargs)[source]#

Plot flux points.

Parameters
axAxes, optional

Axis object to plot on. Default is None.

sed_type{“dnde”, “flux”, “eflux”, “e2dnde”}, optional

SED type. Default is None.

energy_powerfloat, optional

Power of energy to multiply flux axis with. Default is 0.

time_format{“iso”, “mjd”}

Used time format is a time axis is present. Default is “iso”.

**kwargsdict, optional

Keyword arguments passed to plot.

Returns
axAxes

Axis object.

plot_ts_profiles(ax=None, sed_type=None, add_cbar=True, **kwargs)[source]#

Plot fit statistic SED profiles as a density plot.

Parameters
axAxes, optional

Axis object to plot on. Default is None.

sed_type{“dnde”, “flux”, “eflux”, “e2dnde”}, optional

SED type. Default is None.

add_cbarbool, optional

Whether to add a colorbar to the plot. Default is True.

**kwargsdict, optional

Keyword arguments passed to pcolormesh.

Returns
axAxes

Axis object.

classmethod read(filename, sed_type=None, format=None, reference_model=None, checksum=False, **kwargs)[source]#

Read precomputed flux points.

Parameters
filenamestr

Filename.

sed_type{“dnde”, “flux”, “eflux”, “e2dnde”, “likelihood”}

SED type.

format{“gadf-sed”, “lightcurve”, “profile”}, optional

Format string. If None, the format is extracted from the input. Default is None.

reference_modelSpectralModel

Reference spectral model.

checksumbool

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

**kwargsdict, optional

Keyword arguments passed to astropy.table.Table.read.

Returns
flux_pointsFluxPoints

Flux points.

recompute_ul(n_sigma_ul=2, **kwargs)[source]#

Recompute upper limits corresponding to the given value. The pre-computed statistic profiles must exist for the re-computation.

Parameters
n_sigma_ulint

Number of sigma to use for upper limit computation. Default is 2.

**kwargsdict, optional

Keyword arguments passed to brentq.

Returns
flux_pointsFluxPoints

A new FluxPoints object with modified upper limits.

Examples

>>> from gammapy.estimators import FluxPoints
>>> filename = '$GAMMAPY_DATA/tests/spectrum/flux_points/binlike.fits'
>>> flux_points = FluxPoints.read(filename)
>>> flux_points_recomputed = flux_points.recompute_ul(n_sigma_ul=3)
>>> print(flux_points.meta["n_sigma_ul"], flux_points.flux_ul.data[0])
2.0 [[3.95451985e-09]]
>>> print(flux_points_recomputed.meta["n_sigma_ul"], flux_points_recomputed.flux_ul.data[0])
3 [[6.22245374e-09]]
resample_axis(axis_new)[source]#

Rebin the flux point object along the new axis. The log-likelihoods profiles in each bin are summed to compute the resultant quantities. Stat profiles must be present on the fluxpoints object for this method to work.

For now, works only for flat fluxpoints.

Parameters
axis_newMapAxis or TimeMapAxis

The new axis to resample along

Returns
flux_pointsFluxPoints

A new FluxPoints object with modified axis.

slice_by_coord(slices)#

Slice flux maps by coordinate values

Parameters
slicesdict

Dict of axes names and astropy.Quantity or astropy.Time 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.

Returns
flux_mapsFluxMaps

Sliced flux maps object.

slice_by_energy(energy_min, energy_max)#

Slice flux maps by coordinate values along the energy axis.

Parameters
energy_min, energy_maxQuantity

Energy bounds used to slice the flux map.

Returns
flux_mapsFluxMaps

Sliced flux maps object.

slice_by_idx(slices)#

Slice flux maps by index.

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.

Returns
flux_mapsFluxMaps

Sliced flux maps object.

slice_by_time(time_min, time_max)#

Slice flux maps by coordinate values along the time axis.

Parameters
time_min, time_maxTime

Time bounds used to slice the flux map.

Returns
flux_mapsFluxMaps

Sliced flux maps object.

to_hdulist(sed_type=None, hdu_bands=None)#

Convert flux map to list of HDUs.

For now, one cannot export the reference model.

Parameters
sed_typestr, optional

SED type to convert to. If None, set to “likelihood”. Default is None.

hdu_bandsstr, optional

Name of the HDU with the BANDS table. Default is ‘BANDS’ If set to None, each map will have its own hdu_band. Default is None.

Returns
hdulistHDUList

Map dataset list of HDUs.

to_maps(sed_type=None)#

Return maps in a given SED type.

Parameters
sed_type{“likelihood”, “dnde”, “e2dnde”, “flux”, “eflux”}, optional

SED type to convert to. If None, set to Likelihood. Default is None.

Returns
mapsMaps

Maps object containing the requested maps.

to_table(sed_type=None, format=None, formatted=False)[source]#

Create table for a given SED type.

Parameters
sed_type{“likelihood”, “dnde”, “e2dnde”, “flux”, “eflux”}

SED type to convert to. Default is likelihood.

format{“gadf-sed”, “lightcurve”, “binned-time-series”, “profile”}, optional

Format specification. The following formats are supported:

  • “gadf-sed”: format for SED flux points see SED for details

  • “lightcurve”: Gammapy internal format to store energy dependent lightcurves. Basically a generalisation of the “gadf” format, but currently there is no detailed documentation available.

  • “binned-time-series”: table format support by Astropy’s BinnedTimeSeries.

  • “profile”: Gammapy internal format to store energy dependent flux profiles. Basically a generalisation of the “gadf” format, but currently there is no detailed documentation available.

If None, the format will be guessed by looking at the axes that are present in the object. Default is None.

formattedbool

Formatted version with column formats applied. Numerical columns are formatted to .3f and .3e respectively. Default is False.

Returns
tableTable

Flux points table.

Examples

This is how to read and plot example flux points:

>>> from gammapy.estimators import FluxPoints
>>> fp = FluxPoints.read("$GAMMAPY_DATA/hawc_crab/HAWC19_flux_points.fits")
>>> table = fp.to_table(sed_type="flux", formatted=True)
>>> print(table[:2])
e_ref e_min e_max     flux      flux_err    flux_ul      ts    sqrt_ts is_ul
 TeV   TeV   TeV  1 / (s cm2) 1 / (s cm2) 1 / (s cm2)
----- ----- ----- ----------- ----------- ----------- -------- ------- -----
1.334 1.000 1.780   1.423e-11   3.135e-13         nan 2734.000  52.288 False
2.372 1.780 3.160   5.780e-12   1.082e-13         nan 4112.000  64.125 False
write(filename, sed_type=None, format=None, overwrite=False, checksum=False)[source]#

Write flux points.

Parameters
filenamestr

Filename.

sed_type{“dnde”, “flux”, “eflux”, “e2dnde”, “likelihood”}, optional

SED type. Default is None.

format{“gadf-sed”, “lightcurve”, “binned-time-series”, “profile”}, optional

Format specification. The following formats are supported:

  • “gadf-sed”: format for SED flux points see SED for details

  • “lightcurve”: Gammapy internal format to store energy dependent lightcurves. Basically a generalisation of the “gadf” format, but currently there is no detailed documentation available.

  • “binned-time-series”: table format support by Astropy’s BinnedTimeSeries.

  • “profile”: Gammapy internal format to store energy dependent flux profiles. Basically a generalisation of the “gadf” format, but currently there is no detailed documentation available.

If None, the format will be guessed by looking at the axes that are present in the object. Default is None.

overwritebool, optional

Overwrite existing file. Default is False.

checksumbool, optional

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