FluxPoints#
- class gammapy.estimators.FluxPoints[source]#
Bases:
FluxMaps
Flux points container.
The supported formats are described here: SED.
In summary, the following formats and minimum required columns are:
Format
dnde
: columnse_ref
anddnde
Format
e2dnde
: columnse_ref
,e2dnde
Format
flux
: columnse_min
,e_max
,flux
Format
eflux
: columnse_min
,e_max
,eflux
- Parameters:
- datadict of
RegionNDMap
The maps dictionary. Expected entries are the following:
norm : the norm factor.
norm_err : optional, the error on the norm factor.
norm_errn : optional, the negative error on the norm factor.
norm_errp : optional, the positive error on the norm factor.
norm_ul : optional, the upper limit on the norm factor.
norm_scan : optional, the norm values of the test statistic scan.
stat_scan : optional, the test statistic scan values.
ts : optional, the delta test statistic associated with the flux value.
sqrt_ts : optional, the square root of the test statistic, when relevant.
success : optional, a boolean tagging the validity of the estimation.
n_dof : optional, the number of degrees of freedom used in TS computation
alpha : optional, normalisation factor to accounts for differences between the test region and the background
acceptance_off : optional, acceptance from the off region
acceptance_on : optional, acceptance from the on region
- reference_model
SkyModel
, optional The 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.
- metadict, optional
Dict of metadata.
- gti
GTI
, optional Maps GTI information.
- filter_success_nanboolean, optional
Set fitted norm and error to NaN when the fit has not succeeded.
- datadict of
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 ofastropy.table.Table
, which contains the required columns, such as'e_ref'
and'dnde'
. The correspondingsed_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
Methods Summary
from_table
(table[, sed_type, format, ...])Create flux points from a table.
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.
to_table
([sed_type, format, formatted])Create table for a given SED type.
write
(filename[, sed_type, format, ...])Write flux points.
Methods Documentation
- 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:
- table
Table
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_model
SpectralModel
, optional Reference spectral model. Default is None.
- gti
GTI
, optional Good time intervals. Default is None.
- table
- Returns:
- flux_points
FluxPoints
Flux points.
- flux_points
- plot(ax=None, sed_type=None, energy_power=0, time_format='iso', **kwargs)[source]#
Plot flux points.
- Parameters:
- ax
Axes
, 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
.
- ax
- Returns:
- ax
Axes
Axis object.
- ax
- plot_ts_profiles(ax=None, sed_type=None, add_cbar=True, **kwargs)[source]#
Plot fit statistic SED profiles as a density plot.
- Parameters:
- ax
Axes
, 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
.
- ax
- Returns:
- ax
Axes
Axis object.
- ax
- classmethod read(filename, sed_type=None, format=None, reference_model=None, checksum=False, table_format='ascii.ecsv', **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_model
SpectralModel
Reference spectral model.
- checksumbool
If True checks both DATASUM and CHECKSUM cards in the file headers. Default is False.
- table_formatstr
Format string for the ~astropy.Table object. Default is “ascii.ecsv”
- **kwargsdict, optional
Keyword arguments passed to
astropy.table.Table.read
.
- Returns:
- flux_points
FluxPoints
Flux points.
- 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_ulfloat, optional
Number of sigma to use for upper limit computation. Must be a positive value. Default is 2.
- **kwargsdict, optional
Keyword arguments passed to
brentq
.
- Returns:
- flux_points
FluxPoints
A new FluxPoints object with modified upper limits.
- flux_points
Examples
>>> from gammapy.estimators import FluxPoints >>> filename = "$GAMMAPY_DATA/estimators/crab_hess_fp/crab_hess_fp.fits" >>> flux_points = FluxPoints.read(filename) >>> flux_points_recomputed = flux_points.recompute_ul(n_sigma_ul=4) >>> print(f'{flux_points.meta["n_sigma_ul"]}, {float(flux_points.flux_ul.data[1][0]):.2e}') 3.0, 3.99e-11 >>> print(f'{flux_points_recomputed.meta["n_sigma_ul"]}, {float(flux_points_recomputed.flux_ul.data[1][0]):.2e}') 4, 4.25e-11
- 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_new
MapAxis
orTimeMapAxis
The new axis to resample along
- axis_new
- Returns:
- flux_points
FluxPoints
A new FluxPoints object with modified axis.
- flux_points
- 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:
- table
Table
Flux points table.
- 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.
- __init__(data, reference_model, meta=None, gti=None, filter_success_nan=True)#
- classmethod __new__(*args, **kwargs)#