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

1D spectrum simulation#

Prerequisites#

Context#

To simulate a specific observation, it is not always necessary to simulate the full photon list. For many uses cases, simulating directly a reduced binned dataset is enough: the IRFs reduced in the correct geometry are combined with a source model to predict an actual number of counts per bin. The latter is then used to simulate a reduced dataset using Poisson probability distribution.

This can be done to check the feasibility of a measurement, to test whether fitted parameters really provide a good fit to the data etc.

Here we will see how to perform a 1D spectral simulation of a CTA observation, in particular, we will generate OFF observations following the template background stored in the CTA IRFs.

Objective: simulate a number of spectral ON-OFF observations of a source with a power-law spectral model with CTA using the CTA 1DC response, fit them with the assumed spectral model and check that the distribution of fitted parameters is consistent with the input values.

Proposed approach#

We will use the following classes:

Setup#

[1]:
%matplotlib inline
import matplotlib.pyplot as plt
[2]:
import numpy as np
import astropy.units as u
from astropy.coordinates import SkyCoord, Angle
from regions import CircleSkyRegion
from gammapy.datasets import SpectrumDatasetOnOff, SpectrumDataset, Datasets
from gammapy.makers import SpectrumDatasetMaker
from gammapy.modeling import Fit
from gammapy.modeling.models import (
    PowerLawSpectralModel,
    SkyModel,
)
from gammapy.irf import load_cta_irfs
from gammapy.data import Observation, observatory_locations
from gammapy.maps import MapAxis, RegionGeom

Simulation of a single spectrum#

To do a simulation, we need to define the observational parameters like the livetime, the offset, the assumed integration radius, the energy range to perform the simulation for and the choice of spectral model. We then use an in-memory observation which is convolved with the IRFs to get the predicted number of counts. This is Poission fluctuated using the fake() to get the simulated counts for each observation.

[3]:
# Define simulation parameters parameters
livetime = 1 * u.h

pointing = SkyCoord(0, 0, unit="deg", frame="galactic")
offset = 0.5 * u.deg

# Reconstructed and true energy axis
energy_axis = MapAxis.from_edges(
    np.logspace(-0.5, 1.0, 10), unit="TeV", name="energy", interp="log"
)
energy_axis_true = MapAxis.from_edges(
    np.logspace(-1.2, 2.0, 31), unit="TeV", name="energy_true", interp="log"
)

on_region_radius = Angle("0.11 deg")

center = pointing.directional_offset_by(
    position_angle=0 * u.deg, separation=offset
)
on_region = CircleSkyRegion(center=center, radius=on_region_radius)
[4]:
# Define spectral model - a simple Power Law in this case
model_simu = PowerLawSpectralModel(
    index=3.0,
    amplitude=2.5e-12 * u.Unit("cm-2 s-1 TeV-1"),
    reference=1 * u.TeV,
)
print(model_simu)
# we set the sky model used in the dataset
model = SkyModel(spectral_model=model_simu, name="source")
PowerLawSpectralModel

  type      name     value         unit        error   min max frozen is_norm link
-------- --------- ---------- -------------- --------- --- --- ------ ------- ----
spectral     index 3.0000e+00                0.000e+00 nan nan  False   False
spectral amplitude 2.5000e-12 cm-2 s-1 TeV-1 0.000e+00 nan nan  False    True
spectral reference 1.0000e+00            TeV 0.000e+00 nan nan   True   False
[5]:
# Load the IRFs
# In this simulation, we use the CTA-1DC irfs shipped with gammapy.
irfs = load_cta_irfs(
    "$GAMMAPY_DATA/cta-1dc/caldb/data/cta/1dc/bcf/South_z20_50h/irf_file.fits"
)
Invalid unit found in background table! Assuming (s-1 MeV-1 sr-1)
[6]:
location = observatory_locations["cta_south"]
obs = Observation.create(
    pointing=pointing,
    livetime=livetime,
    irfs=irfs,
    location=location,
)
print(obs)
Observation

        obs id            : 0
        tstart            : 51544.00
        tstop             : 51544.04
        duration          : 3600.00 s
        pointing (icrs)   : 266.4 deg, -28.9 deg

        deadtime fraction : 0.0%

[7]:
# Make the SpectrumDataset
geom = RegionGeom.create(region=on_region, axes=[energy_axis])

dataset_empty = SpectrumDataset.create(
    geom=geom, energy_axis_true=energy_axis_true, name="obs-0"
)
maker = SpectrumDatasetMaker(selection=["exposure", "edisp", "background"])

dataset = maker.run(dataset_empty, obs)
[8]:
# Set the model on the dataset, and fake
dataset.models = model
dataset.fake(random_state=42)
print(dataset)
SpectrumDataset
---------------

  Name                            : obs-0

  Total counts                    : 298
  Total background counts         : 22.29
  Total excess counts             : 275.71

  Predicted counts                : 303.66
  Predicted background counts     : 22.29
  Predicted excess counts         : 281.37

  Exposure min                    : 2.53e+08 m2 s
  Exposure max                    : 1.77e+10 m2 s

  Number of total bins            : 9
  Number of fit bins              : 9

  Fit statistic type              : cash
  Fit statistic value (-2 log(L)) : -1811.58

  Number of models                : 1
  Number of parameters            : 3
  Number of free parameters       : 2

  Component 0: SkyModel

    Name                      : source
    Datasets names            : None
    Spectral model type       : PowerLawSpectralModel
    Spatial  model type       :
    Temporal model type       :
    Parameters:
      index                         :      3.000   +/-    0.00
      amplitude                     :   2.50e-12   +/- 0.0e+00 1 / (cm2 s TeV)
      reference             (frozen):      1.000       TeV


You can see that background counts are now simulated

On-Off analysis#

To do an on off spectral analysis, which is the usual science case, the standard would be to use SpectrumDatasetOnOff, which uses the acceptance to fake off-counts

[9]:
dataset_on_off = SpectrumDatasetOnOff.from_spectrum_dataset(
    dataset=dataset, acceptance=1, acceptance_off=5
)
dataset_on_off.fake(npred_background=dataset.npred_background())
print(dataset_on_off)
SpectrumDatasetOnOff
--------------------

  Name                            : cUYhcb1h

  Total counts                    : 301
  Total background counts         : 17.00
  Total excess counts             : 284.00

  Predicted counts                : 298.40
  Predicted background counts     : 17.03
  Predicted excess counts         : 281.37

  Exposure min                    : 2.53e+08 m2 s
  Exposure max                    : 1.77e+10 m2 s

  Number of total bins            : 9
  Number of fit bins              : 9

  Fit statistic type              : wstat
  Fit statistic value (-2 log(L)) : 7.27

  Number of models                : 1
  Number of parameters            : 3
  Number of free parameters       : 2

  Component 0: SkyModel

    Name                      : source
    Datasets names            : None
    Spectral model type       : PowerLawSpectralModel
    Spatial  model type       :
    Temporal model type       :
    Parameters:
      index                         :      3.000   +/-    0.00
      amplitude                     :   2.50e-12   +/- 0.0e+00 1 / (cm2 s TeV)
      reference             (frozen):      1.000       TeV

    Total counts_off                : 85
  Acceptance                      : 9
  Acceptance off                  : 45

You can see that off counts are now simulated as well. We now simulate several spectra using the same set of observation conditions.

[10]:
%%time

n_obs = 100
datasets = Datasets()

for idx in range(n_obs):
    dataset_on_off.fake(
        random_state=idx, npred_background=dataset.npred_background()
    )
    dataset_fake = dataset_on_off.copy(name=f"obs-{idx}")
    dataset_fake.meta_table["OBS_ID"] = [idx]
    datasets.append(dataset_fake)
CPU times: user 765 ms, sys: 20.4 ms, total: 785 ms
Wall time: 929 ms
[11]:
table = datasets.info_table()
table
[11]:
Table length=100
namecountsexcesssqrt_tsbackgroundnprednpred_backgroundnpred_signalexposure_minexposure_maxlivetimeontimecounts_ratebackground_rateexcess_raten_binsn_fit_binsstat_typestat_sumcounts_offacceptanceacceptance_offalpha
m2 sm2 sss1 / s1 / s1 / s
str6int64float64float64float64float64float64float64float64float64float64float64float64float64float64int64int64str5float64int64float64float64float64
obs-0317298.600006103515627.0824018081312618.39999961853027368.1666675131354168.16666751313541nan252718176.017719697408.03600.0000000000023600.0000000000020.088055555555555510.00511111100514729560.0829444461398654299wstat738.7245834271496929.045.00.20000000298023224
obs-1275253.023.7678536548728522.064.1666666666666964.16666666666669nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.076388888888888850.0061111111111111080.0702777777777777499wstat575.7795127847381109.045.00.2
obs-2293272.425.1711055540465520.666.066.0nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.081388888888888840.005722222222222220.0756666666666666299wstat645.49938240753031039.045.00.2
obs-3280257.623.98295173740535422.465.3333333333333465.33333333333334nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.077777777777777740.0062222222222222180.0715555555555555399wstat585.92415469858721129.045.000000000000010.19999999999999998
obs-4337316.427.68270994518474720.673.3333333333333473.33333333333334nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.093611111111111060.005722222222222220.0878888888888888499wstat787.3147239494481039.045.00.2
obs-5283258.623.72715478234789524.40000000000000267.5000000000000167.50000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.078611111111111080.0067777777777777750.071833333333333399wstat574.85257377274991229.045.00.2
obs-6330307.626.88918447572786622.40000000000000673.6666666666666773.66666666666667nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.091666666666666620.0062222222222222210.085444444444444499wstat734.24147552859371129.044.999999999999990.20000000000000004
obs-7283257.223.4308797479585325.868.6666666666666768.66666666666667nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.078611111111111080.0071666666666666630.0714444444444444199wstat552.36987794480441299.045.00.2
obs-8308284.625.4204927332833123.40000000000000270.8333333333333470.83333333333334nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.085555555555555510.0064999999999999970.0790555555555555299wstat652.26215725676331179.045.00.2
obs-9299278.625.5708507148686320.466.8333333333333466.83333333333334nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.083055555555555510.0056666666666666640.0773888888888888599wstat666.90626702607861029.045.000000000000010.19999999999999998
obs-10310294.827.48897216135677415.264.3333333333333464.33333333333334nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.086111111111111070.004222222222222220.0818888888888888499wstat768.5404492529331769.045.000000000000010.19999999999999998
obs-11285261.023.93383374545468524.00000000000000467.5000000000000167.50000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.079166666666666620.00666666666666666450.0724999999999999799wstat583.4487532927331209.044.999999999999990.20000000000000004
obs-12299278.025.432526389511721.00000000000000467.3333333333333667.33333333333336nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.083055555555555510.0058333333333333310.0772222222222221899wstat655.27682464803571059.044.999999999999990.20000000000000004
obs-13309287.425.87740623501261821.669.5000000000000169.50000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.085833333333333290.00599999999999999750.0798333333333332899wstat672.3557482743361089.045.00.2
obs-14320297.426.2828260281925522.59999999999999872.1666666666666772.16666666666667nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.088888888888888850.0062777777777777740.0826111111111110799wstat697.59497454151691139.045.000000000000010.19999999999999998
obs-15283261.024.2540897930413722.065.5000000000000165.50000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.078611111111111080.0061111111111111080.0724999999999999799wstat592.87545637398891109.045.00.2
obs-16298273.824.66421820169735224.20000000000000369.8333333333333669.83333333333336nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.082777777777777740.006722222222222220.0760555555555555299wstat619.25841883868051219.045.00.2
obs-17301272.824.0118015192522728.20000000000000373.6666666666666773.66666666666667nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.083611111111111070.007833333333333330.0757777777777777599wstat588.49763795040561419.045.00.2
obs-18311289.225.94747697694981521.870.070.0nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.086388888888888850.0060555555555555530.0803333333333332899wstat682.62503087509841099.045.00.2
obs-19280257.623.98295173740537622.40000000000000265.3333333333333465.33333333333334nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.077777777777777740.0062222222222222190.0715555555555555399wstat590.99956498143291129.045.00.2
obs-20327304.226.6332428682802422.80000000000000473.5000000000000373.50000000000003nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.09083333333333330.00633333333333333140.0844999999999999599wstat720.48401558232791149.044.999999999999990.20000000000000004
obs-21286266.825.0833207306816419.263.6666666666666863.66666666666668nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.07944444444444440.00533333333333333060.0741111111111110799wstat637.0646310076313969.045.000000000000010.19999999999999998
.....................................................................
obs-77304285.026.19215711358046419.066.5000000000000166.50000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.08444444444444440.0052777777777777750.0791666666666666299wstat696.82356000151959.045.00.2
obs-78337314.427.2333118017933222.675.0000000000000175.00000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.093611111111111060.0062777777777777750.0873333333333332899wstat750.40262937493821139.045.00.2
obs-79315294.626.495723628343920.469.569.5nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.087499999999999950.0056666666666666640.081833333333333399wstat711.05460461397391029.045.000000000000010.19999999999999998
obs-80304282.825.6786874038582321.20000000000000368.3333333333333668.33333333333336nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.08444444444444440.0058888888888888870.0785555555555555299wstat669.96543897555061069.044.999999999999990.20000000000000004
obs-81301281.625.92234741483462219.40000000000000266.3333333333333466.33333333333334nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.083611111111111070.0053888888888888870.0782222222222221999wstat678.9813616402096979.045.00.2
obs-82280258.024.07258704408423722.00000000000000465.065.0nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.077777777777777740.0061111111111111090.0716666666666666399wstat581.86100831608771109.044.999999999999990.20000000000000004
obs-83303280.825.39492828229910522.20000000000000369.0000000000000169.00000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.084166666666666630.0061666666666666640.0779999999999999699wstat654.46508964703431119.044.999999999999990.20000000000000004
obs-84269247.823.58024714002698721.20000000000000362.50000000000001462.500000000000014nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.074722222222222190.0058888888888888870.068833333333333399wstat577.23628917662691069.044.999999999999990.20000000000000004
obs-85297274.624.99893584573557322.468.1666666666666768.16666666666667nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.082499999999999960.0062222222222222180.0762777777777777599wstat631.62983611980531129.045.000000000000010.19999999999999998
obs-86286268.225.4234089695154617.80000000000000462.5000000000000162.50000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.07944444444444440.0049444444444444430.0744999999999999699wstat655.2082535886425899.044.999999999999990.20000000000000004
obs-87333313.227.64685145114721719.872.0000000000000172.00000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.092499999999999960.0054999999999999970.0869999999999999599wstat770.1313412321031999.045.00.2
obs-88315296.827.017205925109718.20000000000000367.6666666666666967.66666666666669nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.087499999999999950.0050555555555555540.082444444444444499wstat734.6476098262428919.044.999999999999990.20000000000000004
obs-89287265.024.4945655868051522.00000000000000466.1666666666666766.16666666666667nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.079722222222222180.0061111111111111090.0736111111111110799wstat611.23340636850891109.044.999999999999990.20000000000000004
obs-90286267.025.13122188704343819.00000000000000463.5000000000000163.50000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.07944444444444440.0052777777777777760.0741666666666666399wstat635.8304116110406959.044.999999999999990.20000000000000004
obs-91285259.823.6775459193106925.20000000000000368.5000000000000168.50000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.079166666666666620.00699999999999999750.0721666666666666399wstat573.35647429633811269.045.00.2
obs-92313289.425.66493542019417623.671.8333333333333671.83333333333336nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.08694444444444440.0065555555555555520.0803888888888888499wstat669.44961428251961189.045.00.2
obs-93302283.226.12386752260549718.866.0000000000000166.00000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.083888888888888850.005222222222222220.0786666666666666299wstat700.5066738339824949.045.00.2
obs-94322300.226.57402975747320221.871.8333333333333471.83333333333334nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.08944444444444440.0060555555555555530.0833888888888888599wstat715.25704980371471099.045.00.2
obs-95305280.625.03073126049315224.40000000000000671.1666666666666771.16666666666667nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.084722222222222180.0067777777777777760.0779444444444444199wstat643.28468794047931229.044.999999999999990.20000000000000004
obs-96301277.424.96984596942142823.669.8333333333333469.83333333333334nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.083611111111111070.0065555555555555520.0770555555555555299wstat626.73968424286241189.045.00.2
obs-97290271.225.41798219445482618.864.0000000000000164.00000000000001nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.080555555555555520.005222222222222220.075333333333333399wstat650.9762484493856949.045.00.2
obs-98301280.625.68783296467500720.40000000000000267.1666666666666767.16666666666667nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.083611111111111070.00566666666666666450.0779444444444444199wstat664.56200994049521029.045.00.2
obs-99323302.226.8570784237687120.871.1666666666666971.16666666666669nan252718170.9728769417719697919.599293600.0000000000023600.0000000000020.089722222222222180.0057777777777777750.083944444444444499wstat733.15765987680171049.045.00.2

Before moving on to the fit let’s have a look at the simulated observations.

[12]:
fix, axes = plt.subplots(1, 3, figsize=(12, 4))
axes[0].hist(table["counts"])
axes[0].set_xlabel("Counts")
axes[1].hist(table["counts_off"])
axes[1].set_xlabel("Counts Off")
axes[2].hist(table["excess"])
axes[2].set_xlabel("excess");
../../../_images/tutorials_analysis_1D_spectrum_simulation_19_0.png

Now, we fit each simulated spectrum individually

[13]:
%%time
results = []

fit = Fit()

for dataset in datasets:
    dataset.models = model.copy()
    result = fit.optimize(dataset)
    results.append(
        {
            "index": result.parameters["index"].value,
            "amplitude": result.parameters["amplitude"].value,
        }
    )
CPU times: user 15.5 s, sys: 247 ms, total: 15.7 s
Wall time: 19.9 s

We take a look at the distribution of the fitted indices. This matches very well with the spectrum that we initially injected.

[14]:
index = np.array([_["index"] for _ in results])
plt.hist(index, bins=10, alpha=0.5)
plt.axvline(x=model_simu.parameters["index"].value, color="red")
print(f"index: {index.mean()} += {index.std()}")
index: 3.00369255504466 += 0.08081469527083394
../../../_images/tutorials_analysis_1D_spectrum_simulation_23_1.png

Exercises#

  • Change the observation time to something longer or shorter. Does the observation and spectrum results change as you expected?

  • Change the spectral model, e.g. add a cutoff at 5 TeV, or put a steep-spectrum source with spectral index of 4.0

  • Simulate spectra with the spectral model we just defined. How much observation duration do you need to get back the injected parameters?

[ ]: