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

Spectral analysis of extended sources

Prerequisites:

Context

Many VHE sources in the Galaxy are extended. Studying them with a 1D spectral analysis is more complex than studying point sources. One often has to use complex (i.e. non circular) regions and more importantly, one has to take into account the fact that the instrument response is non uniform over the selectred region. A typical example is given by the supernova remnant RX J1713-3935 which is nearly 1 degree in diameter. See the following article.

Objective: Measure the spectrum of RX J1713-3945 in a 1 degree region fully enclosing it.

Proposed approach:

We have seen in the general presentation of the spectrum extraction for point sources, see the corresponding notebook, that Gammapy uses specific datasets makers to first produce reduced spectral data and then to extract OFF measurements with reflected background techniques: the gammapy.makers.SpectrumDatasetMaker and the gammapy.makers.ReflectedRegionsBackgroundMaker. The former simply computes the reduced IRF at the center of the ON region (assumed to be circular).

This is no longer valid for extended sources. To be able to compute average responses in the ON region, Gammapy relies on the creation of a cube enclosing it (i.e. a gammapy.datasets.MapDataset) which can be reduced to a simple spectrum (i.e. a gammapy.datasets.SpectrumDataset). We can then proceed with the OFF extraction as the standard point source case.

In summary, we have to:

Here, we will use the RX J1713-3945 observations from the H.E.S.S. first public test data release. The tutorial is implemented with the intermediate level API.

Setup

As usual, we’ll start with some general imports…

[1]:
%matplotlib inline
import matplotlib.pyplot as plt
[2]:
import astropy.units as u
from astropy.coordinates import SkyCoord, Angle
from regions import CircleSkyRegion
from gammapy.maps import Map, MapAxis, WcsGeom
from gammapy.modeling import Fit
from gammapy.data import DataStore
from gammapy.modeling.models import PowerLawSpectralModel, SkyModel
from gammapy.datasets import Datasets, MapDataset
from gammapy.makers import (
    SafeMaskMaker,
    MapDatasetMaker,
    ReflectedRegionsBackgroundMaker,
)

Select the data

We first set the datastore and retrieve a few observations from our source.

[3]:
datastore = DataStore.from_dir("$GAMMAPY_DATA/hess-dl3-dr1/")
obs_ids = [20326, 20327, 20349, 20350, 20396, 20397]
# In case you want to use all RX J1713 data in the HESS DR1
# other_ids=[20421, 20422, 20517, 20518, 20519, 20521, 20898, 20899, 20900]

observations = datastore.get_observations(obs_ids)

Prepare the datasets creation

Select the ON region

Here we take a simple 1 degree circular region because it fits well with the morphology of RX J1713-3945. More complex regions could be used e.g. ~regions.EllipseSkyRegion or ~regions.RectangleSkyRegion.

[4]:
target_position = SkyCoord(347.3, -0.5, unit="deg", frame="galactic")
radius = Angle("0.5 deg")
on_region = CircleSkyRegion(target_position, radius)

Define the geometries

This part is especially important. - We have to define first energy axes. They define the axes of the resulting gammapy.datasets.SpectrumDatasetOnOff. In particular, we have to be careful to the true energy axis: it has to cover a larger range than the reconstructed energy one. - Then we define the geometry itself. It does not need to be very finely binned and should enclose all the ON region. To limit CPU and memory usage, one should avoid using a much larger region.

[5]:
# The binning of the final spectrum is defined here.
energy_axis = MapAxis.from_energy_bounds(0.3, 40.0, 10, unit="TeV")

# Reduced IRFs are defined in true energy (i.e. not measured energy).
energy_axis_true = MapAxis.from_energy_bounds(
    0.05, 100, 30, unit="TeV", name="energy_true"
)

# Here we use 1.5 degree which is slightly larger than needed.
geom = WcsGeom.create(
    skydir=target_position,
    binsz=0.04,
    width=(1.5, 1.5),
    frame="galactic",
    proj="CAR",
    axes=[energy_axis],
)

Create the makers

First we instantiate the target gammapy.datasets.MapDataset.

[6]:
stacked = MapDataset.create(
    geom=geom, energy_axis_true=energy_axis_true, name="rxj-stacked"
)

Now we create its associated maker. Here we need to produce, counts, exposure and edisp (energy dispersion) entries. PSF and IRF background are not needed, therefore we don’t compute them.

[7]:
maker = MapDatasetMaker(selection=["counts", "exposure", "edisp"])

Now we create the OFF background maker for the spectra. If we have an exclusion region, we have to pass it here. We also define the safe range maker.

[8]:
bkg_maker = ReflectedRegionsBackgroundMaker()
safe_mask_maker = SafeMaskMaker(
    methods=["aeff-default", "aeff-max"], aeff_percent=10
)

Perform the data reduction loop.

We can now run over selected observations. For each of them, we: - create the map dataset and stack it on our target dataset. - squeeze the map dataset to a spectral dataset in the ON region - Compute the OFF and create a gammapy.datasets.SpectrumDatasetOnOff object - Run the safe mask maker on it - Add the gammapy.datasets.SpectrumDatasetOnOff to the list.

[9]:
%%time
datasets = Datasets()

for obs in observations:
    # A MapDataset is filled in this geometry
    dataset = maker.run(stacked, obs)
    # To make images, the resulting dataset cutout is stacked onto the final one
    stacked.stack(dataset)

    # Extract 1D spectrum
    spectrum_dataset = dataset.to_spectrum_dataset(
        on_region, name=f"obs-{obs.obs_id}"
    )
    # Compute OFF
    spectrum_dataset = bkg_maker.run(spectrum_dataset, obs)
    # Define safe mask
    spectrum_dataset = safe_mask_maker.run(spectrum_dataset, obs)
    # Append dataset to the list
    datasets.append(spectrum_dataset)
CPU times: user 3.25 s, sys: 208 ms, total: 3.46 s
Wall time: 3.49 s
[10]:
datasets.meta_table
[10]:
Table length=6
NAMETYPETELESCOPOBS_IDRA_PNTDEC_PNT
degdeg
str9str20str4int64float64float64
obs-20326SpectrumDatasetOnOffHESS20326259.29852294921875-39.76222229003906
obs-20327SpectrumDatasetOnOffHESS20327257.4773254394531-39.76222229003906
obs-20349SpectrumDatasetOnOffHESS20349259.29852294921875-39.76222229003906
obs-20350SpectrumDatasetOnOffHESS20350257.4773254394531-39.76222229003906
obs-20396SpectrumDatasetOnOffHESS20396258.3879089355469-39.06222152709961
obs-20397SpectrumDatasetOnOffHESS20397258.3879089355469-40.462223052978516

Explore the results

First let’s look at the data to see if our region is correct. We plot it over the excess. To do so we convert it to a pixel region using the WCS information stored on the geom.

[11]:
stacked.counts.sum_over_axes().smooth(width="0.05 deg").plot()
on_region.to_pixel(stacked.counts.geom.wcs).plot()
[11]:
<matplotlib.axes._subplots.WCSAxesSubplot at 0x122954518>
../_images/tutorials_extended_source_spectral_analysis_22_1.png

We now turn to the spectral datasets. We can peek at their content:

[12]:
datasets[0].peek()
[12]:
(<matplotlib.axes._subplots.AxesSubplot at 0x122db2be0>,
 <matplotlib.axes._subplots.AxesSubplot at 0x122f04d68>,
 <matplotlib.axes._subplots.AxesSubplot at 0x126f50f98>)
../_images/tutorials_extended_source_spectral_analysis_24_1.png

Cumulative excess and signficance

Finally, we can look at cumulative significance and number of excesses. This is done with the info_table method of gammapy.datasets.Datasets.

[13]:
info_table = datasets.info_table(cumulative=True)
[14]:
info_table
[14]:
Table length=6
namecountsbackgroundexcesssqrt_tsnprednpred_backgroundnpred_signalexposure_minexposure_maxlivetimeontimecounts_ratebackground_rateexcess_raten_binsn_fit_binsstat_typestat_sumcounts_offacceptanceacceptance_offalpha
m2 sm2 sss1 / s1 / s1 / s
str9float32float64float64float64float64float64float64float64float64float64float64float64float64float64int64int64str5float64float32float64float64float64
obs-20326450.0343.5106.54.405517626186096378.99999999999994378.99999999999994nan4309258.042057533422238157.512381551683.01683.00.267379679144385040.204099821746880570.063279857397504451010wstat41.619587630493804687.010.020.00.5
obs-20326845.0670.5174.55.215405763158092728.6666666666667728.6666666666667nan14662719.327748926832059034.4789453366.03366.00.25103980986333930.199197860962566840.051841948900772431010wstat68.997588727764551341.010.020.00.5
obs-203261288.0983.0305.07.4578748831235571084.66666666666671084.6666666666667nan25046841.7695577441241412951.51141645048.05048.00.25515055467511890.194730586370839930.060419968304278921010wstat110.361929295469341966.010.020.00.5
obs-203261725.01341.0384.08.0752575056283821468.99999999999981468.9999999999998nan29517317.734697961662189705.54993496730.06730.00.25631500742942050.199257057949479950.057057949479940561010wstat150.007960427822752682.010.020.00.5
obs-203262198.01653.0001220703125544.999877929687510.4001660373092091823.14663019796331823.1466301979633nan39224336.772897542071146470.36662868413.08413.00.26126233210507550.19648165007373260.064780682031342861010wstat193.548339870705833618.010.021.8874759674072270.44348058104515076
obs-203262626.02001.5001220703125624.499877929687510.8376845948849032198.5933688830692198.593368883069nan41782883.403462152500415797.96830310095.010095.00.260128776622090160.19826648064094230.061862295981147851010wstat213.441659457942454315.010.021.558830261230470.4485357403755188
[15]:
fig = plt.figure(figsize=(10, 6))
ax = fig.add_subplot(121)
ax.plot(
    info_table["livetime"].to("h"), info_table["excess"], marker="o", ls="none"
)
plt.xlabel("Livetime [h]")
plt.ylabel("Excess events")

ax = fig.add_subplot(122)
ax.plot(
    info_table["livetime"].to("h"),
    info_table["sqrt_ts"],
    marker="o",
    ls="none",
)
plt.xlabel("Livetime [h]")
plt.ylabel("Sqrt(TS)");
../_images/tutorials_extended_source_spectral_analysis_28_0.png

Perform spectral model fitting

Here we perform a joint fit.

We first create the model, here a simple powerlaw, and assign it to every dataset in the gammapy.datasets.Datasets.

[16]:
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, name="RXJ 1713")

datasets.models = [model]

Now we can run the fit

[17]:
fit_joint = Fit(datasets)
result_joint = fit_joint.run()
print(result_joint)
OptimizeResult

        backend    : minuit
        method     : minuit
        success    : True
        message    : Optimization terminated successfully.
        nfev       : 39
        total stat : 72.95

Explore the fit results

First the fitted parameters values and their errors.

[18]:
result_joint.parameters.to_table()
[18]:
Table length=3
namevalueunitminmaxfrozenerror
str9float64str14float64float64boolfloat64
index2.1021e+00nannanFalse6.657e-02
amplitude1.2869e-11cm-2 s-1 TeV-1nannanFalse1.023e-12
reference1.0000e+00TeVnannanTrue0.000e+00

Then plot the fit result to compare measured and expected counts. Rather than plotting them for each individual dataset, we stack all datasets and plot the fit result on the result.

[19]:
# First stack them all
reduced = datasets.stack_reduce()
# Assign the fitted model
reduced.models = model
# Plot the result
reduced.plot_fit();
../_images/tutorials_extended_source_spectral_analysis_36_0.png
[ ]: