-
Notifications
You must be signed in to change notification settings - Fork 187
/
spectrum_simulation.py
255 lines (198 loc) · 7.75 KB
/
spectrum_simulation.py
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
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
"""
1D spectrum simulation
======================
Simulate a number of spectral on-off observations of a source with a power-law spectral
model using the CTA 1DC response and fit them with the assumed spectral model.
Prerequisites
-------------
- Knowledge of spectral extraction and datasets used in gammapy, see
for instance the :doc:`spectral analysis
tutorial </tutorials/analysis-1d/spectral_analysis`
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:
- `~gammapy.datasets.SpectrumDatasetOnOff`
- `~gammapy.datasets.SpectrumDataset`
- `~gammapy.irf.load_cta_irfs`
- `~gammapy.modeling.models.PowerLawSpectralModel`
"""
import numpy as np
import astropy.units as u
from astropy.coordinates import Angle, SkyCoord
from regions import CircleSkyRegion
# %matplotlib inline
import matplotlib.pyplot as plt
######################################################################
# Setup
# -----
#
from IPython.display import display
from gammapy.data import Observation, observatory_locations
from gammapy.datasets import Datasets, SpectrumDataset, SpectrumDatasetOnOff
from gammapy.irf import load_cta_irfs
from gammapy.makers import SpectrumDatasetMaker
from gammapy.maps import MapAxis, RegionGeom
from gammapy.modeling import Fit
from gammapy.modeling.models import PowerLawSpectralModel, SkyModel
######################################################################
# Check setup
# -----------
from gammapy.utils.check import check_tutorials_setup
check_tutorials_setup()
######################################################################
# 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.
#
# 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)
# 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")
######################################################################
# 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"
)
location = observatory_locations["cta_south"]
obs = Observation.create(
pointing=pointing,
livetime=livetime,
irfs=irfs,
location=location,
)
print(obs)
######################################################################
# Simulate a spectra
#
# 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)
# Set the model on the dataset, and fake
dataset.models = model
dataset.fake(random_state=42)
print(dataset)
######################################################################
# 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
#
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)
######################################################################
# You can see that off counts are now simulated as well. We now simulate
# several spectra using the same set of observation conditions.
#
# %%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)
table = datasets.info_table()
display(table)
######################################################################
# Before moving on to the fit let’s have a look at the simulated
# observations.
#
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")
######################################################################
# Now, we fit each simulated spectrum individually
#
# %%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,
}
)
######################################################################
# We take a look at the distribution of the fitted indices. This matches
# very well with the spectrum that we initially injected.
#
fig, ax = plt.subplots()
index = np.array([_["index"] for _ in results])
ax.hist(index, bins=10, alpha=0.5)
ax.axvline(x=model_simu.parameters["index"].value, color="red")
ax.set_xlabel("Reconstructed spectral index")
print(f"index: {index.mean()} += {index.std()}")
plt.show()
######################################################################
# 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?
#