FluxPoints¶
-
class
gammapy.spectrum.
FluxPoints
(table)[source]¶ Bases:
object
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: - table :
Table
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.spectrum import FluxPoints filename = '$GAMMAPY_DATA/tests/spectrum/flux_points/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:from astropy import units as u from astropy.table import Table from gammapy.spectrum import FluxPoints from gammapy.spectrum.models import PowerLaw table = Table() pwl = PowerLaw() e_ref = np.logspace(0, 2, 7) * u.TeV table['e_ref'] = e_ref table['dnde'] = pwl(e_ref) table.meta['SED_TYPE'] = 'dnde' flux_points = FluxPoints(table) flux_points.plot()
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.spectrum import FluxPoints table = Table.read('$GAMMAPY_DATA/tests/spectrum/flux_points/flux_points_ctb_37b.txt', format='ascii.csv', delimiter=' ', comment='#') table.meta['SED_TYPE'] = 'dnde' 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(table) flux_points.plot()
Attributes: - table :
Table
Table with flux point data
Attributes Summary
e_edges
Edges of the energy bin. e_max
Upper bound of energy bin. e_min
Lower bound of energy bin. e_ref
Reference energy. is_ul
sed_type
SED type (str). table_formatted
Return formatted version of the flux points table. Methods Summary
drop_ul
(self)Drop upper limit flux points. plot
(self[, ax, energy_unit, flux_unit, …])Plot flux points. plot_likelihood
(self[, ax, energy_unit, …])Plot likelihood SED profiles as a density plot.. read
(filename, \*\*kwargs)Read flux points. stack
(flux_points)Create flux points by stacking list of flux points. to_sed_type
(self, sed_type[, method, model, …])Convert to a different SED type (return new FluxPoints
).write
(self, filename, \*\*kwargs)Write flux points. Attributes Documentation
-
e_max
¶ Upper bound of energy bin.
Defined by
e_max
column intable
.Returns: - e_max :
Quantity
Upper bound of energy bin.
- e_max :
-
e_min
¶ Lower bound of energy bin.
Defined by
e_min
column inFluxPoints.table
.Returns: - e_min :
Quantity
Lower bound of energy bin.
- e_min :
-
e_ref
¶ Reference energy.
Defined by
e_ref
column inFluxPoints.table
or computed as log center, ife_min
ande_max
columns are present inFluxPoints.table
.Returns: - e_ref :
Quantity
Reference energy.
- e_ref :
-
is_ul
¶
-
sed_type
¶ SED type (str).
One of: {‘dnde’, ‘e2dnde’, ‘flux’, ‘eflux’}
-
table_formatted
¶ Return formatted version of the flux points table. Used for pretty printing
Methods Documentation
-
drop_ul
(self)[source]¶ Drop upper limit flux points.
Returns: - flux_points :
FluxPoints
Flux points with upper limit points removed.
Examples
>>> from gammapy.spectrum import FluxPoints >>> filename = '$GAMMAPY_DATA/tests/spectrum/flux_points/flux_points.fits' >>> flux_points = FluxPoints.read(filename) >>> print(flux_points) FluxPoints(sed_type="flux", n_points=24) >>> print(flux_points.drop_ul()) FluxPoints(sed_type="flux", n_points=19)
- flux_points :
-
plot
(self, ax=None, energy_unit='TeV', flux_unit=None, energy_power=0, **kwargs)[source]¶ Plot flux points.
Parameters: - ax :
Axes
Axis object to plot on.
- energy_unit : str,
Unit
, optional Unit of the energy axis
- flux_unit : str,
Unit
, optional Unit of the flux axis
- energy_power : int
Power of energy to multiply y axis with
- kwargs : dict
Keyword arguments passed to
matplotlib.pyplot.errorbar()
Returns: - ax :
Axes
Axis object
- ax :
-
plot_likelihood
(self, ax=None, energy_unit='TeV', add_cbar=True, y_values=None, y_unit=None, **kwargs)[source]¶ Plot likelihood SED profiles as a density plot..
Parameters: - ax :
Axes
Axis object to plot on.
- energy_unit : str,
Unit
, optional Unit of the energy axis
- y_values :
astropy.units.Quantity
Array of y-values to use for the likelihood profile evaluation.
- y_unit : str or
astropy.units.Unit
Unit to use for the y-axis.
- add_cbar : bool
Whether to add a colorbar to the plot.
- kwargs : dict
Keyword arguments passed to
matplotlib.pyplot.pcolormesh()
Returns: - ax :
Axes
Axis object
- ax :
-
classmethod
read
(filename, **kwargs)[source]¶ Read flux points.
Parameters: - filename : str
Filename
- kwargs : dict
Keyword arguments passed to
astropy.table.Table.read
.
-
classmethod
stack
(flux_points)[source]¶ 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: - flux_points : list of
FluxPoints
List of flux points to stack.
Returns: - flux_points :
FluxPoints
Flux points without upper limit points.
- flux_points : list of
-
to_sed_type
(self, sed_type, method='log_center', model=None, pwl_approx=False)[source]¶ Convert to a different SED type (return new
FluxPoints
).See: https://ui.adsabs.harvard.edu/abs/1995NIMPA.355..541L for details on the
'lafferty'
method.Parameters: - sed_type : {‘dnde’}
SED type to convert to.
- model :
SpectralModel
Spectral model assumption. Note that the value of the amplitude parameter does not matter. Still it is recommended to use something with the right scale and units. E.g.
amplitude = 1e-12 * u.Unit('cm-2 s-1 TeV-1')
- method : {‘lafferty’, ‘log_center’, ‘table’}
Flux points
e_ref
estimation method:'laferty'
Lafferty & Wyatt model-based e_ref'log_center'
log bin center e_ref'table'
using column ‘e_ref’ from input flux_points
- pwl_approx : bool
Use local power law appoximation at e_ref to compute differential flux from the integral flux. This method is used by the Fermi-LAT catalogs.
Returns: - flux_points :
FluxPoints
Flux points including differential quantity columns
dnde
anddnde_err
(optional),dnde_ul
(optional).
Examples
>>> from gammapy.spectrum import FluxPoints >>> from gammapy.spectrum.models import PowerLaw >>> filename = '$GAMMAPY_DATA/tests/spectrum/flux_points/flux_points.fits' >>> flux_points = FluxPoints.read(filename) >>> model = PowerLaw(index=2.2) >>> flux_points_dnde = flux_points.to_sed_type('dnde', model=model)
-
write
(self, filename, **kwargs)[source]¶ Write flux points.
Parameters: - filename : str
Filename
- kwargs : dict
Keyword arguments passed to
astropy.table.Table.write
.
- Format