This is a fixed-text formatted version of a Jupyter notebook.

Light curves

Introduction

This tutorial explain how to compute a light curve with Gammapy.

We will use the four Crab nebula observations from the H.E.S.S. first public test data release and compute per-observation fluxes. The Crab nebula is not known to be variable at TeV energies, so we expect constant brightness within statistical and systematic errors.

The main classes we will use are:

Setup

As usual, we’ll start with some setup…

In [1]:
%matplotlib inline
import matplotlib.pyplot as plt
import astropy.units as u
from astropy.coordinates import SkyCoord, Angle
from regions import CircleSkyRegion
from gammapy.utils.energy import EnergyBounds
from gammapy.data import DataStore
from gammapy.spectrum import SpectrumExtraction
from gammapy.spectrum.models import PowerLaw
from gammapy.background import ReflectedRegionsBackgroundEstimator
from gammapy.time import LightCurve, LightCurveEstimator

Spectrum

The LightCurveEstimator is based on a 1d spectral analysis within each time bin. So before we can make the light curve, we have to extract 1d spectra.

In [2]:
data_store = DataStore.from_dir("$GAMMAPY_DATA/hess-dl3-dr1/")
In [3]:
mask = data_store.obs_table["TARGET_NAME"] == "Crab"
obs_ids = data_store.obs_table["OBS_ID"][mask].data
observations = data_store.get_observations(obs_ids)
print(observations)
Observations
Number of observations: 4
Info for OBS_ID = 23523
- Start time: 53343.92
- Pointing pos: RA 83.63 deg / Dec 21.51 deg
- Observation duration: 1687.0 s
- Dead-time fraction: 6.240 %
Info for OBS_ID = 23526
- Start time: 53343.95
- Pointing pos: RA 83.63 deg / Dec 22.51 deg
- Observation duration: 1683.0 s
- Dead-time fraction: 6.555 %
Info for OBS_ID = 23559
- Start time: 53345.96
- Pointing pos: RA 85.25 deg / Dec 22.01 deg
- Observation duration: 1686.0 s
- Dead-time fraction: 6.398 %
Info for OBS_ID = 23592
- Start time: 53347.91
- Pointing pos: RA 82.01 deg / Dec 22.01 deg
- Observation duration: 1686.0 s
- Dead-time fraction: 6.212 %

In [4]:
# Target definition
target_position = SkyCoord(ra=83.63308, dec=22.01450, unit="deg")
on_region_radius = Angle("0.2 deg")
on_region = CircleSkyRegion(center=target_position, radius=on_region_radius)
In [5]:
%%time
bkg_estimator = ReflectedRegionsBackgroundEstimator(
    on_region=on_region, observations=observations
)
bkg_estimator.run()
CPU times: user 1.75 s, sys: 19.5 ms, total: 1.77 s
Wall time: 1.77 s
In [6]:
%%time
ebounds = EnergyBounds.equal_log_spacing(0.7, 100, 50, unit="TeV")
extraction = SpectrumExtraction(
    observations=observations,
    bkg_estimate=bkg_estimator.result,
    containment_correction=False,
    e_reco=ebounds,
    e_true=ebounds,
)
extraction.run()
spectrum_observations = extraction.spectrum_observations
CPU times: user 961 ms, sys: 18.6 ms, total: 979 ms
Wall time: 979 ms

Light curve estimation

OK, so now that we have prepared 1D spectra (not spectral models, just the 1D counts and exposure and 2D energy dispersion matrix), we can compute a lightcurve.

To compute the light curve, a spectral model shape has to be assumed, and an energy band chosen. The method is then to adjust the amplitude parameter of the spectral model in each time bin to the data, resulting in a flux measurement in each time bin.

In [7]:
# Creat list of time bin intervals
# Here we do one time bin per observation
def time_intervals_per_obs(observations):
    for obs in observations:
        yield obs.tstart, obs.tstop


time_intervals = list(time_intervals_per_obs(observations))
In [8]:
# Assumed spectral model
spectral_model = PowerLaw(
    index=2, amplitude=2.0e-11 * u.Unit("1 / (cm2 s TeV)"), reference=1 * u.TeV
)
In [9]:
energy_range = [1, 100] * u.TeV
In [10]:
%%time
lc_estimator = LightCurveEstimator(extraction)
lc = lc_estimator.light_curve(
    time_intervals=time_intervals,
    spectral_model=spectral_model,
    energy_range=energy_range,
)
CPU times: user 1.05 s, sys: 17.3 ms, total: 1.07 s
Wall time: 1.07 s

Results

The light curve measurement result is stored in a table. Let’s have a look at the results:

In [11]:
print(lc.table.colnames)
['time_min', 'time_max', 'flux', 'flux_err', 'flux_ul', 'is_ul', 'livetime', 'alpha', 'n_on', 'n_off', 'measured_excess', 'expected_excess']
In [12]:
lc.table["time_min", "time_max", "flux", "flux_err"]
Out[12]:
Table length=4
time_mintime_maxfluxflux_err
1 / (cm2 s)1 / (cm2 s)
float64float64float64float64
53343.9223400925953343.941865555561.839834552434914e-111.8993229068363982e-12
53343.9542150925953343.973694259261.9975853809663366e-112.0202859435674313e-12
53345.9619812962953345.981495185182.193293401955398e-112.5526997556899865e-12
53347.91319657407453347.932710462962.317987907365536e-112.5670181770343347e-12
In [13]:
lc.plot();
../_images/notebooks_light_curve_17_0.png
In [14]:
# Let's compare to the expected flux of this source
from gammapy.spectrum import CrabSpectrum

crab_spec = CrabSpectrum().model
crab_flux = crab_spec.integral(*energy_range).to("cm-2 s-1")
crab_flux
Out[14]:
$$2.0740149 \times 10^{-11} \; \mathrm{\frac{1}{s\,cm^{2}}}$$
In [15]:
ax = lc.plot(marker="o", lw=2)
ax.hlines(
    crab_flux.value,
    xmin=lc.table["time_min"].min(),
    xmax=lc.table["time_max"].max(),
);
../_images/notebooks_light_curve_19_0.png

Exercises

  • Change the assumed spectral model shape (e.g. to a steeper power-law), and see how the integral flux estimate for the lightcurve changes.
  • Try a time binning where you split the observation time for every run into two time bins.
  • Try to analyse the PKS 2155 flare data from the H.E.S.S. first public test data release. Start with per-observation fluxes, and then try fluxes within 5 minute time bins for one or two of the observations where the source was very bright.
In [16]: