File: hess.py

package info (click to toggle)
gammapy 2.0.1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 10,800 kB
  • sloc: python: 81,999; makefile: 211; sh: 11; javascript: 10
file content (135 lines) | stat: -rw-r--r-- 4,449 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
"""
H.E.S.S. with Gammapy
=====================

Explore H.E.S.S. event lists and IRFs.


Introduction
------------

`H.E.S.S. <https://www.mpi-hd.mpg.de/hfm/HESS/>`__ is an array of
gamma-ray telescopes located in Namibia. Gammapy is regularly used and
fully supports H.E.S.S. high level data analysis, after export to the
current `open data level 3
format <https://gamma-astro-data-formats.readthedocs.io/>`__.

The H.E.S.S. data is private, and H.E.S.S. analysis is mostly documented
and discussed in the internal Wiki pages and in
H.E.S.S.-internal communication channels. However, in 2018, a small
sub-set of archival H.E.S.S. data was publicly released, called the
`H.E.S.S. DL3
DR1 <https://www.mpi-hd.mpg.de/hfm/HESS/pages/dl3-dr1/>`__, the data
level 3, data release number 1. This dataset is 50 MB in size and is
used in many Gammapy analysis tutorials, and can be downloaded via
`gammapy
download <https://docs.gammapy.org/dev/getting-started/index.html#quickstart-setup>`__.

This notebook is a quick introduction to this specific DR1 release. It
briefly describes H.E.S.S. data and instrument responses and show a
simple exploration of the data with the creation of theta-squared plot.

H.E.S.S. members can find details on the DL3 FITS production on this
`Confluence
page <https://cchesswiki.in2p3.fr/en/hess/working_groups/analysis_and_reconstruction_working_group/ar_active_tasks/hess_fits_data>`__
and access more detailed tutorials in the BitBucket `hess-open-tools` repository.

DL3 DR1
-------

This is how to access data and IRFs from the H.E.S.S. data level 3, data
release 1.

"""

import astropy.units as u

# %matplotlib inline
import matplotlib.pyplot as plt
from IPython.display import display
from gammapy.data import DataStore


######################################################################
# A useful way to organize the relevant files are the index tables. The
# observation index table contains information on each particular run,
# such as the pointing, or the run ID. The HDU index table has a row per
# relevant file (i.e., events, effective area, psf…) and contains the path
# to said file. Together they can be loaded into a Datastore by indicating
# the directory in which they can be found, in this case
# `$GAMMAPY_DATA/hess-dl3-dr1`:
#

######################################################################
# Create and get info on the data store

data_store = DataStore.from_dir("$GAMMAPY_DATA/hess-dl3-dr1")

data_store.info()

######################################################################
# Preview an excerpt from the observation table

display(data_store.obs_table[:2][["OBS_ID", "DATE-OBS", "RA_PNT", "DEC_PNT", "OBJECT"]])

######################################################################
# Get a single observation

obs = data_store.obs(23523)

######################################################################
# Select and peek events

obs.events.select_offset([0, 2.5] * u.deg).peek()
plt.show()

######################################################################
# Peek the effective area

obs.aeff.peek()
plt.show()

######################################################################
# Peek the energy dispersion

obs.edisp.peek()
plt.show()

######################################################################
# Peek the psf
obs.psf.peek()
plt.show()

######################################################################
# Peek the background rate
obs.bkg.to_2d().plot()
plt.show()


######################################################################
# Exercises
# ---------
#
# -  Find the `OBS_ID` for the runs of the Crab nebula
# -  Compute the expected number of background events in the whole RoI for
#    `OBS_ID=23523` in the 1 TeV to 3 TeV energy band, from the
#    background IRF.
#


######################################################################
# Next steps
# ----------
#
# Now you know how to access and work with H.E.S.S. data. All other
# tutorials and documentation apply to H.E.S.S. and CTAO or any other IACT
# that provides DL3 data and IRFs in the standard format.
#
# You can see the following tutorials for more detailed analysis using H.E.S.S. data
#
# -  :doc:`/tutorials/starting/analysis_1`
# -  :doc:`/tutorials/starting/analysis_2`
# -  :doc:`/tutorials/analysis-1d/spectral_analysis`
# -  :doc:`/tutorials/analysis-1d/extended_source_spectral_analysis`
# -  :doc:`/tutorials/analysis-2d/ring_background`
#