Note
Go to the end to download the full example code or to run this example in your browser via Binder
Source catalogs#
Access and explore thew most common gamma-ray source catalogs.
Introduction#
catalog
provides convenient access to common gamma-ray
source catalogs. This module is mostly independent from the rest of
Gammapy. Typically you use it to compare new analyses against catalog
results, e.g. overplot the spectral model, or compare the source
position.
Moreover as creating a source model and flux points for a given catalog
from the FITS table is tedious, catalog
has this already
implemented. So you can create initial source models for your analyses.
This is very common for Fermi-LAT, to start with a catalog model. For
TeV analysis, especially in crowded Galactic regions, using the HGPS,
gamma-cat or 2HWC catalog in this way can also be useful.
In this tutorial you will learn how to:
List available catalogs
Load a catalog
Access the source catalog table data
Select a catalog subset or a single source
Get source spectral and spatial models
Get flux points (if available)
Get lightcurves (if available)
Access the source catalog table data
Pretty-print the source information
In this tutorial we will show examples using the following catalogs:
All catalog and source classes work the same, as long as some
information is available. E.g. trying to access a lightcurve from a
catalog and source that doesn’t have that information will return
None
.
Further information is available at catalog
.
import numpy as np
import astropy.units as u
# %matplotlib inline
import matplotlib.pyplot as plt
from IPython.display import display
from gammapy.catalog import CATALOG_REGISTRY
Check setup#
from gammapy.utils.check import check_tutorials_setup
check_tutorials_setup()
System:
python_executable : /home/runner/work/gammapy-docs/gammapy-docs/gammapy/.tox/build_docs/bin/python
python_version : 3.9.18
machine : x86_64
system : Linux
Gammapy package:
version : 1.0.2
path : /home/runner/work/gammapy-docs/gammapy-docs/gammapy/.tox/build_docs/lib/python3.9/site-packages/gammapy
Other packages:
numpy : 1.26.2
scipy : 1.11.4
astropy : 5.2.2
regions : 0.8
click : 8.1.7
yaml : 6.0.1
IPython : 8.18.1
jupyterlab : not installed
matplotlib : 3.8.2
pandas : not installed
healpy : 1.16.6
iminuit : 2.24.0
sherpa : 4.16.0
naima : 0.10.0
emcee : 3.1.4
corner : 2.2.2
Gammapy environment variables:
GAMMAPY_DATA : /home/runner/work/gammapy-docs/gammapy-docs/gammapy-datasets/1.0.2
List available catalogs#
catalog
contains a catalog registry CATALOG_REGISTRY
,
which maps catalog names (e.g. “3fhl”) to catalog classes
(e.g. SourceCatalog3FHL
).
print(CATALOG_REGISTRY)
Registry
--------
SourceCatalogGammaCat: gamma-cat
SourceCatalogHGPS : hgps
SourceCatalog2HWC : 2hwc
SourceCatalog3FGL : 3fgl
SourceCatalog4FGL : 4fgl
SourceCatalog2FHL : 2fhl
SourceCatalog3FHL : 3fhl
SourceCatalog3HWC : 3hwc
Load catalogs#
If you have run gammapy download datasets
or
gammapy download tutorials
, you have a copy of the catalogs as FITS
files in $GAMMAPY_DATA/catalogs
, and that is the default location
where catalog
loads from.
# # # !ls -1 $GAMMAPY_DATA/catalogs
# # # !ls -1 $GAMMAPY_DATA/catalogs/fermi
So a catalog can be loaded directly from its corresponding class
from gammapy.catalog import SourceCatalog4FGL
catalog = SourceCatalog4FGL()
print("Number of sources :", len(catalog.table))
Number of sources : 6659
Note that it loads the default catalog from $GAMMAPY_DATA/catalogs
,
you could pass a different filename
when creating the catalog. For
example here we load an older version of 4FGL catalog:
catalog = SourceCatalog4FGL("$GAMMAPY_DATA/catalogs/fermi/gll_psc_v20.fit.gz")
print("Number of sources :", len(catalog.table))
/home/runner/work/gammapy-docs/gammapy-docs/gammapy/.tox/build_docs/lib/python3.9/site-packages/astropy/units/core.py:2097: UnitsWarning: '_' did not parse as fits unit: At col 0, Unit '_' not supported by the FITS standard. If this is meant to be a custom unit, define it with 'u.def_unit'. To have it recognized inside a file reader or other code, enable it with 'u.add_enabled_units'. For details, see https://docs.astropy.org/en/latest/units/combining_and_defining.html
warnings.warn(msg, UnitsWarning)
Number of sources : 5066
Alternatively you can load a catalog by name via
CATALOG_REGISTRY.get_cls(name)()
(note the ()
to instantiate a
catalog object from the catalog class - only this will load the catalog
and be useful), or by importing the catalog class
(e.g. SourceCatalog3FGL
) directly. The two ways are equivalent, the
result will be the same.
# FITS file is loaded
catalog = CATALOG_REGISTRY.get_cls("3fgl")()
print(catalog)
SourceCatalog3FGL:
name: 3fgl
description: LAT 4-year point source catalog
sources: 3034
Let’s load the source catalogs we will use throughout this tutorial
catalog_gammacat = CATALOG_REGISTRY.get_cls("gamma-cat")()
catalog_3fhl = CATALOG_REGISTRY.get_cls("3fhl")()
catalog_4fgl = CATALOG_REGISTRY.get_cls("4fgl")()
catalog_hgps = CATALOG_REGISTRY.get_cls("hgps")()
Catalog table#
Source catalogs are given as FITS
files that contain one or multiple
tables.
However, you can also access the underlying astropy.table.Table
for
a catalog, and the row data as a Python dict
. This can be useful if
you want to do something that is not pre-scripted by the
SourceCatalog
classes, such as e.g. selecting sources by sky
position or association class, or accessing special source information.
print(type(catalog_3fhl.table))
<class 'astropy.table.table.Table'>
print(len(catalog_3fhl.table))
1556
display(catalog_3fhl.table[:3][["Source_Name", "RAJ2000", "DEJ2000"]])
Source_Name RAJ2000 DEJ2000
deg deg
------------------ -------- --------
3FHL J0001.2-0748 0.3107 -7.8075
3FHL J0001.9-4155 0.4849 -41.9303
3FHL J0002.1-6728 0.5283 -67.4825
Note that the catalogs object include a helper property that gives
directly the sources positions as a SkyCoord
object (we will show an
usage example in the following).
print(catalog_3fhl.positions[:3])
<SkyCoord (ICRS): (ra, dec) in deg
[(0.31067517, -7.8075185), (0.4848653 , -41.93026 ),
(0.52826166, -67.48248 )]>
Source object#
Select a source#
The catalog entries for a single source are represented by a
SourceCatalogObject
. In order to select a source object index into
the catalog using []
, with a catalog table row index (zero-based,
first row is [0]
), or a source name. If a name is given, catalog
table columns with source names and association names (“ASSOC1” in the
example below) are searched top to bottom. There is no name resolution
web query.
source = catalog_4fgl[49]
print(source)
*** Basic info ***
Catalog row index (zero-based) : 49
Source name : 4FGL J0009.3+5030
Extended name :
Associations : NVSS J000922+503028, 3FGL J0009.3+5030, 3FHL J0009.4+5030
ASSOC_PROB_BAY : 1.000
ASSOC_PROB_LR : 0.928
Class1 : bll
Class2 :
TeVCat flag : N
*** Other info ***
Significance (100 MeV - 1 TeV) : 61.802
Npred : 3263.4
Other flags : 0
*** Position info ***
RA : 2.347 deg
DEC : 50.512 deg
GLON : 116.124 deg
GLAT : -11.804 deg
Semimajor (68%) : 0.0084 deg
Semiminor (68%) : 0.0083 deg
Position angle (68%) : 36.45 deg
Semimajor (95%) : 0.0137 deg
Semiminor (95%) : 0.0135 deg
Position angle (95%) : 36.45 deg
ROI number : 1272
*** Spectral info ***
Spectrum type : LogParabola
Detection significance (100 MeV - 1 TeV) : 61.802
beta : 0.0695 +- 0.01443
Significance curvature : 5.3
Pivot energy : 1379 MeV
Spectral index : 1.884 +- 0.033
Flux Density at pivot energy : 1.33e-12 +- 4.38e-14 cm-2 MeV-1 s-1
Integral flux (1 - 100 GeV) : 2.49e-09 +- 7.72e-11 cm-2 s-1
Energy flux (100 MeV - 100 GeV) : 2.29e-11 +- 9.18e-13 erg cm-2 s-1
*** Spectral points ***
e_min e_max flux flux_errn flux_errp e2dnde e2dnde_errn e2dnde_errp is_ul flux_ul e2dnde_ul sqrt_ts
MeV MeV 1 / (cm2 s) 1 / (cm2 s) 1 / (cm2 s) erg / (cm2 s) erg / (cm2 s) erg / (cm2 s) 1 / (cm2 s) erg / (cm2 s)
---------- ----------- ----------- ----------- ----------- ------------- ------------- ------------- ----- ----------- ------------- -------
50.000 100.000 1.968e-08 1.752e-08 1.359e-08 3.220e-12 2.866e-12 2.224e-12 False nan nan 1.121
100.000 300.000 1.225e-08 1.876e-09 1.903e-09 3.065e-12 4.693e-13 4.761e-13 False nan nan 6.894
300.000 1000.000 4.892e-09 2.759e-10 2.759e-10 3.458e-12 1.950e-13 1.950e-13 False nan nan 22.124
1000.000 3000.000 1.614e-09 7.652e-11 7.652e-11 3.911e-12 1.854e-13 1.854e-13 False nan nan 33.513
3000.000 10000.000 6.074e-10 3.607e-11 3.607e-11 4.133e-12 2.455e-13 2.455e-13 False nan nan 36.197
10000.000 30000.000 1.562e-10 1.598e-11 1.710e-11 3.668e-12 3.753e-13 4.015e-13 False nan nan 25.276
30000.000 100000.000 4.397e-11 8.076e-12 9.104e-12 2.884e-12 5.296e-13 5.970e-13 False nan nan 15.104
100000.000 1000000.000 1.139e-17 nan 2.399e-12 1.602e-18 nan 3.374e-13 True 4.798e-12 6.749e-13 0.000
*** Lightcurve info ***
Lightcurve measured in the energy band: 100 MeV - 100 GeV
Variability index : 204.735
Significance peak (100 MeV - 100 GeV) : 30.253
Integral flux peak (100 MeV - 100 GeV) : 3.3e-08 +- 2.01e-09 cm^-2 s^-1
Time peak : 3.18e+08 s (Mission elapsed time)
Peak interval : 3.65e+02 day
print(source.row_index, source.name)
49 4FGL J0009.3+5030
source = catalog_4fgl["4FGL J0010.8-2154"]
print(source.row_index, source.name)
60 4FGL J0010.8-2154
print(source.data["ASSOC1"])
PKS 0008-222
source = catalog_4fgl["PKS 0008-222"]
print(source.row_index, source.name)
60 4FGL J0010.8-2154
Note that you can also do a for source in catalog
loop, to find or
process sources of interest.
Source information#
The source objects have a data
property that contains the
information of the catalog row corresponding to the source.
print(source.data["Npred"])
333.74478
60.28118133544922 deg -79.40050506591797 deg
As for the catalog object, the source object has a position
property.
print(source.position.galactic)
<SkyCoord (Galactic): (l, b) in deg
(60.28120079, -79.40051035)>
Select a catalog subset#
The catalog objects support selection using boolean arrays (of the same length), so one can create a new catalog as a subset of the main catalog that verify a set of conditions.
In the next example we selection only few of the brightest sources brightest sources in the 100 to 200 GeV energy band.
mask_bright = np.zeros(len(catalog_3fhl.table), dtype=bool)
for k, source in enumerate(catalog_3fhl):
flux = source.spectral_model().integral(100 * u.GeV, 200 * u.GeV).to("cm-2 s-1")
if flux > 1e-10 * u.Unit("cm-2 s-1"):
mask_bright[k] = True
print(f"{source.row_index:<7d} {source.name:20s} {flux:.3g}")
352 3FHL J0534.5+2201 2.99e-10 1 / (cm2 s)
553 3FHL J0851.9-4620e 1.24e-10 1 / (cm2 s)
654 3FHL J1036.3-5833e 1.57e-10 1 / (cm2 s)
691 3FHL J1104.4+3812 3.34e-10 1 / (cm2 s)
1111 3FHL J1653.8+3945 1.27e-10 1 / (cm2 s)
1219 3FHL J1824.5-1351e 1.77e-10 1 / (cm2 s)
1361 3FHL J2028.6+4110e 1.75e-10 1 / (cm2 s)
SourceCatalog3FHL:
name: 3fhl
description: LAT third high-energy source catalog
sources: 7
print(catalog_3fhl_bright.table["Source_Name"])
Source_Name
------------------
3FHL J0534.5+2201
3FHL J0851.9-4620e
3FHL J1036.3-5833e
3FHL J1104.4+3812
3FHL J1653.8+3945
3FHL J1824.5-1351e
3FHL J2028.6+4110e
Similarly we can select only sources within a region of interest. Here
for example we use the position
property of the catalog object to
select sources within 5 degrees from “PKS 0008-222”:
source = catalog_4fgl["PKS 0008-222"]
mask_roi = source.position.separation(catalog_4fgl.positions) < 5 * u.deg
catalog_4fgl_roi = catalog_4fgl[mask_roi]
print("Number of sources :", len(catalog_4fgl_roi.table))
Number of sources : 15
Source models#
The SourceCatalogObject
classes have a
sky_model()
model which creates a
SkyModel
object, with model parameter values
and parameter errors from the catalog filled in.
In most cases, the spectral_model()
method provides the
SpectralModel
part of the sky model, and the
spatial_model()
method the SpatialModel
part individually.
We use the SourceCatalog3FHL
for the examples in
this section.
source = catalog_4fgl["PKS 2155-304"]
model = source.sky_model()
print(model)
SkyModel
Name : 4FGL J2158.8-3013
Datasets names : None
Spectral model type : LogParabolaSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
amplitude : 1.34e-11 +/- 1.5e-13 1 / (cm2 MeV s)
reference (frozen): 1146.894 MeV
alpha : 1.767 +/- 0.01
beta : 0.041 +/- 0.00
lon_0 : 329.714 +/- 0.00 deg
lat_0 : -30.225 +/- 0.00 deg
print(model)
SkyModel
Name : 4FGL J2158.8-3013
Datasets names : None
Spectral model type : LogParabolaSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
amplitude : 1.34e-11 +/- 1.5e-13 1 / (cm2 MeV s)
reference (frozen): 1146.894 MeV
alpha : 1.767 +/- 0.01
beta : 0.041 +/- 0.00
lon_0 : 329.714 +/- 0.00 deg
lat_0 : -30.225 +/- 0.00 deg
print(model.spatial_model)
PointSpatialModel
type name value unit error ... max frozen is_norm link
------- ----- ----------- ---- --------- ... --------- ------ ------- ----
spatial lon_0 3.2971e+02 deg 3.735e-03 ... nan False False
spatial lat_0 -3.0225e+01 deg 3.227e-03 ... 9.000e+01 False False
print(model.spectral_model)
LogParabolaSpectralModel
type name value unit ... max frozen is_norm link
-------- --------- ---------- -------------- ... --- ------ ------- ----
spectral amplitude 1.3413e-11 cm-2 MeV-1 s-1 ... nan False True
spectral reference 1.1469e+03 MeV ... nan True False
spectral alpha 1.7672e+00 ... nan False False
spectral beta 4.0590e-02 ... nan False False
plt.figure()
energy_bounds = (100 * u.MeV, 100 * u.GeV)
opts = dict(sed_type="e2dnde", yunits=u.Unit("TeV cm-2 s-1"))
model.spectral_model.plot(energy_bounds, **opts)
model.spectral_model.plot_error(energy_bounds, **opts)
<Axes: xlabel='Energy [GeV]', ylabel='e2dnde [TeV / (cm2 s)]'>
You can create initial source models for your analyses using the
to_models()
method of the catalog objects. Here for example we
create a Models
object from the 4FGL catalog subset we previously
defined:
Models
Component 0: SkyModel
Name : 4FGL J0001.8-2153
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 1.877 +/- 0.17
amplitude : 4.94e-15 +/- 1.3e-15 1 / (cm2 MeV s)
reference (frozen): 4429.934 MeV
lon_0 : 0.465 +/- 0.04 deg
lat_0 : -21.886 +/- 0.05 deg
Component 1: SkyModel
Name : 4FGL J0003.3-1928
Datasets names : None
Spectral model type : LogParabolaSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
amplitude : 5.61e-13 +/- 5.6e-14 1 / (cm2 MeV s)
reference (frozen): 1021.799 MeV
alpha : 2.103 +/- 0.14
beta : 0.234 +/- 0.09
lon_0 : 0.846 +/- 0.03 deg
lat_0 : -19.468 +/- 0.03 deg
Component 2: SkyModel
Name : 4FGL J0006.3-1813
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.182 +/- 0.22
amplitude : 1.83e-14 +/- 4.9e-15 1 / (cm2 MeV s)
reference (frozen): 2248.979 MeV
lon_0 : 1.578 +/- 0.03 deg
lat_0 : -18.229 +/- 0.03 deg
Component 3: SkyModel
Name : 4FGL J0008.4-2339
Datasets names : None
Spectral model type : LogParabolaSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
amplitude : 1.30e-14 +/- 2.4e-15 1 / (cm2 MeV s)
reference (frozen): 4584.337 MeV
alpha : 1.463 +/- 0.21
beta : 0.199 +/- 0.11
lon_0 : 2.111 +/- 0.02 deg
lat_0 : -23.651 +/- 0.02 deg
Component 4: SkyModel
Name : 4FGL J0010.2-2431
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.139 +/- 0.18
amplitude : 1.78e-14 +/- 4.5e-15 1 / (cm2 MeV s)
reference (frozen): 2286.647 MeV
lon_0 : 2.560 +/- 0.03 deg
lat_0 : -24.524 +/- 0.04 deg
Component 5: SkyModel
Name : 4FGL J0010.8-2154
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.379 +/- 0.15
amplitude : 1.08e-13 +/- 2.0e-14 1 / (cm2 MeV s)
reference (frozen): 1205.260 MeV
lon_0 : 2.717 +/- 0.05 deg
lat_0 : -21.900 +/- 0.05 deg
Component 6: SkyModel
Name : 4FGL J0013.9-1854
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 1.934 +/- 0.09
amplitude : 2.99e-14 +/- 3.6e-15 1 / (cm2 MeV s)
reference (frozen): 3038.329 MeV
lon_0 : 3.480 +/- 0.02 deg
lat_0 : -18.912 +/- 0.02 deg
Component 7: SkyModel
Name : 4FGL J0021.5-2552
Datasets names : None
Spectral model type : LogParabolaSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
amplitude : 9.49e-13 +/- 6.1e-14 1 / (cm2 MeV s)
reference (frozen): 935.760 MeV
alpha : 2.035 +/- 0.06
beta : 0.055 +/- 0.03
lon_0 : 5.391 +/- 0.01 deg
lat_0 : -25.868 +/- 0.01 deg
Component 8: SkyModel
Name : 4FGL J0022.1-1854
Datasets names : None
Spectral model type : LogParabolaSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
amplitude : 3.97e-13 +/- 2.7e-14 1 / (cm2 MeV s)
reference (frozen): 1451.120 MeV
alpha : 1.804 +/- 0.07
beta : 0.087 +/- 0.03
lon_0 : 5.535 +/- 0.01 deg
lat_0 : -18.907 +/- 0.01 deg
Component 9: SkyModel
Name : 4FGL J0025.0-2320
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.330 +/- 0.22
amplitude : 2.48e-14 +/- 6.4e-15 1 / (cm2 MeV s)
reference (frozen): 1975.588 MeV
lon_0 : 6.268 +/- 0.04 deg
lat_0 : -23.338 +/- 0.05 deg
Component 10: SkyModel
Name : 4FGL J0025.2-2231
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.395 +/- 0.20
amplitude : 3.90e-14 +/- 9.1e-15 1 / (cm2 MeV s)
reference (frozen): 1708.668 MeV
lon_0 : 6.303 +/- 0.03 deg
lat_0 : -22.533 +/- 0.03 deg
Component 11: SkyModel
Name : 4FGL J0031.0-2327
Datasets names : None
Spectral model type : LogParabolaSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
amplitude : 1.23e-13 +/- 2.3e-14 1 / (cm2 MeV s)
reference (frozen): 1515.482 MeV
alpha : 1.615 +/- 0.40
beta : 0.659 +/- 0.36
lon_0 : 7.756 +/- 0.04 deg
lat_0 : -23.458 +/- 0.05 deg
Component 12: SkyModel
Name : 4FGL J2357.7-1937
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.259 +/- 0.17
amplitude : 3.55e-14 +/- 8.2e-15 1 / (cm2 MeV s)
reference (frozen): 1816.197 MeV
lon_0 : 359.430 +/- 0.04 deg
lat_0 : -19.619 +/- 0.04 deg
Component 13: SkyModel
Name : 4FGL J2358.5-1808
Datasets names : None
Spectral model type : LogParabolaSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
amplitude : 1.28e-13 +/- 1.1e-14 1 / (cm2 MeV s)
reference (frozen): 2145.829 MeV
alpha : 1.732 +/- 0.10
beta : 0.086 +/- 0.04
lon_0 : 359.639 +/- 0.01 deg
lat_0 : -18.141 +/- 0.01 deg
Component 14: SkyModel
Name : 4FGL J2359.3-2049
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 1.984 +/- 0.08
amplitude : 5.20e-14 +/- 5.2e-15 1 / (cm2 MeV s)
reference (frozen): 2594.250 MeV
lon_0 : 359.836 +/- 0.02 deg
lat_0 : -20.819 +/- 0.02 deg
Specificities of the HGPS catalog#
Using the to_models()
method for the
SourceCatalogHGPS
will return only the models
components of the sources retained in the main catalog, several
candidate objects appears only in the Gaussian components table (see
section 4.9 of the HGPS paper, https://arxiv.org/abs/1804.02432). To
access these components you can do the following:
discarded_ind = np.where(
["Discarded" in _ for _ in catalog_hgps.table_components["Component_Class"]]
)[0]
discarded_table = catalog_hgps.table_components[discarded_ind]
There is no spectral model available for these components but you can access their spatial models:
discarded_spatial = [
catalog_hgps.gaussian_component(idx).spatial_model() for idx in discarded_ind
]
In addition to the source components the HGPS catalog include a large
scale diffuse component built by fitting a gaussian model in a sliding
window along the Galactic plane. Information on this model can be
accessed via the properties table_large_scale_component
and
large_scale_component
of SourceCatalogHGPS
.
# here we show the 5 first elements of the table
display(catalog_hgps.table_large_scale_component[:5])
# you can also try :
# help(catalog_hgps.large_scale_component)
GLON GLAT GLAT_Err ... Surface_Brightness_Err Width Width_Err
deg deg deg ... 1 / (cm2 s sr) deg deg
---------- --------- -------- ... ---------------------- -------- ---------
270.000000 0.205357 0.251932 ... 4.064108e-10 0.269385 0.137990
272.959198 -0.120154 0.058201 ... 7.346488e-10 0.088742 0.041882
275.918365 -0.095232 0.089881 ... 6.117877e-10 0.167219 0.111797
278.877563 -0.257642 0.065071 ... 5.230542e-10 0.156525 0.056130
281.836731 -0.283487 0.066442 ... 4.336444e-10 0.205192 0.049676
Flux points#
The flux points are available via the flux_points
property as a
FluxPoints
object.
source = catalog_4fgl["PKS 2155-304"]
flux_points = source.flux_points
print(flux_points)
FluxPoints
----------
geom : RegionGeom
axes : ['lon', 'lat', 'energy']
shape : (1, 1, 8)
quantities : ['norm', 'norm_errp', 'norm_errn', 'norm_ul', 'sqrt_ts', 'is_ul']
ref. model : lp
n_sigma : 1
n_sigma_ul : 2
sqrt_ts_threshold_ul : 1
sed type init : flux
display(flux_points.to_table(sed_type="flux"))
e_ref e_min e_max ... sqrt_ts is_ul
MeV MeV MeV ...
------------------ ------------------ ------------------ ... --------- -----
70.71067811865478 49.99999999999999 100.00000000000004 ... 2.961729 False
173.20508075688775 100.00000000000004 299.99999999999994 ... 25.350155 False
547.722557505166 299.99999999999994 999.9999999999998 ... 93.123116 False
1732.0508075688763 999.9999999999998 2999.9999999999977 ... 134.43538 False
5477.225575051666 2999.9999999999977 10000.00000000001 ... 129.33427 False
17320.50807568877 10000.00000000001 30000.000000000007 ... 91.40714 False
54772.255750516626 30000.000000000007 100000.00000000001 ... 59.852314 False
316227.7660168382 100000.00000000001 999999.9999999995 ... 29.674353 False
plt.figure()
flux_points.plot(sed_type="e2dnde")
<Axes: xlabel='Energy [MeV]', ylabel='e2dnde [erg / (cm2 s)]'>
Lightcurves#
The Fermi catalogs contain lightcurves for each source. It is available
via the source.lightcurve
method as a
FluxPoints
object with a time axis.
lightcurve = catalog_4fgl["4FGL J0349.8-2103"].lightcurve()
print(lightcurve)
FluxPoints
----------
geom : RegionGeom
axes : ['lon', 'lat', 'energy', 'time']
shape : (1, 1, 1, 12)
quantities : ['norm', 'norm_errp', 'norm_errn', 'norm_ul', 'ts']
ref. model : lp
n_sigma : 1
n_sigma_ul : 2
sqrt_ts_threshold_ul : 1
sed type init : flux
display(lightcurve.to_table(format="lightcurve", sed_type="flux"))
time_min time_max e_ref ... sqrt_ts is_ul
MeV ...
------------------ ------------------ ------------------ ... --------- -----
54682.65603794185 55045.301668796295 3872.9833462074166 ... 30.332916 False
55045.301668796295 55410.57944657407 3872.9833462074166 ... 13.702817 False
55410.57944657407 55775.85722435185 3872.9833462074166 ... 6.40809 False
55775.85722435185 56141.13500212963 3872.9833462074166 ... 8.410365 False
56141.13500212963 56506.412779907405 3872.9833462074166 ... 8.440284 False
56506.412779907405 56871.690557685186 3872.9833462074166 ... 10.094164 False
56871.690557685186 57236.96833546296 3872.9833462074166 ... 13.103473 False
57236.96833546296 57602.24611324074 3872.9833462074166 ... 5.9435167 False
57602.24611324074 57967.523891018514 3872.9833462074166 ... 7.861615 False
57967.523891018514 58332.801668796295 3872.9833462074166 ... 9.68207 False
58332.801668796295 58698.07944657407 3872.9833462074166 ... 4.562178 False
58698.07944657407 59063.35722435185 3872.9833462074166 ... 5.4070954 False
<Axes: xlabel='Time [iso]', ylabel='flux [1 / (cm2 s)]'>
Pretty-print source information#
A source object has a nice string representation that you can print.
source = catalog_hgps["MSH 15-52"]
print(source)
*** Basic info ***
Catalog row index (zero-based) : 18
Source name : HESS J1514-591
Analysis reference : HGPS
Source class : PWN
Identified object : MSH 15-52
Gamma-Cat id : 79
*** Info from map analysis ***
RA : 228.499 deg = 15h14m00s
DEC : -59.161 deg = -59d09m41s
GLON : 320.315 +/- 0.008 deg
GLAT : -1.188 +/- 0.007 deg
Position Error (68%) : 0.020 deg
Position Error (95%) : 0.033 deg
ROI number : 13
Spatial model : 3-Gaussian
Spatial components : HGPSC 023, HGPSC 024, HGPSC 025
TS : 1763.4
sqrt(TS) : 42.0
Size : 0.145 +/- 0.026 (UL: 0.000) deg
R70 : 0.215 deg
RSpec : 0.215 deg
Total model excess : 3502.8
Excess in RSpec : 2440.5
Model Excess in RSpec : 2414.3
Background in RSpec : 1052.5
Livetime : 41.4 hours
Energy threshold : 0.61 TeV
Source flux (>1 TeV) : (6.434 +/- 0.211) x 10^-12 cm^-2 s^-1 = (28.47 +/- 0.94) % Crab
Fluxes in RSpec (> 1 TeV):
Map measurement : 4.552 x 10^-12 cm^-2 s^-1 = 20.14 % Crab
Source model : 4.505 x 10^-12 cm^-2 s^-1 = 19.94 % Crab
Other component model : 0.000 x 10^-12 cm^-2 s^-1 = 0.00 % Crab
Large scale component model : 0.000 x 10^-12 cm^-2 s^-1 = 0.00 % Crab
Total model : 4.505 x 10^-12 cm^-2 s^-1 = 19.94 % Crab
Containment in RSpec : 70.0 %
Contamination in RSpec : 0.0 %
Flux correction (RSpec -> Total) : 142.8 %
Flux correction (Total -> RSpec) : 70.0 %
*** Info from spectral analysis ***
Livetime : 13.7 hours
Energy range: : 0.38 to 61.90 TeV
Background : 1825.9
Excess : 2061.1
Spectral model : ECPL
TS ECPL over PL : 10.2
Best-fit model flux(> 1 TeV) : (5.720 +/- 0.417) x 10^-12 cm^-2 s^-1 = (25.31 +/- 1.84) % Crab
Best-fit model energy flux(1 to 10 TeV) : (20.779 +/- 1.878) x 10^-12 erg cm^-2 s^-1
Pivot energy : 1.54 TeV
Flux at pivot energy : (2.579 +/- 0.083) x 10^-12 cm^-2 s^-1 TeV^-1 = (11.41 +/- 0.24) % Crab
PL Flux(> 1 TeV) : (5.437 +/- 0.186) x 10^-12 cm^-2 s^-1 = (24.06 +/- 0.82) % Crab
PL Flux(@ 1 TeV) : (6.868 +/- 0.241) x 10^-12 cm^-2 s^-1 TeV^-1 = (30.39 +/- 0.69) % Crab
PL Index : 2.26 +/- 0.03
ECPL Flux(@ 1 TeV) : (6.860 +/- 0.252) x 10^-12 cm^-2 s^-1 TeV^-1 = (30.35 +/- 0.72) % Crab
ECPL Flux(> 1 TeV) : (5.720 +/- 0.417) x 10^-12 cm^-2 s^-1 = (25.31 +/- 1.84) % Crab
ECPL Index : 2.05 +/- 0.06
ECPL Lambda : 0.052 +/- 0.014 TeV^-1
ECPL E_cut : 19.20 +/- 5.01 TeV
*** Flux points info ***
Number of flux points: 6
Flux points table:
e_ref e_min e_max dnde dnde_errn dnde_errp dnde_ul is_ul
TeV TeV TeV 1 / (cm2 s TeV) 1 / (cm2 s TeV) 1 / (cm2 s TeV) 1 / (cm2 s TeV)
------ ------ ------ --------------- --------------- --------------- --------------- -----
0.562 0.383 0.825 2.439e-11 1.419e-12 1.509e-12 2.732e-11 False
1.212 0.825 1.778 4.439e-12 2.489e-13 2.654e-13 4.970e-12 False
2.738 1.778 4.217 7.295e-13 4.788e-14 4.898e-14 8.302e-13 False
6.190 4.217 9.085 1.305e-13 1.220e-14 1.282e-14 1.571e-13 False
13.991 9.085 21.544 1.994e-14 2.723e-15 2.858e-15 2.588e-14 False
31.623 21.544 46.416 9.474e-16 3.480e-16 4.329e-16 1.919e-15 False
*** Gaussian component info ***
Number of components: 3
Spatial components : HGPSC 023, HGPSC 024, HGPSC 025
Component HGPSC 023:
GLON : 320.303 +/- 0.005 deg
GLAT : -1.124 +/- 0.007 deg
Size : 0.057 +/- 0.005 deg
Flux (>1 TeV) : (2.01 +/- 0.23) x 10^-12 cm^-2 s^-1 = (8.9 +/- 1.0) % Crab
Component HGPSC 024:
GLON : 320.298 +/- 0.020 deg
GLAT : -1.168 +/- 0.021 deg
Size : 0.206 +/- 0.020 deg
Flux (>1 TeV) : (2.54 +/- 0.29) x 10^-12 cm^-2 s^-1 = (11.2 +/- 1.3) % Crab
Component HGPSC 025:
GLON : 320.351 +/- 0.005 deg
GLAT : -1.284 +/- 0.007 deg
Size : 0.055 +/- 0.005 deg
Flux (>1 TeV) : (1.88 +/- 0.22) x 10^-12 cm^-2 s^-1 = (8.3 +/- 1.0) % Crab
*** Source associations info ***
Source_Name Association_Catalog Association_Name Separation
deg
---------------- ------------------- --------------------- ----------
HESS J1514-591 2FHL 2FHL J1514.0-5915e 0.098903
HESS J1514-591 3FGL 3FGL J1513.9-5908 0.026914
HESS J1514-591 3FGL 3FGL J1514.0-5915e 0.094834
HESS J1514-591 COMP G320.4-1.2 0.070483
HESS J1514-591 PSR B1509-58 0.026891
*** Source identification info ***
Source_Name: HESS J1514-591
Identified_Object: MSH 15-52
Class: PWN
Evidence: Morphology
Reference: 2005A%26A...435L..17A
Distance_Reference: SNRCat
Distance: 5.199999809265137 kpc
Distance_Min: 3.799999952316284 kpc
Distance_Max: 6.599999904632568 kpc
You can also call source.info()
instead and pass as an option what
information to print. The options available depend on the catalog, you
can learn about them using help()
help(source.info)
Help on method info in module gammapy.catalog.hess:
info(info='all') method of gammapy.catalog.hess.SourceCatalogObjectHGPS instance
Info string.
Parameters
----------
info : {'all', 'basic', 'map', 'spec', 'flux_points', 'components', 'associations', 'id'}
Comma separated list of options
print(source.info("associations"))
*** Source associations info ***
Source_Name Association_Catalog Association_Name Separation
deg
---------------- ------------------- --------------------- ----------
HESS J1514-591 2FHL 2FHL J1514.0-5915e 0.098903
HESS J1514-591 3FGL 3FGL J1513.9-5908 0.026914
HESS J1514-591 3FGL 3FGL J1514.0-5915e 0.094834
HESS J1514-591 COMP G320.4-1.2 0.070483
HESS J1514-591 PSR B1509-58 0.026891
plt.show()