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

Spectral analysis with Gammapy

Prerequisites

  • Understanding how spectral extraction is performed in Cherenkov astronomy, in particular regarding OFF background measurements.

  • Understanding the basics data reduction and modeling/fitting process with the gammapy library API as shown in the first gammapy analysis with the gammapy library API tutorial

Context

While 3D analysis allows in principle to deal with complex situations such as overlapping sources, in many cases, it is not required to extract the spectrum of a source. Spectral analysis, where all data inside a ON region are binned into 1D datasets, provides a nice alternative.

In classical Cherenkov astronomy, it is used with a specific background estimation technique that relies on OFF measurements taken in the field-of-view in regions where the background rate is assumed to be equal to the one in the ON region.

This allows to use a specific fit statistics for ON-OFF measurements, the wstat (see gammapy.stats.fit_statistics), where no background model is assumed. Background is treated as a set of nuisance parameters. This removes some systematic effects connected to the choice or the quality of the background model. But this comes at the expense of larger statistical uncertainties on the fitted model parameters.

Objective: perform a full region based spectral analysis of 4 Crab observations of H.E.S.S. data release 1 and fit the resulting datasets.

Introduction

Here, as usual, we use the gammapy.data.DataStore to retrieve a list of selected observations (gammapy.data.Observations). Then, we define the ON region containing the source and the geometry of the gammapy.cube.SpectrumDataset object we want to produce. We then create the corresponding dataset Maker.

We have to define the Maker object that will extract the OFF counts from reflected regions in the field-of-view. To ensure we use data in an energy range where the quality of the IRFs is good enough we also create a safe range Maker.

We can then proceed with data reduction with a loop over all selected observations to produce datasets in the relevant geometry.

We can then explore the resulting datasets and look at the cumulative signal and significance of our source. We finally proceed with model fitting.

In practice, we have to: - Create a gammapy.data.DataStore poiting to the relevant data - Apply an observation selection to produce a list of observations, a gammapy.data.Observations object. - Define a geometry of the spectrum we want to produce: - Create a ~regions.CircleSkyRegion for the ON extraction region - Create a gammapy.maps.MapAxis for the energy binnings: one for the reconstructed (i.e. measured) energy, the other for the true energy (i.e. the one used by IRFs and models) - Create the necessary makers : - the spectrum dataset maker : gammapy.cube.SpectrumDatasetMaker - the OFF background maker, here a gammapy.cube.ReflectedRegionsBackgroundMaker - and the safe range maker : gammapy.cube.SafeRangeMaker - Perform the data reduction loop. And for every observation: - Apply the makers sequentially to produce a gammapy.maps.SpectrumDatasetOnOff - Append it to list of datasets - Define the gammapy.modeling.models.SkyModel to apply to the dataset. - Create a gammapy.modeling.Fit object and run it to fit the model parameters - Apply a gammapy.spectrum.FluxPointsEstimator to compute flux points for the spectral part of the fit.

Setup

As usual, we’ll start with some setup …

[1]:
%matplotlib inline
import matplotlib.pyplot as plt
[2]:
# Check package versions
import gammapy
import numpy as np
import astropy
import regions

print("gammapy:", gammapy.__version__)
print("numpy:", np.__version__)
print("astropy", astropy.__version__)
print("regions", regions.__version__)
gammapy: 0.16
numpy: 1.18.1
astropy 4.1.dev27293
regions 0.4
[3]:
from pathlib import Path
import astropy.units as u
from astropy.coordinates import SkyCoord, Angle
from regions import CircleSkyRegion
from gammapy.maps import Map
from gammapy.modeling import Fit, Datasets
from gammapy.data import DataStore
from gammapy.modeling.models import (
    PowerLawSpectralModel,
    create_crab_spectral_model,
    SkyModel,
)
from gammapy.cube import SafeMaskMaker
from gammapy.spectrum import (
    SpectrumDatasetMaker,
    SpectrumDatasetOnOff,
    SpectrumDataset,
    FluxPointsEstimator,
    FluxPointsDataset,
    ReflectedRegionsBackgroundMaker,
    plot_spectrum_datasets_off_regions,
)

Load Data

First, we select and load some H.E.S.S. observations of the Crab nebula (simulated events for now).

We will access the events, effective area, energy dispersion, livetime and PSF for containement correction.

[4]:
datastore = DataStore.from_dir("$GAMMAPY_DATA/hess-dl3-dr1/")
obs_ids = [23523, 23526, 23559, 23592]
observations = datastore.get_observations(obs_ids)

Define Target Region

The next step is to define a signal extraction region, also known as on region. In the simplest case this is just a CircleSkyRegion, but here we will use the Target class in gammapy that is useful for book-keeping if you run several analysis in a script.

[5]:
target_position = SkyCoord(ra=83.63, dec=22.01, unit="deg", frame="icrs")
on_region_radius = Angle("0.11 deg")
on_region = CircleSkyRegion(center=target_position, radius=on_region_radius)

Create exclusion mask

We will use the reflected regions method to place off regions to estimate the background level in the on region. To make sure the off regions don’t contain gamma-ray emission, we create an exclusion mask.

Using http://gamma-sky.net/ we find that there’s only one known gamma-ray source near the Crab nebula: the AGN called RGB J0521+212 at GLON = 183.604 deg and GLAT = -8.708 deg.

[6]:
exclusion_region = CircleSkyRegion(
    center=SkyCoord(183.604, -8.708, unit="deg", frame="galactic"),
    radius=0.5 * u.deg,
)

skydir = target_position.galactic
exclusion_mask = Map.create(
    npix=(150, 150), binsz=0.05, skydir=skydir, proj="TAN", frame="icrs"
)

mask = exclusion_mask.geom.region_mask([exclusion_region], inside=False)
exclusion_mask.data = mask
exclusion_mask.plot();
../_images/notebooks_spectrum_analysis_12_0.png

Run data reduction chain

We begin with the configuration of the maker classes:

[7]:
e_reco = np.logspace(-1, np.log10(40), 40) * u.TeV
e_true = np.logspace(np.log10(0.05), 2, 200) * u.TeV
dataset_empty = SpectrumDataset.create(
    e_reco=e_reco, e_true=e_true, region=on_region
)
[8]:
dataset_maker = SpectrumDatasetMaker(
    containment_correction=False, selection=["counts", "aeff", "edisp"]
)
bkg_maker = ReflectedRegionsBackgroundMaker(exclusion_mask=exclusion_mask)
safe_mask_masker = SafeMaskMaker(methods=["aeff-max"], aeff_percent=10)
[9]:
%%time
datasets = []

for obs_id, observation in zip(obs_ids, observations):
    dataset = dataset_maker.run(
        dataset_empty.copy(name=str(obs_id)), observation
    )
    dataset_on_off = bkg_maker.run(dataset, observation)
    dataset_on_off = safe_mask_masker.run(dataset_on_off, observation)
    datasets.append(dataset_on_off)
CPU times: user 2.85 s, sys: 134 ms, total: 2.98 s
Wall time: 3.36 s

Plot off regions

[10]:
plt.figure(figsize=(8, 8))
_, ax, _ = exclusion_mask.plot()
on_region.to_pixel(ax.wcs).plot(ax=ax, edgecolor="k")
plot_spectrum_datasets_off_regions(ax=ax, datasets=datasets)
../_images/notebooks_spectrum_analysis_18_0.png

Source statistic

Next we’re going to look at the overall source statistics in our signal region.

[11]:
datasets_all = Datasets(datasets)
[12]:
info_table = datasets_all.info_table(cumulative=True)
[13]:
info_table
[13]:
Table length=4
namelivetimen_onbackgroundexcesssignificancebackground_rategamma_ratea_onn_offa_offalpha
s1 / s1 / s
str7float64int64float64float64float64float64float64float64int64float64float64
stacked1581.736758410930617213.5158.521.1081393207378380.0085349220900465020.100206307501657091.016212.00.08333333333333333
stacked3154.423482418060336531.999999999999993333.029.9338166900580160.010144484460745270.105566041419630481.038412.00.08333333333333333
stacked4732.54699993133551241.90243902439024470.0975609756098437.371273295857450.0088540988657900710.099332887973945211.079018.8533178114086160.05304106205619018
stacked6313.81164062023263652.54132791327913583.458672086720841.988553624767070.0083216495682658320.092409895210210171.0117322.3252827171796650.0447922659107239
[14]:
plt.plot(
    info_table["livetime"].to("h"), info_table["excess"], marker="o", ls="none"
)
plt.xlabel("Livetime [h]")
plt.ylabel("Excess");
../_images/notebooks_spectrum_analysis_23_0.png
[15]:
plt.plot(
    info_table["livetime"].to("h"),
    info_table["significance"],
    marker="o",
    ls="none",
)
plt.xlabel("Livetime [h]")
plt.ylabel("Significance");
../_images/notebooks_spectrum_analysis_24_0.png
[16]:
datasets[0].peek()
../_images/notebooks_spectrum_analysis_25_0.png

Finally you can write the extrated datasets to disk using the OGIP format (PHA, ARF, RMF, BKG, see here for details):

[17]:
path = Path("spectrum_analysis")
path.mkdir(exist_ok=True)
[18]:
for dataset in datasets:
    dataset.to_ogip_files(outdir=path, overwrite=True)

If you want to read back the datasets from disk you can use:

[19]:
datasets = []
for obs_id in obs_ids:
    filename = path / f"pha_obs{obs_id}.fits"
    datasets.append(SpectrumDatasetOnOff.from_ogip_files(filename))

Fit spectrum

Now we’ll fit a global model to the spectrum. First we do a joint likelihood fit to all observations. If you want to stack the observations see below. We will also produce a debug plot in order to show how the global fit matches one of the individual observations.

[20]:
spectral_model = PowerLawSpectralModel(
    index=2, amplitude=2e-11 * u.Unit("cm-2 s-1 TeV-1"), reference=1 * u.TeV
)
model = SkyModel(spectral_model=spectral_model)

for dataset in datasets:
    dataset.models = model

fit_joint = Fit(datasets)
result_joint = fit_joint.run()

# we make a copy here to compare it later
model_best_joint = model.copy()
model_best_joint.spectral_model.parameters.covariance = (
    result_joint.parameters.covariance
)
[21]:
print(result_joint)
OptimizeResult

        backend    : minuit
        method     : minuit
        success    : True
        message    : Optimization terminated successfully.
        nfev       : 48
        total stat : 122.01

[22]:
plt.figure(figsize=(8, 6))
ax_spectrum, ax_residual = datasets[0].plot_fit()
ax_spectrum.set_ylim(0.1, 40)
[22]:
(0.1, 40)
../_images/notebooks_spectrum_analysis_34_1.png

Compute Flux Points

To round up our analysis we can compute flux points by fitting the norm of the global model in energy bands. We’ll use a fixed energy binning for now:

[23]:
e_min, e_max = 0.7, 30
e_edges = np.logspace(np.log10(e_min), np.log10(e_max), 11) * u.TeV

Now we create an instance of the gammapy.spectrum.FluxPointsEstimator, by passing the dataset and the energy binning:

[24]:
fpe = FluxPointsEstimator(datasets=datasets, e_edges=e_edges)
flux_points = fpe.run()
/Users/adonath/github/adonath/astropy/astropy/units/quantity.py:481: RuntimeWarning: invalid value encountered in true_divide
  result = super().__array_ufunc__(function, method, *arrays, **kwargs)

Here is a the table of the resulting flux points:

[25]:
flux_points.table_formatted
[25]:
Table length=10
e_refe_mine_maxref_dnderef_fluxref_efluxref_e2dndenormstatnorm_errcounts [4]norm_errpnorm_errnnorm_ulsqrt_tstsnorm_scan [11]stat_scan [11]dndednde_uldnde_errdnde_errpdnde_errn
TeVTeVTeV1 / (cm2 s TeV)1 / (cm2 s)TeV / (cm2 s)TeV / (cm2 s)1 / (cm2 s TeV)1 / (cm2 s TeV)1 / (cm2 s TeV)1 / (cm2 s TeV)1 / (cm2 s TeV)
float64float64float64float64float64float64float64float64float64float64int64float64float64float64float64float64float64float64float64float64float64float64float64
0.8590.7371.0024.040e-111.077e-119.177e-122.982e-110.96918.3940.08749 .. 310.0890.0841.15220.687427.9480.200 .. 5.000188.879 .. 679.4413.915e-114.655e-113.502e-123.607e-123.400e-12
1.2611.0021.5881.488e-118.850e-121.095e-112.368e-110.97212.5850.09131 .. 360.0930.0881.16420.513420.7860.200 .. 5.000172.382 .. 615.3051.447e-111.733e-111.349e-121.391e-121.308e-12
1.8521.5882.1605.482e-123.151e-125.786e-121.881e-111.1969.9480.14927 .. 120.1560.1431.51915.286233.6710.200 .. 5.000115.964 .. 242.7346.559e-128.329e-128.190e-138.527e-137.861e-13
2.5182.1602.9362.466e-121.927e-124.811e-121.564e-111.2519.1750.17710 .. 110.1850.1691.63813.945194.4710.200 .. 5.00096.552 .. 175.1893.085e-124.038e-124.361e-134.573e-134.156e-13
3.6972.9364.6569.082e-131.583e-125.741e-121.242e-110.93914.9410.15817 .. 110.1670.1501.28910.703114.5480.200 .. 5.00060.497 .. 213.6798.529e-131.171e-121.436e-131.517e-131.359e-13
5.4294.6566.3303.345e-135.637e-133.034e-129.860e-121.1368.7590.2669 .. 50.2860.2461.7498.14166.2710.200 .. 5.00038.262 .. 80.9483.799e-135.852e-138.896e-149.574e-148.246e-14
7.9716.33010.0371.232e-134.630e-133.620e-127.829e-121.04012.4160.2745 .. 40.2970.2521.6806.85546.9950.200 .. 5.00032.845 .. 80.5281.281e-132.070e-133.378e-143.659e-143.110e-14
11.70310.03713.6474.539e-141.649e-131.913e-126.217e-120.90210.5000.4300 .. 10.4910.3742.0093.36911.3510.200 .. 5.00015.516 .. 38.2724.094e-149.118e-141.951e-142.227e-141.697e-14
17.18313.64721.6361.672e-141.354e-132.282e-124.936e-120.3107.2060.2811 .. 00.3580.3101.1850.9530.9080.200 .. 5.0007.390 .. 42.9775.189e-151.980e-144.705e-155.977e-155.189e-15
25.22921.63629.4196.158e-154.822e-141.206e-123.920e-120.0000.5240.0010 .. 00.2770.0001.1080.0030.0000.200 .. 5.0001.246 .. 18.5691.099e-206.825e-158.662e-181.706e-151.099e-20

Now we plot the flux points and their likelihood profiles. For the plotting of upper limits we choose a threshold of TS < 4.

[26]:
plt.figure(figsize=(8, 5))
flux_points.table["is_ul"] = flux_points.table["ts"] < 4
ax = flux_points.plot(
    energy_power=2, flux_unit="erg-1 cm-2 s-1", color="darkorange"
)
flux_points.to_sed_type("e2dnde").plot_ts_profiles(ax=ax)
[26]:
<matplotlib.axes._subplots.AxesSubplot at 0x1c177bee80>
../_images/notebooks_spectrum_analysis_42_1.png

The final plot with the best fit model, flux points and residuals can be quickly made like this:

[27]:
flux_points_dataset = FluxPointsDataset(
    data=flux_points, models=model_best_joint
)
[28]:
plt.figure(figsize=(8, 6))
flux_points_dataset.peek();
../_images/notebooks_spectrum_analysis_45_0.png

Stack observations

And alternative approach to fitting the spectrum is stacking all observations first and the fitting a model. For this we first stack the individual datasets:

[29]:
dataset_stacked = Datasets(datasets).stack_reduce()

Again we set the model on the dataset we would like to fit (in this case it’s only a single one) and pass it to the gammapy.modeling.Fit object:

[30]:
dataset_stacked.models = model
stacked_fit = Fit([dataset_stacked])
result_stacked = stacked_fit.run()

# make a copy to compare later
model_best_stacked = model.copy()
model_best_stacked.spectral_model.parameters.covariance = (
    result_stacked.parameters.covariance
)
[31]:
print(result_stacked)
OptimizeResult

        backend    : minuit
        method     : minuit
        success    : True
        message    : Optimization terminated successfully.
        nfev       : 26
        total stat : 29.05

[32]:
model_best_joint.parameters.to_table()
[32]:
Table length=3
namevalueerrorunitminmaxfrozen
str9float64float64str14float64float64bool
index2.600e+005.888e-02nannanFalse
amplitude2.723e-111.240e-12cm-2 s-1 TeV-1nannanFalse
reference1.000e+000.000e+00TeVnannanTrue
[33]:
model_best_stacked.parameters.to_table()
[33]:
Table length=3
namevalueerrorunitminmaxfrozen
str9float64float64str14float64float64bool
index2.601e+005.887e-02nannanFalse
amplitude2.723e-111.240e-12cm-2 s-1 TeV-1nannanFalse
reference1.000e+000.000e+00TeVnannanTrue

Finally, we compare the results of our stacked analysis to a previously published Crab Nebula Spectrum for reference. This is available in gammapy.spectrum.

[34]:
plot_kwargs = {
    "energy_range": [0.1, 30] * u.TeV,
    "energy_power": 2,
    "flux_unit": "erg-1 cm-2 s-1",
}

# plot stacked model
model_best_stacked.spectral_model.plot(
    **plot_kwargs, label="Stacked analysis result"
)
model_best_stacked.spectral_model.plot_error(**plot_kwargs)

# plot joint model
model_best_joint.spectral_model.plot(
    **plot_kwargs, label="Joint analysis result", ls="--"
)
model_best_joint.spectral_model.plot_error(**plot_kwargs)

create_crab_spectral_model("hess_pl").plot(
    **plot_kwargs, label="Crab reference"
)
plt.legend()
[34]:
<matplotlib.legend.Legend at 0x1c16a66f98>
../_images/notebooks_spectrum_analysis_54_1.png

Exercises

Now you have learned the basics of a spectral analysis with Gammapy. To practice you can continue with the following exercises:

What next?

The methods shown in this tutorial is valid for point-like or midly extended sources where we can assume that the IRF taken at the region center is valid over the whole region. If one wants to extract the 1D spectrum of a large source and properly average the response over the extraction region, one has to use a different approach explained in the extended source spectral analysis tutorial.

[ ]: