Note
Click here to download the full example code or to run this example in your browser via Binder
Simulating and fitting a time varying source#
Simulate and fit a time decaying light curve of a source using the CTA 1DC response.
Prerequisites#
To understand how a single binned simulation works, please refer to 1D spectrum simulation tutorial and 3D map simulation tutorial for 1D and 3D simulations respectively.
For details of light curve extraction using gammapy, refer to the two tutorials Light curves and Light curves for flares tutorial.
Context#
Frequently, studies of variable sources (eg: decaying GRB light curves, AGN flares, etc) require time variable simulations. For most use cases, generating an event list is an overkill, and it suffices to use binned simulations using a temporal model.
Objective: Simulate and fit a time decaying light curve of a source with CTA using the CTA 1DC response
Proposed approach#
We will simulate 10 spectral datasets within given time intervals (Good Time Intervals) following a given spectral (a power law) and temporal profile (an exponential decay, with a decay time of 6 hr ). These are then analysed using the light curve estimator to obtain flux points.
Modelling and fitting of lightcurves can be done either - directly on
the output of the LighCurveEstimator
(at the DL5 level) - fit the
simulated datasets (at the DL4 level)
In summary, necessary steps are:
Choose observation parameters including a list of
gammapy.data.GTI
Define temporal and spectral models from :ref:model-gallery as per science case
Perform the simulation (in 1D or 3D)
Extract the light curve from the reduced dataset as shown in Light curves tutorial.
Optionally, we show here how to fit the simulated datasets using a source model
Setup#
As usual, we’ll start with some general imports…
Setup#
from IPython.display import display
log = logging.getLogger(__name__)
And some gammapy specific imports
from gammapy.data import Observation, observatory_locations
from gammapy.datasets import Datasets, FluxPointsDataset, SpectrumDataset
from gammapy.estimators import LightCurveEstimator
from gammapy.irf import load_cta_irfs
from gammapy.makers import SpectrumDatasetMaker
from gammapy.maps import MapAxis, RegionGeom, TimeMapAxis
from gammapy.modeling import Fit
from gammapy.modeling.models import (
ExpDecayTemporalModel,
PowerLawSpectralModel,
SkyModel,
)
Check setup#
from gammapy.utils.check import check_tutorials_setup
check_tutorials_setup()
System:
python_executable : /Users/terrier/Code/anaconda3/envs/gammapy-dev/bin/python
python_version : 3.8.13
machine : x86_64
system : Darwin
Gammapy package:
version : 1.0rc2
path : /Users/terrier/Code/gammapy-dev/gammapy/gammapy
Other packages:
numpy : 1.22.4
scipy : 1.9.3
astropy : 5.1
regions : 0.6
click : 8.1.3
yaml : 6.0
IPython : 8.4.0
jupyterlab : 3.4.8
matplotlib : 3.5.3
pandas : 1.5.0
healpy : 1.16.1
iminuit : 2.17.0
sherpa : 4.15.0
naima : 0.10.0
emcee : 3.1.3
corner : 2.2.1
Gammapy environment variables:
GAMMAPY_DATA : /Users/terrier/Code/gammapy-dev/gammapy-data
We first define our preferred time format:
TimeMapAxis.time_format = "iso"
Simulating a light curve#
We will simulate 10 spectra between 300 GeV and 10 TeV using an
PowerLawSpectralModel
and a ExpDecayTemporalModel
. The important
thing to note here is how to attach a different GTI
to each dataset.
Since we use spectrum datasets here, we will use a RegionGeom
.
# Loading IRFs
irfs = load_cta_irfs(
"$GAMMAPY_DATA/cta-1dc/caldb/data/cta/1dc/bcf/South_z20_50h/irf_file.fits"
)
# 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"
)
geom = RegionGeom.create("galactic;circle(0, 0, 0.11)", axes=[energy_axis])
# Pointing position
pointing = SkyCoord(0.5, 0.5, unit="deg", frame="galactic")
/Users/terrier/Code/anaconda3/envs/gammapy-dev/lib/python3.8/site-packages/astropy/units/core.py:2042: UnitsWarning: '1/s/MeV/sr' did not parse as fits unit: Numeric factor not supported by FITS If this is meant to be a custom unit, define it with 'u.def_unit'. To have it recognized inside a file reader or other code, enable it with 'u.add_enabled_units'. For details, see https://docs.astropy.org/en/latest/units/combining_and_defining.html
warnings.warn(msg, UnitsWarning)
Note that observations are usually conducted in Wobble mode, in which the source is not in the center of the camera. This allows to have a symmetrical sky position from which background can be estimated.
# Define the source model: A combination of spectral and temporal model
gti_t0 = Time("2020-03-01")
spectral_model = PowerLawSpectralModel(
index=3, amplitude="1e-11 cm-2 s-1 TeV-1", reference="1 TeV"
)
temporal_model = ExpDecayTemporalModel(t0="6 h", t_ref=gti_t0.mjd * u.d)
model_simu = SkyModel(
spectral_model=spectral_model,
temporal_model=temporal_model,
name="model-simu",
)
# Look at the model
display(model_simu.parameters.to_table())
/Users/terrier/Code/anaconda3/envs/gammapy-dev/lib/python3.8/site-packages/astropy/units/quantity.py:611: RuntimeWarning: overflow encountered in exp
result = super().__array_ufunc__(function, method, *arrays, **kwargs)
/Users/terrier/Code/anaconda3/envs/gammapy-dev/lib/python3.8/site-packages/astropy/units/quantity.py:611: RuntimeWarning: invalid value encountered in subtract
result = super().__array_ufunc__(function, method, *arrays, **kwargs)
type name value unit ... max frozen is_norm link
-------- --------- ---------- -------------- ... --- ------ ------- ----
spectral index 3.0000e+00 ... nan False False
spectral amplitude 1.0000e-11 cm-2 s-1 TeV-1 ... nan False True
spectral reference 1.0000e+00 TeV ... nan True False
temporal t0 6.0000e+00 h ... nan False False
temporal t_ref 5.8909e+04 d ... nan True False
Now, define the start and observation livetime wrt to the reference
time, gti_t0
Now perform the simulations
datasets = Datasets()
empty = SpectrumDataset.create(
geom=geom, energy_axis_true=energy_axis_true, name="empty"
)
maker = SpectrumDatasetMaker(selection=["exposure", "background", "edisp"])
for idx in range(n_obs):
obs = Observation.create(
pointing=pointing,
livetime=lvtm[idx],
tstart=tstart[idx],
irfs=irfs,
reference_time=gti_t0,
obs_id=idx,
location=observatory_locations["cta_south"],
)
empty_i = empty.copy(name=f"dataset-{idx}")
dataset = maker.run(empty_i, obs)
dataset.models = model_simu
dataset.fake()
datasets.append(dataset)
The reduced datasets have been successfully simulated. Let’s take a quick look into our datasets.
display(datasets.info_table())
name counts excess ... n_fit_bins stat_type stat_sum
...
--------- ------ ------------------ ... ---------- --------- -------------------
dataset-0 789 768.6812133789062 ... 9 cash nan
dataset-1 299 289.7641959176923 ... 9 cash -1825.4990055104624
dataset-2 278 268.3947637544 ... 9 cash -1663.5292807438886
dataset-3 302 287.22271346830775 ... 9 cash -1799.3614644502222
dataset-4 196 181.22271346830775 ... 9 cash -1025.081000048927
dataset-5 203 184.52839183538467 ... 9 cash -1097.4875261091647
dataset-6 40 25.22271346830774 ... 9 cash -79.10846210737706
dataset-7 43 23.789527508800063 ... 9 cash -110.38846690151904
dataset-8 33 17.114416978430825 ... 9 cash -54.39811202595201
dataset-9 41 23.6366883252616 ... 9 cash -79.5182330516443
Extract the lightcurve#
This section uses standard light curve estimation tools for a 1D extraction. Only a spectral model needs to be defined in this case. Since the estimator returns the integrated flux separately for each time bin, the temporal model need not be accounted for at this stage. We extract the lightcurve in 3 energy binsç
# Define the model:
spectral_model = PowerLawSpectralModel(
index=3, amplitude="1e-11 cm-2 s-1 TeV-1", reference="1 TeV"
)
model_fit = SkyModel(spectral_model=spectral_model, name="model-fit")
# Attach model to all datasets
datasets.models = model_fit
lc_maker_1d = LightCurveEstimator(
energy_edges=[0.3, 0.6, 1.0, 10] * u.TeV,
source="model-fit",
selection_optional=["ul"],
)
lc_1d = lc_maker_1d.run(datasets)
fig, ax = plt.subplots(
figsize=(8, 6),
gridspec_kw={"left": 0.16, "bottom": 0.2, "top": 0.98, "right": 0.98},
)
lc_1d.plot(ax=ax, marker="o", axis_name="time", sed_type="flux")
<AxesSubplot:xlabel='Time [iso]', ylabel='flux [1 / (cm2 s)]'>
Fitting temporal models#
We have the reconstructed lightcurve at this point. Now we want to fit a profile to the obtained light curves, using a joint fitting across the different datasets, while simultaneously minimising across the temporal model parameters as well. The temporal models can be applied
directly on the obtained lightcurve
on the simulated datasets
Fitting the obtained light curve#
The fitting will proceed through a joint fit of the flux points. First,
we need to obtain a set of FluxPointDatasets
, one for each time bin
# Create the datasets by iterating over the returned lightcurve
datasets = Datasets()
for idx, fp in enumerate(lc_1d.iter_by_axis(axis_name="time")):
dataset = FluxPointsDataset(data=fp, name=f"time-bin-{idx}")
datasets.append(dataset)
We will fit the amplitude, spectral index and the decay time scale. Note
that t_ref
should be fixed by default for the
ExpDecayTemporalModel
.
# Define the model:
spectral_model1 = PowerLawSpectralModel(
index=2.0, amplitude="1e-12 cm-2 s-1 TeV-1", reference="1 TeV"
)
temporal_model1 = ExpDecayTemporalModel(t0="10 h", t_ref=gti_t0.mjd * u.d)
model = SkyModel(
spectral_model=spectral_model1,
temporal_model=temporal_model1,
name="model-test",
)
datasets.models = model
/Users/terrier/Code/anaconda3/envs/gammapy-dev/lib/python3.8/site-packages/astropy/units/quantity.py:611: RuntimeWarning: overflow encountered in exp
result = super().__array_ufunc__(function, method, *arrays, **kwargs)
/Users/terrier/Code/anaconda3/envs/gammapy-dev/lib/python3.8/site-packages/astropy/units/quantity.py:611: RuntimeWarning: invalid value encountered in subtract
result = super().__array_ufunc__(function, method, *arrays, **kwargs)
Do a joint fit
/Users/terrier/Code/anaconda3/envs/gammapy-dev/lib/python3.8/site-packages/astropy/units/quantity.py:611: RuntimeWarning: overflow encountered in exp
result = super().__array_ufunc__(function, method, *arrays, **kwargs)
/Users/terrier/Code/anaconda3/envs/gammapy-dev/lib/python3.8/site-packages/astropy/units/quantity.py:611: RuntimeWarning: invalid value encountered in subtract
result = super().__array_ufunc__(function, method, *arrays, **kwargs)
Now let’s plot model and data. We plot only the normalisation of the temporal model in relative units for one particular energy range
fig, ax = plt.subplots(
figsize=(8, 6),
gridspec_kw={"left": 0.16, "bottom": 0.2, "top": 0.98, "right": 0.98},
)
lc_1TeV_10TeV = lc_1d.slice_by_idx({"energy": slice(2, 3)})
lc_1TeV_10TeV.plot(ax=ax, sed_type="norm", axis_name="time")
time_range = lc_1TeV_10TeV.geom.axes["time"].time_bounds
temporal_model1.plot(ax=ax, time_range=time_range, label="Best fit model")
ax.set_yscale("linear")
ax.legend()
<matplotlib.legend.Legend object at 0x142f39340>
Fit the datasets#
Here, we apply the models directly to the simulated datasets.
For modelling and fitting more complex flares, you should attach the
relevant model to each group of datasets
. The parameters of a model
in a given group of dataset will be tied. For more details on joint
fitting in Gammapy, see 3D detailed analysis
# Define the model:
spectral_model2 = PowerLawSpectralModel(
index=2.0, amplitude="1e-12 cm-2 s-1 TeV-1", reference="1 TeV"
)
temporal_model2 = ExpDecayTemporalModel(t0="10 h", t_ref=gti_t0.mjd * u.d)
model2 = SkyModel(
spectral_model=spectral_model2,
temporal_model=temporal_model2,
name="model-test2",
)
display(model2.parameters.to_table())
datasets.models = model2
/Users/terrier/Code/anaconda3/envs/gammapy-dev/lib/python3.8/site-packages/astropy/units/quantity.py:611: RuntimeWarning: overflow encountered in exp
result = super().__array_ufunc__(function, method, *arrays, **kwargs)
/Users/terrier/Code/anaconda3/envs/gammapy-dev/lib/python3.8/site-packages/astropy/units/quantity.py:611: RuntimeWarning: invalid value encountered in subtract
result = super().__array_ufunc__(function, method, *arrays, **kwargs)
type name value unit ... max frozen is_norm link
-------- --------- ---------- -------------- ... --- ------ ------- ----
spectral index 2.0000e+00 ... nan False False
spectral amplitude 1.0000e-12 cm-2 s-1 TeV-1 ... nan False True
spectral reference 1.0000e+00 TeV ... nan True False
temporal t0 1.0000e+01 h ... nan False False
temporal t_ref 5.8909e+04 d ... nan True False
Do a joint fit
/Users/terrier/Code/anaconda3/envs/gammapy-dev/lib/python3.8/site-packages/astropy/units/quantity.py:611: RuntimeWarning: overflow encountered in exp
result = super().__array_ufunc__(function, method, *arrays, **kwargs)
/Users/terrier/Code/anaconda3/envs/gammapy-dev/lib/python3.8/site-packages/astropy/units/quantity.py:611: RuntimeWarning: invalid value encountered in subtract
result = super().__array_ufunc__(function, method, *arrays, **kwargs)
type name value unit ... max frozen is_norm link
-------- --------- ---------- -------------- ... --- ------ ------- ----
spectral index 3.0296e+00 ... nan False False
spectral amplitude 9.0456e-12 cm-2 s-1 TeV-1 ... nan False True
spectral reference 1.0000e+00 TeV ... nan True False
temporal t0 6.4523e+00 h ... nan False False
temporal t_ref 5.8909e+04 d ... nan True False
We see that the fitted parameters are consistent between fitting flux points and datasets, and match well with the simulated ones
Exercises#
Re-do the analysis with
MapDataset
instead ofSpectralDataset
Model the flare of PKS 2155-304 which you obtained using the
Light curves for flares tutorial. Use a combination of a Gaussian and Exponential flare profiles.
Do a joint fitting of the datasets.
Total running time of the script: ( 0 minutes 34.601 seconds)