{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "
\n", "\n", "**This is a fixed-text formatted version of a Jupyter notebook**\n", "\n", "- Try online [![Binder](https://static.mybinder.org/badge.svg)](https://mybinder.org/v2/gh/gammapy/gammapy-webpage/v0.18.2?urlpath=lab/tree/spectrum_simulation.ipynb)\n", "- You can contribute with your own notebooks in this\n", "[GitHub repository](https://github.com/gammapy/gammapy/tree/master/docs/tutorials).\n", "- **Source files:**\n", "[spectrum_simulation.ipynb](../_static/notebooks/spectrum_simulation.ipynb) |\n", "[spectrum_simulation.py](../_static/notebooks/spectrum_simulation.py)\n", "
\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# Spectrum simulation\n", "\n", "## Prerequisites\n", "\n", "- Knowledge of spectral extraction and datasets used in gammapy, see for instance the [spectral analysis tutorial](spectrum_analysis.ipynb)\n", "\n", "## Context\n", "\n", "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.\n", "\n", "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.\n", "\n", "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.\n", "\n", "**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.**\n", "\n", "## Proposed approach:\n", "\n", "We will use the following classes:\n", "\n", "* `~gammapy.datasets.SpectrumDatasetOnOff`\n", "* `~gammapy.datasets.SpectrumDataset`\n", "* `~gammapy.irf.load_cta_irfs`\n", "* `~gammapy.modeling.models.PowerLawSpectralModel`" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Setup\n" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "%matplotlib inline\n", "import matplotlib.pyplot as plt" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "import numpy as np\n", "import astropy.units as u\n", "from astropy.coordinates import SkyCoord, Angle\n", "from regions import CircleSkyRegion\n", "from gammapy.datasets import SpectrumDatasetOnOff, SpectrumDataset, Datasets\n", "from gammapy.makers import SpectrumDatasetMaker\n", "from gammapy.modeling import Fit\n", "from gammapy.modeling.models import (\n", " PowerLawSpectralModel,\n", " SkyModel,\n", ")\n", "from gammapy.irf import load_cta_irfs\n", "from gammapy.data import Observation\n", "from gammapy.maps import MapAxis" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Simulation of a single spectrum\n", "\n", "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. " ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "# Define simulation parameters parameters\n", "livetime = 1 * u.h\n", "\n", "pointing = SkyCoord(0, 0, unit=\"deg\", frame=\"galactic\")\n", "offset = 0.5 * u.deg\n", "\n", "# Reconstructed and true energy axis\n", "energy_axis = MapAxis.from_edges(\n", " np.logspace(-0.5, 1.0, 10), unit=\"TeV\", name=\"energy\", interp=\"log\"\n", ")\n", "energy_axis_true = MapAxis.from_edges(\n", " np.logspace(-1.2, 2.0, 31), unit=\"TeV\", name=\"energy_true\", interp=\"log\"\n", ")\n", "\n", "on_region_radius = Angle(\"0.11 deg\")\n", "\n", "center = pointing.directional_offset_by(\n", " position_angle=0 * u.deg, separation=offset\n", ")\n", "on_region = CircleSkyRegion(center=center, radius=on_region_radius)" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "PowerLawSpectralModel\n", "\n", " name value unit min max frozen error \n", "--------- ---------- -------------- --- --- ------ ---------\n", " index 3.0000e+00 nan nan False 0.000e+00\n", "amplitude 2.5000e-12 cm-2 s-1 TeV-1 nan nan False 0.000e+00\n", "reference 1.0000e+00 TeV nan nan True 0.000e+00\n" ] } ], "source": [ "# Define spectral model - a simple Power Law in this case\n", "model_simu = PowerLawSpectralModel(\n", " index=3.0,\n", " amplitude=2.5e-12 * u.Unit(\"cm-2 s-1 TeV-1\"),\n", " reference=1 * u.TeV,\n", ")\n", "print(model_simu)\n", "# we set the sky model used in the dataset\n", "model = SkyModel(spectral_model=model_simu, name=\"source\")" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "Invalid unit found in background table! Assuming (s-1 MeV-1 sr-1)\n" ] } ], "source": [ "# Load the IRFs\n", "# In this simulation, we use the CTA-1DC irfs shipped with gammapy.\n", "irfs = load_cta_irfs(\n", " \"$GAMMAPY_DATA/cta-1dc/caldb/data/cta/1dc/bcf/South_z20_50h/irf_file.fits\"\n", ")" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Observation\n", "\n", "\tobs id : 0 \n", " \ttstart : 51544.00\n", "\ttstop : 51544.04\n", "\tduration : 3600.00 s\n", "\tpointing (icrs) : 266.4 deg, -28.9 deg\n", "\n", "\tdeadtime fraction : 0.0%\n", "\n" ] } ], "source": [ "obs = Observation.create(pointing=pointing, livetime=livetime, irfs=irfs)\n", "print(obs)" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [], "source": [ "# Make the SpectrumDataset\n", "dataset_empty = SpectrumDataset.create(\n", " e_reco=energy_axis, e_true=energy_axis_true, region=on_region, name=\"obs-0\"\n", ")\n", "maker = SpectrumDatasetMaker(selection=[\"exposure\", \"edisp\", \"background\"])\n", "\n", "dataset = maker.run(dataset_empty, obs)" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "SpectrumDataset\n", "---------------\n", "\n", " Name : obs-0 \n", "\n", " Total counts : 298 \n", " Total background counts : 22.32\n", " Total excess counts : 275.68\n", "\n", " Predicted counts : 303.69\n", " Predicted background counts : 22.32\n", " Predicted excess counts : 281.37\n", "\n", " Exposure min : 2.53e+08 m2 s\n", " Exposure max : 1.77e+10 m2 s\n", "\n", " Number of total bins : 9 \n", " Number of fit bins : 9 \n", "\n", " Fit statistic type : cash\n", " Fit statistic value (-2 log(L)) : -1811.58\n", "\n", " Number of models : 1 \n", " Number of parameters : 3\n", " Number of free parameters : 2\n", "\n", " Component 0: SkyModel\n", " \n", " Name : source\n", " Datasets names : None\n", " Spectral model type : PowerLawSpectralModel\n", " Spatial model type : \n", " Temporal model type : \n", " Parameters:\n", " index : 3.000 \n", " amplitude : 2.50e-12 1 / (cm2 s TeV)\n", " reference (frozen) : 1.000 TeV \n", " \n", " \n" ] } ], "source": [ "# Set the model on the dataset, and fake\n", "dataset.models = model\n", "dataset.fake(random_state=42)\n", "print(dataset)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "You can see that backgound counts are now simulated" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### On-Off analysis\n", "\n", "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 " ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "SpectrumDatasetOnOff\n", "--------------------\n", "\n", " Name : obs-0 \n", "\n", " Total counts : 285 \n", " Total off counts : 112.00\n", " Total background counts : 22.40\n", " Total excess counts : 262.60\n", "\n", " Predicted counts : 303.54\n", " Predicted background counts : 22.17\n", " Predicted excess counts : 281.37\n", "\n", " Exposure min : 2.53e+08 m2 s\n", " Exposure max : 1.77e+10 m2 s\n", "\n", " Acceptance mean : 1.000\n", " Acceptance off : 45.000\n", "\n", " Number of total bins : 9 \n", " Number of fit bins : 9 \n", "\n", " Fit statistic type : wstat\n", " Fit statistic value (-2 log(L)) : 4.48\n", "\n", " Number of models : 1 \n", " Number of parameters : 3\n", " Number of free parameters : 2\n", "\n", " Component 0: SkyModel\n", " \n", " Name : source\n", " Datasets names : None\n", " Spectral model type : PowerLawSpectralModel\n", " Spatial model type : \n", " Temporal model type : \n", " Parameters:\n", " index : 3.000 \n", " amplitude : 2.50e-12 1 / (cm2 s TeV)\n", " reference (frozen) : 1.000 TeV \n", " \n", " \n" ] } ], "source": [ "dataset_on_off = SpectrumDatasetOnOff.from_spectrum_dataset(\n", " dataset=dataset, acceptance=1, acceptance_off=5\n", ")\n", "dataset_on_off.fake(npred_background=dataset.npred_background())\n", "print(dataset_on_off)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "You can see that off counts are now simulated as well. We now simulate several spectra using the same set of observation conditions." ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "CPU times: user 607 ms, sys: 13 ms, total: 620 ms\n", "Wall time: 621 ms\n" ] } ], "source": [ "%%time\n", "\n", "n_obs = 100\n", "datasets = Datasets()\n", "\n", "for idx in range(n_obs):\n", " dataset_on_off.fake(\n", " random_state=idx, npred_background=dataset.npred_background()\n", " )\n", " dataset_fake = dataset_on_off.copy(name=f\"obs-{idx}\")\n", " dataset_fake.meta_table[\"OBS_ID\"] = [idx]\n", " datasets.append(dataset_fake)" ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "data": { "text/html": [ "Table length=100\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
namecountsbackgroundexcesssqrt_tsnprednpred_backgroundnpred_signalexposure_minexposure_maxlivetimeontimecounts_ratebackground_rateexcess_raten_binsn_fit_binsstat_typestat_sumcounts_offacceptanceacceptance_offalpha
m2 sm2 sss1 / s1 / s1 / s
str6float64float64float64float64float64float64float64float64float64float64float64float64float64float64int64int64str5float64int64float64float64float64
obs-0317.018.400000000000002298.627.08240194504324300.023721792529318.653546302364763281.37017549016446252718170.9728751517719697919.599263600.03600.00.088055555555555550.0051111111111111110.0829444444444444599wstat9.917001109717983929.044.999999999999990.2
obs-1275.022.0253.023.76785365487285302.9527632060656421.582587715901177281.37017549016446252718170.9728751517719697919.599263600.03600.00.07638888888888890.0061111111111111110.0702777777777777799wstat8.818660828556191109.045.00.2
obs-2293.020.6272.425.17110555404655301.839734872925520.469559382761037281.37017549016446252718170.9728751517719697919.599263600.03600.00.081388888888888890.0057222222222222220.0756666666666666699wstat7.3591110029987461039.045.00.2
obs-3280.022.4257.623.982951737405376303.455713902603522.085538412438968281.37017549016446252718170.9728751517719697919.599263600.03600.00.077777777777777780.0062222222222222220.0715555555555555799wstat11.8329834891773161129.045.00.2
obs-4337.020.6316.427.682709945184747302.281284162580720.911108672416212281.37017549016446252718170.9728751517719697919.599263600.03600.00.093611111111111120.0057222222222222220.0878888888888888899wstat17.8671037271713121039.045.00.2
obs-5283.024.400000000000002258.623.727154782347895305.383813845902624.013638355738106281.37017549016446252718170.9728751517719697919.599263600.03600.00.078611111111111120.0067777777777777780.0718333333333333599wstat8.0051961804382581229.044.999999999999990.2
obs-6330.022.400000000000006307.626.889184475727866304.171682312679122.8015068225146281.37017549016446252718170.9728751517719697919.599263600.03600.00.091666666666666660.0062222222222222240.0854444444444444599wstat10.0987664444295581129.044.9999999999999860.2
obs-7283.026.000000000000004257.023.389178133235443307.0218633524487525.651687862284245281.37017549016446252718170.9728751517719697919.599263600.03600.00.078611111111111120.0072222222222222240.0713888888888888999wstat4.3108997791891191309.044.999999999999990.2
obs-8308.023.400000000000002284.625.42049273328333304.835178286128323.46500279596388281.37017549016446252718170.9728751517719697919.599263600.03600.00.085555555555555550.0065000000000000010.0790555555555555699wstat4.1429472135169961179.044.999999999999990.2
.....................................................................
obs-90286.019.000000000000004267.025.131221887043438300.2660027111190318.895827220954608281.37017549016446252718170.9728751517719697919.599263600.03600.00.079444444444444440.0052777777777777790.0741666666666666799wstat6.724388228060566959.044.999999999999990.2
obs-91285.025.200000000000003259.823.67754591931069306.219148597141124.848973106976626281.37017549016446252718170.9728751517719697919.599263600.03600.00.079166666666666660.0070000000000000010.0721666666666666799wstat14.5779408268717591269.044.999999999999990.2
obs-92313.023.6289.425.664935420194176305.05832734099223.688151850827467281.37017549016446252718170.9728751517719697919.599263600.03600.00.086944444444444450.0065555555555555560.0803888888888888899wstat6.3046465683080981189.045.00.2
obs-93302.018.8283.226.123867522605497300.1259513871509618.75577589698646281.37017549016446252718170.9728751517719697919.599263600.03600.00.083888888888888890.0052222222222222230.0786666666666666699wstat5.86503429029057949.045.00.2
obs-94322.022.0300.026.5292481657426303.6642079652339722.294032475069507281.37017549016446252718170.9728751517719697919.599263600.03600.00.089444444444444440.0061111111111111110.0833333333333333399wstat10.0058837608730521109.045.00.2
obs-95305.024.600000000000005280.424.98804632088198305.893020446573724.522844956409255281.37017549016446252718170.9728751517719697919.599263600.03600.00.084722222222222230.00683333333333333450.0778888888888888899wstat5.6277991710462881239.044.999999999999990.2
obs-96301.023.6277.424.969845969421428305.039959264549323.669783774384797281.37017549016446252718170.9728751517719697919.599263600.03600.00.083611111111111110.0065555555555555560.0770555555555555499wstat5.511772031205331189.045.00.2
obs-97290.018.8271.225.417982194454826300.093377655459318.723202165294804281.37017549016446252718170.9728751517719697919.599263600.03600.00.080555555555555560.0052222222222222230.0753333333333333499wstat5.667231272491558949.045.00.2
obs-98301.020.400000000000002280.625.687832964675007301.75614855910520.38597306894053281.37017549016446252718170.9728751517719697919.599263600.03600.00.083611111111111110.0056666666666666670.0779444444444444699wstat7.142141153018121029.044.999999999999990.2
obs-99323.020.8302.226.85707842376871302.476691765042521.106516274878086281.37017549016446252718170.9728751517719697919.599263600.03600.00.089722222222222220.0057777777777777780.0839444444444444599wstat5.0596239298516051049.045.00.2
" ], "text/plain": [ "\n", " name counts background ... acceptance acceptance_off alpha \n", " ... \n", " str6 float64 float64 ... float64 float64 float64\n", "------ ------- ------------------ ... ---------- ------------------ -------\n", " obs-0 317.0 18.400000000000002 ... 9.0 44.99999999999999 0.2\n", " obs-1 275.0 22.0 ... 9.0 45.0 0.2\n", " obs-2 293.0 20.6 ... 9.0 45.0 0.2\n", " obs-3 280.0 22.4 ... 9.0 45.0 0.2\n", " obs-4 337.0 20.6 ... 9.0 45.0 0.2\n", " obs-5 283.0 24.400000000000002 ... 9.0 44.99999999999999 0.2\n", " obs-6 330.0 22.400000000000006 ... 9.0 44.999999999999986 0.2\n", " obs-7 283.0 26.000000000000004 ... 9.0 44.99999999999999 0.2\n", " obs-8 308.0 23.400000000000002 ... 9.0 44.99999999999999 0.2\n", " ... ... ... ... ... ... ...\n", "obs-90 286.0 19.000000000000004 ... 9.0 44.99999999999999 0.2\n", "obs-91 285.0 25.200000000000003 ... 9.0 44.99999999999999 0.2\n", "obs-92 313.0 23.6 ... 9.0 45.0 0.2\n", "obs-93 302.0 18.8 ... 9.0 45.0 0.2\n", "obs-94 322.0 22.0 ... 9.0 45.0 0.2\n", "obs-95 305.0 24.600000000000005 ... 9.0 44.99999999999999 0.2\n", "obs-96 301.0 23.6 ... 9.0 45.0 0.2\n", "obs-97 290.0 18.8 ... 9.0 45.0 0.2\n", "obs-98 301.0 20.400000000000002 ... 9.0 44.99999999999999 0.2\n", "obs-99 323.0 20.8 ... 9.0 45.0 0.2" ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" } ], "source": [ "table = datasets.info_table()\n", "table" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Before moving on to the fit let's have a look at the simulated observations." ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAskAAAEGCAYAAACXYwgRAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4yLjIsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy+WH4yJAAAb4klEQVR4nO3df7AldXnn8fcngD8ibAC5ICLjGAtJkNURR9RFDQQ1MFCiKRNhTSTR7IgRI5auO2pCTLJbNf62FFfESIErAU2UyIZRIURFXQEHHGAQEWTHdWSWgRAFSlcdfPaP0wOX9px7L/f8vvf9qjp1ur/97e6n585zz3P7dPc3VYUkSZKkB/zKuAOQJEmSJo1FsiRJktRikSxJkiS1WCRLkiRJLRbJkiRJUsuu4w6gm3322adWrlw57jCkiXH11VffWVUz446jG/NVejDzVZoec+XrRBbJK1euZOPGjeMOQ5oYSb437hh6MV+lBzNfpekxV756uYUkSZLUYpEsSZIktVgkS5IkSS0WyZIkSVKLRbIkSZLUYpEsSZIktVgkS5IkSS0WyZIkSVKLRbIkSZLUMpEj7qk/K9ddPJDtbFl/3EC2I6k381Vansz9yeeZZEmSJKnFIlmSJElqsUiWJEmSWiySJUmSpBaLZEmSJKnFIlmSJElqsUiWJEmSWiySJUmSpBaLZEmSJKnFIlmSJElqsUiWJEmSWiySJUmSpJZd5+uQ5GzgeGB7VR3atH0SOLjpsifww6pa1WXdLcA9wH3AjqpaPaC4JUmSpKGZt0gGzgHOAD6+s6GqXrZzOsl7gB/Nsf5RVXXnYgOUJEmSRm3eyy2q6nLgrm7LkgT4feD8AcclaYiSHJjki0luTHJDktc37XsnuTTJzc37XuOOVVruzFdpPPq9Jvm5wO1VdXOP5QVckuTqJGvn2lCStUk2Jtl4xx139BmWpHnsAN5YVb8JPAt4bZJDgHXAZVV1EHBZMy9pvMxXaQz6LZJPYu6zyEdU1WHAsXSS+nm9OlbVWVW1uqpWz8zM9BmWpLlU1baquqaZvge4ETgAOAE4t+l2LvDi8UQoaSfzVRqPRRfJSXYFfhf4ZK8+VXVb874duBA4fLH7kzQcSVYCTwOuBParqm3Q+WAG9h1fZJLazFdpdPo5k/x84NtVtbXbwiSPSrLHzmnghcDmPvYnacCS7A58Gjitqu5+COt5eZQ0YuarNFrzFslJzge+DhycZGuSVzWLTqR1qUWSxybZ0MzuB3w1ybXAVcDFVfX5wYUuqR9JdqPzgXteVX2mab49yf7N8v2B7d3W9fIoabTMV2n05n0EXFWd1KP9j7q03QasaaZvBZ7aZ3yShqB5Ms3HgBur6r2zFl0EnAysb94/O4bwJM1ivkrjsZDnJGuZWrnu4oFsZ8v64wayHQ3UEcAfAtcn2dS0vZXOh+2nmm+M/g/we2OKT9IDzFdpDCySpWWoqr4KpMfio0cZi6S5ma/SePT7CDhJkiRpybFIliRJkloskiVJkqQWi2RJkiSpxSJZkiRJarFIliRJkloskiVJkqQWi2RJkiSpxSJZkiRJarFIliRJkloskiVJkqQWi2RJkiSpxSJZkiRJarFIliRJkloskiVJkqQWi2RJkiSpZd4iOcnZSbYn2Tyr7e1JfpBkU/Na02PdY5LclOSWJOsGGbgkSZI0LAs5k3wOcEyX9vdV1armtaG9MMkuwIeAY4FDgJOSHNJPsJIkSdIozFskV9XlwF2L2PbhwC1VdWtV/Qy4ADhhEduRJEmSRmrXPtY9NckrgI3AG6vq31rLDwC+P2t+K/DMXhtLshZYC7BixYo+whq9lesuHsh2tqw/biDbkSRJUn8We+Peh4EnAquAbcB7uvRJl7bqtcGqOquqVlfV6pmZmUWGJUmSJPVvUUVyVd1eVfdV1S+Aj9K5tKJtK3DgrPnHAbctZn+SJEnSKC2qSE6y/6zZlwCbu3T7BnBQkickeRhwInDRYvYnSZIkjdK81yQnOR84EtgnyVbgL4Ejk6yic/nEFuDVTd/HAn9bVWuqakeSU4EvALsAZ1fVDUM5CkmSJGmA5i2Sq+qkLs0f69H3NmDNrPkNwC89Hk6SJEmaZI64J0mSJLVYJEuSJEktFsmSJElSi0WyJEmS1GKRLEmSJLX0Myy1tCAO2y1JkqaNZ5IlSZKkFs8kS9IiDepbkkHwGxstVZP2f3uS8l7D5ZlkSZIkqcUiWZIkSWqxSJYkSZJaLJIlSZKkFotkSZIkqcUiWVqGkpydZHuSzbPa3p7kB0k2Na8144xR0gPMWWn0LJKl5ekc4Jgu7e+rqlXNa8OIY5LU2zmYs9JIWSRLy1BVXQ7cNe44JC2MOSuNnkWypNlOTXJd89XuXuMORtK8zFlpSCySJe30YeCJwCpgG/CeXh2TrE2yMcnGO+64Y1TxSXqwBeWs+SotzrxFco+bBd6V5NvNX68XJtmzx7pbklzf3FCwcZCBSxqsqrq9qu6rql8AHwUOn6PvWVW1uqpWz8zMjC5ISfdbaM6ar9LiLORM8jn88s0ClwKHVtVTgO8Ab5lj/aOaGwpWLy5ESaOQZP9Zsy8BNvfqK2n8zFlpuHadr0NVXZ5kZavtklmzVwAvHWxYkoYpyfnAkcA+SbYCfwkcmWQVUMAW4NVjC1DSg5iz0ujNWyQvwCuBT/ZYVsAlSQr4SFWdNYD9SepTVZ3UpfljIw9E0oKYs9Lo9VUkJ3kbsAM4r0eXI6rqtiT7Apcm+XbzGJtu21oLrAVYsWJFP2FJkiRJfVn00y2SnAwcD7y8qqpbn6q6rXnfDlyINwJJkiRpCiyqSE5yDPBfgBdV1Y979HlUkj12TgMvxJsKJEmSNAUW8gi484GvAwcn2ZrkVcAZwB50LqHYlOTMpu9jk+wcFnM/4KtJrgWuAi6uqs8P5SgkSZKkAVrI0y0WfLNAc3nFmmb6VuCpfUUnSUOwct3F4w5BkjThHHFPkiRJarFIliRJkloskiVJkqQWi2RJkiSpxSJZkiRJarFIliRJkloskiVJkqQWi2RJkiSpZd7BRCRJy8egBlrZsv64gWxH0tzM2eHxTLIkSZLUYpEsSZIktSzryy0G9RWFJEmSlhbPJEuSJEktFsmSJElSi0WyJEmS1GKRLEmSJLVYJEuSJEktFsmSJElSi0WyJEmS1DJvkZzk7CTbk2ye1bZ3kkuT3Ny879Vj3WOS3JTkliTrBhm4JEmSNCwLOZN8DnBMq20dcFlVHQRc1sw/SJJdgA8BxwKHACclOaSvaCVJkqQRmLdIrqrLgbtazScA5zbT5wIv7rLq4cAtVXVrVf0MuKBZT5IkSZpoix2Wer+q2gZQVduS7NulzwHA92fNbwWe2WuDSdYCawFWrFixyLC0lA1qGPEt648byHYkSdLSNcwb99KlrXp1rqqzqmp1Va2emZkZYliSJEnS3BZbJN+eZH+A5n17lz5bgQNnzT8OuG2R+5MkSZJGZrFF8kXAyc30ycBnu/T5BnBQkickeRhwYrOeJEmSNNEW8gi484GvAwcn2ZrkVcB64AVJbgZe0MyT5LFJNgBU1Q7gVOALwI3Ap6rqhuEchiRJkjQ48964V1Un9Vh0dJe+twFrZs1vADYsOjpJkiRpDBxxT5IkSWqxSJYkSZJaLJIlSZKkFotkSZIkqcUiWZIkSWqxSJYkSZJaLJKlZSjJ2Um2J9k8q23vJJcmubl532ucMUp6gDkrjZ5FsrQ8nQMc02pbB1xWVQcBlzXzkibDOZiz0khZJEvLUFVdDtzVaj4BOLeZPhd48UiDktSTOSuN3rwj7klaNvarqm0AVbUtyb69OiZZC6wFWLFixYjCk9SyoJw1X7UQK9dd3Pc2tqw/bgCRTA7PJEt6yKrqrKpaXVWrZ2Zmxh2OpDmYr9LiWCRL2un2JPsDNO/bxxyPpLmZs9IQWSRL2uki4ORm+mTgs2OMRdL8zFlpiCySpWUoyfnA14GDk2xN8ipgPfCCJDcDL2jmJU0Ac1YaPW/ck5ahqjqpx6KjRxqIpAUxZ6XR80yyJEmS1GKRLEmSJLVYJEuSJEktiy6SkxycZNOs191JTmv1OTLJj2b1Ob3/kCVJkqThWvSNe1V1E7AKIMkuwA+AC7t0/UpVHb/Y/UiSJEmjNqjLLY4GvltV3xvQ9iRJkqSxGVSRfCJwfo9lz05ybZLPJXlyrw0kWZtkY5KNd9xxx4DCkiRJkh66vovkJA8DXgT8fZfF1wCPr6qnAh8E/rHXdhxbXpIkSZNiEGeSjwWuqarb2wuq6u6qureZ3gDslmSfAexTkiRJGppBFMkn0eNSiySPSZJm+vBmf/86gH1KkiRJQ9PXsNRJfpXOePGvntV2CkBVnQm8FHhNkh3AT4ATq6r62ackSZI0bH0VyVX1Y+DRrbYzZ02fAZzRzz4kSZKkUXPEPUmSJKnFIlmSJElqsUiWJEmSWiySJUmSpBaLZEmSJKnFIlmSJElqsUiWJEmSWiySJUmSpBaLZEmSJKnFIlmSJElq6WtYakmSpGmwct3F4w5BU8YzyZIkSVKLRbIkSZLUYpEsSZIktVgkS5IkSS0WyZIkSVKLRbIkSZLUYpEsSZIktfRVJCfZkuT6JJuSbOyyPEk+kOSWJNclOayf/UmSJEmjMIjBRI6qqjt7LDsWOKh5PRP4cPMuSQ+ZgwFI08N81bQb9uUWJwAfr44rgD2T7D/kfUqSJEl96fdMcgGXJCngI1V1Vmv5AcD3Z81vbdq2tTeUZC2wFmDFihV9hiX1NqizG1vWHzeQ7UiSpMnT75nkI6rqMDqXVbw2yfNay9Nlneq2oao6q6pWV9XqmZmZPsOSJEmSFq+vIrmqbmvetwMXAoe3umwFDpw1/zjgtn72KUmSJA3boovkJI9KssfOaeCFwOZWt4uAVzRPuXgW8KOq+qVLLSRJkqRJ0s81yfsBFybZuZ2/q6rPJzkFoKrOBDYAa4BbgB8Df9xfuJIkSdLwLbpIrqpbgad2aT9z1nQBr13sPiSNXpItwD3AfcCOqlo93ogkzcWclYZjEM9JlrT0zPX8c0mTx5yVBsxhqSVJkqQWzyRLapvv+ec+11zz8nnkIzVnzpqv0uJ4JllS23zPP/e55tJkmTNnzVdpcSySJT3IAp5/LmmCmLPScFgkS7rfAp9/LmlCmLPS8HhNsqTZuj7/fLwhSZqDOSsNiUWypPv1ev65pMlkzkrD4+UWkiRJUotFsiRJktRikSxJkiS1eE2yJGliOSiJND0Gla+D0m/eeyZZkiRJarFIliRJkloskiVJkqQWi2RJkiSpZWpv3Ju0i8MHYSke01LmDUWSJC1dnkmWJEmSWhZdJCc5MMkXk9yY5IYkr+/S58gkP0qyqXmd3l+4kiRJ0vD1c7nFDuCNVXVNkj2Aq5NcWlXfavX7SlUd38d+JEmSpJFa9JnkqtpWVdc00/cANwIHDCowSZIkaVwGcuNekpXA04Aruyx+dpJrgduAN1XVDYPYpyRJC+WNtpIeqr6L5CS7A58GTququ1uLrwEeX1X3JlkD/CNwUI/trAXWAqxYsaLfsCRJkqRF6+vpFkl2o1Mgn1dVn2kvr6q7q+reZnoDsFuSfbptq6rOqqrVVbV6Zmamn7AkSZKkvvTzdIsAHwNurKr39ujzmKYfSQ5v9vevi92nJEmSNAr9XG5xBPCHwPVJNjVtbwVWAFTVmcBLgdck2QH8BDixqqqPfUqSJElDt+giuaq+CmSePmcAZyx2H5IkSdI4OOKeJEmS1GKRLEmSJLVYJEuSJEktAxlMRJIkLQ2DGnhFmnaeSZYkSZJaLJIlSZKkFotkSZIkqcUiWZIkSWqxSJYkSZJaLJIlSZKkFotkSZIkqcUiWZIkSWpxMBFJQ+fgBJKkaeOZZEmSJKnFIlmSJElqsUiWJEmSWiySJUmSpBaLZEmSJKnFIlmSJElq6atITnJMkpuS3JJkXZflSfKBZvl1SQ7rZ3+Shm++vJY0OcxXaXgWXSQn2QX4EHAscAhwUpJDWt2OBQ5qXmuBDy92f5KGb4F5LWkCmK/ScPVzJvlw4JaqurWqfgZcAJzQ6nMC8PHquALYM8n+fexT0nAtJK8lTQbzVRqifkbcOwD4/qz5rcAzF9DnAGBbe2NJ1tI52wxwb5Kb+ohtVPYB7hx3EH2Y5vinOXaYFX/esaD+jx9mMLMsJK+nNV8Hadr//y3Wsj/uZZKv0/ZzNt7hmrZ4oYm533ztp0hOl7ZaRJ9OY9VZwFl9xDNySTZW1epxx7FY0xz/NMcOEx3/gnJ2GvN1kCb45zdUHvfEGUq+TvDxdmW8wzVt8cLgYu7ncoutwIGz5h8H3LaIPpImhzkrTQ/zVRqiforkbwAHJXlCkocBJwIXtfpcBLyiecrFs4AfVdUvXWohaWIsJK8lTQbzVRqiRV9uUVU7kpwKfAHYBTi7qm5Ickqz/ExgA7AGuAX4MfDH/Yc8Uab96+Zpjn+aY4cJjb9XXo85rEk0kT+/EfC4J8gQ83Uij3cOxjtc0xYvDCjmVHW9RFiSJElathxxT5IkSWqxSJYkSZJaLJJ7SHJgki8muTHJDUle37SvSnJFkk1JNiY5fNY6b2mGBr0pye+ML3pI8ogkVyW5ton/r5r2vZNcmuTm5n2vWetMQ/zvSvLtZpjzC5PsOWudiYi/V+yzlr8pSSXZZ1bbRMSujiRnJ9meZPOstqnInX70OO6Jz7l+dDvmWcuWXK72+mybtXyijnmueJO8ronphiTvnIR454p5UuuHaasXRlofVJWvLi9gf+CwZnoP4Dt0hv28BDi2aV8DfKmZPgS4Fng48ATgu8AuY4w/wO7N9G7AlcCzgHcC65r2dcA7piz+FwK7Nu3vmMT4e8XezB9I5yab7wH7TFrsvu7/GT4POAzYPKttKnJnCMc98Tk36GNu2pdkrvb6bJvUY+4VL3AU8M/Aw5tl+05CvPPEPJH1Q6/PrEn9nTdHvAP/XeWZ5B6qaltVXdNM3wPcSGd0owL+XdPt13jgmZQnABdU1U+r6n/TeaLH4YxJddzbzO7WvIpOnOc27ecCL26mpyL+qrqkqnY07VfQeS4oTFD8c/zbA7wPeDMPfuD/xMSujqq6HLir1TwVudOPbsc9DTnXjx4/a1iiuTrHZxtM4DHPEe9rgPVV9dNm2fZJiHeemCeyfpi2emGU9YFF8gIkWQk8jc5fK6cB70ryfeDdwFuabr2G4B6bJLsk2QRsBy6tqiuB/ap5VnXzvm/TfVrin+2VwOea6YmKv1vsSV4E/KCqrm11n6jY1dPU5M4QTWzODdJyydXZn23TcMytz+InAc9NcmWSLyd5RtNtYuKF6akfpq1eGFV9YJE8jyS7A58GTququ+n89fqGqjoQeAPwsZ1du6w+1ufrVdV9VbWKzl9Thyc5dI7uUxV/krcBO4DzdjZ128Two+yuS+xPAd4GnN6l+0TFrodsWfz8Jj3nBiXJr7IMcnX2Zxudn+tEH3OXz+Jdgb3ofM3+n4FPJQkTEi9MV/0wbfXCqOoDi+Q5JNmNzn/w86rqM03zycDO6b/ngVP2Ezs8aFX9EPgScAxwe5L9AZr3nV9RTUv8JDkZOB54eTUXHDGh8c+K/QQ610Jdm2QLnfiuSfIYJjR2/ZKpy51BmaacG4AnssRztctn20Qfc4/P4q3AZ5qv3q8CfgHsMwnxwvTWD9NWLwy9Puh2obKv+y8M/zjw/lb7jcCRzfTRwNXN9JN58IXhtzLeG99mgD2b6UcCX2n+47yLB1+I/84pi/8Y4FvATKv/xMTfK/ZWny08cGPMxMTu60E/o5U8+Aa2qcidIRz3xOfcoI+5tWxJ5Wqvz7ZJPeZe8QKnAH/dTD+JztfpGXe888Q8kfVDr8+sSf2dN0e8A/9dNbL/NNP2Ap5D53T8dcCm5rWmab+6+Qe/Enj6rHXeRueuyZto7mAdY/xPAb7ZxL8ZOL1pfzRwGXBz8773lMV/S/PLcOfP5MxJi79X7K0+W2g+hCYpdl/3/zzOB7YBP6dzFuJV05I7Qzjuic+5QR9za/mSytVen22Tesy94gUeBnyi+R17DfDbkxDvPDFPZP3Q6zNrUn/nzRHvwH9XOSy1JEmS1OI1yZIkSVKLRbIkSZLUYpEsSZIktVgkS5IkSS0WyZIkSVKLRfISluQxSS5I8t0k30qyIcmTBrj9I5P8h0FtT1qqhp2LzT4ecj4meXKSf0nynSQ3J/mLZtQykjw8yT8n2ZTkZUmem+SGZv6Rg4xdkiaRRfIS1XzQXQh8qaqeWFWHAG8F9hvgbo4ELJKlOYwoF+Eh5mNT6F4ErK+qJwFPbdb/06bL04DdqmpVVX0SeDnw7mb+JwONXJImkEXy0nUU8POqOnNnQ1VtAr6a5F1JNie5PsnL4P6zUP+0s2+SM5L8UTO9JclfJbmmWec3kqykM+LRG5ozS89N8nvNdq9NcvkIj1WaZF1zsaq+ko5x5eN/BL5WVZc0Mf0YOBVYl2RfOgM1rGq292rg94HTk5w3+H8iafok+YMkVzU58pEkz0xyXZJHJHlU883LoUl2SfLuJl+vS/K6Zv2nJ/lykquTfGHWENB/1nzjdF2SC5q232r2synJN5PsMc5jXy52HXcAGppD6Yzs0/a7wCo6Z432Ab6xwIL2zqo6LMmfAm+qqj9JciZwb1W9GyDJ9cDvVNUPkuw5mMOQpl6vXITx5uOT23FV1XeT7A78P+BPmm0f32zv2cA/VdU/LCA+aUlL8pvAy4AjqurnSf47cDCdb2f+K53hkj9RVZuTvIbOcMhPq6odSfZOshvwQeCEqrqj+QP5vwGvpDME9BOq6qezcvdNwGur6muzclRD5pnk5ec5wPlVdV9V3Q58GXjGAtb7TPN+NbCyR5+vAeck+U/ALv0GKi0D48zH0Bk6txuHYpXmdjTwdDp/2G5q5n8d+GvgBcBq4J1N3+fTGSJ5B0BV3UWnoD4UuLRZ/8+BxzX9rwPOS/IHwI6m7WvAe5P8GbDnzm1puCySl64b6CRwW3r038GD/z88orX8p837ffT4BqKqTqGT6AcCm5I8esHRSktXr1yE8ebjDXQ+yB8IJvl1Omej7+kRl6SOAOc21+ivqqqDq+rtwN7A7sAePJC33f4gDXDDrPX/fVW9sFl2HPAhOr83rk6ya1Wtp/PtziOBK5L8xlCPToBF8lL2L8DDm7NIACR5BvBvwMuaa6RmgOcBVwHfAw5p7mj/NTp/Fc/nHjq/CHZu/4lVdWVVnQ7cSefDWVruuuZikt8CLmd8+Xge8Jwkz2/6PxL4AA+c/ZLU22XAS5vr92kuoXg8cBbwF3Ty6x1N30uAU5LsurMvcBMw01zGRJLd0nnazK8AB1bVF4E3A3sCuzf5fH1VvQPYCFgkj4DXJC9RVVVJXgK8P8k6OtcvbQFOo/NX7rV0/rJ9c1X9X4Akn6LzNc/NwDcXsJv/CfxDkhOA19G5aeggOn8hX9bsQ1rW5snFy4FnM4Z8rKqfNH0/mORDdC7J+B/AGf0dsbT0VdW3kvw5cElT2P4c+Cywo6r+LskuwP9K8tvA3wJPAq5L8nPgo1V1RpKXAh9o/hDeFXg/8B3gE01bgPdV1Q+T/E2So+h8e/Qt4HMjPuRlKVVeeiZJkiTN5uUWkiRJUotFsiRJktRikSxJkiS1WCRLkiRJLRbJkiRJUotFsiRJktRikSxJkiS1/H9TKYyhzV/B1AAAAABJRU5ErkJggg==\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "fix, axes = plt.subplots(1, 3, figsize=(12, 4))\n", "axes[0].hist(table[\"counts\"])\n", "axes[0].set_xlabel(\"Counts\")\n", "axes[1].hist(table[\"counts_off\"])\n", "axes[1].set_xlabel(\"Counts Off\")\n", "axes[2].hist(table[\"excess\"])\n", "axes[2].set_xlabel(\"excess\");" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now, we fit each simulated spectrum individually " ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "CPU times: user 11.5 s, sys: 86.2 ms, total: 11.6 s\n", "Wall time: 11.7 s\n" ] } ], "source": [ "%%time\n", "results = []\n", "\n", "for dataset in datasets:\n", " dataset.models = model.copy()\n", " fit = Fit([dataset])\n", " result = fit.optimize()\n", " results.append(\n", " {\n", " \"index\": result.parameters[\"index\"].value,\n", " \"amplitude\": result.parameters[\"amplitude\"].value,\n", " }\n", " )" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We take a look at the distribution of the fitted indices. This matches very well with the spectrum that we initially injected." ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "index: 3.0036666673409944 += 0.08075110145690628\n" ] }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAXAAAAD4CAYAAAD1jb0+AAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4yLjIsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy+WH4yJAAAL4UlEQVR4nO3df4zk9V3H8eerQKMRKjS34LWCZxqgvRhLyYlojaHBq9B/rv5KJKYQgjlNWkNpm0j4p2pigjGWxNRozkKgsdY0FlJM/HW5lKBpi10qpeCVg9AfYi/cISgQjXr07R/zvea6zN3M7czO3Hvv+Ug28+u7O+9Pdu+5s9+Z71yqCklSP69Z9gCSpPUx4JLUlAGXpKYMuCQ1ZcAlqakzF3lnW7ZsqW3bti3yLqX5eOKJ0emlly53Dp2WHn744eeqamXt9QsN+LZt21hdXV3kXUrzcdVVo9MHHljmFDpNJfnGuOvdhSJJTRlwSWrKgEtSUwZckpoy4JLUlAGXpKYMuCQ1ZcAlqSkDLklNLfRITOlUdcfeAye8/Rdf+G8A/nLCdutxy85L5v41dXrwEbgkNWXAJakpAy5JTRlwSWrKgEtSUwZckpoy4JLUlAGXpKYMuCQ1ZcAlqSkDLklNGXBJasqAS1JTBlySmjLgktSUAZekpgy4JDVlwCWpKQMuSU0ZcElqyoBLUlMGXJKaMuCS1NTEgCe5MMlnk+xP8niSm4frX59kb5Inh9PzNn5cSdJR0zwCPwJ8sKreAlwJvDfJduBWYF9VXQzsGy5LkhZkYsCr6mBVfWk4/xKwH3gjsAu4Z9jsHuDdGzWkJOnVTmofeJJtwNuAh4ALquogjCIPnH+cz9mdZDXJ6uHDh2ebVpL0HVMHPMnZwKeB91fVi9N+XlXtqaodVbVjZWVlPTNKksaYKuBJzmIU709U1b3D1c8m2TrcvhU4tDEjSpLGmeZVKAHuBPZX1UeOuel+4Ibh/A3AZ+Y/niTpeM6cYpu3A+8BvpLkkeG624DbgU8luQn4JvBLGzOiJGmciQGvqn8Ecpybr57vOJKkaXkkpiQ1ZcAlqSkDLklNTfMkprQwd+w9sOwRFm5Za75l5yVLuV/Nj4/AJakpAy5JTRlwSWrKgEtSUwZckpoy4JLUlAGXpKYMuCQ1ZcAlqSkDLklNGXBJasqAS1JTBlySmjLgktSUAZekpgy4JDVlwCWpKQMuSU0ZcElqyoBLUlMGXJKaMuCS1JQBl6SmDLgkNWXAJakpAy5JTRlwSWrKgEtSUwZckpoy4JLUlAGXpKYMuCQ1NTHgSe5KcijJY8dc91tJ/i3JI8PHuzZ2TEnSWtM8Ar8buGbM9XdU1WXDx1/PdyxJ0iQTA15VDwLPL2AWSdJJOHOGz31fkuuBVeCDVfXCuI2S7AZ2A1x00UUz3J0W5Y69B5Y9gqQprPdJzD8G3gRcBhwE/uB4G1bVnqraUVU7VlZW1nl3kqS11hXwqnq2ql6pqm8DfwpcMd+xJEmTrCvgSbYec/HngMeOt60kaWNM3Aee5JPAVcCWJM8AHwauSnIZUMDXgV/bwBklSWNMDHhVXTfm6js3YBZJ0knwSExJasqAS1JTBlySmprlQB5JjS3zgK1bdl6ytPveTHwELklNGXBJasqAS1JTBlySmjLgktSUAZekpgy4JDVlwCWpKQMuSU0ZcElqyoBLUlMGXJKaMuCS1JQBl6SmDLgkNWXAJakpAy5JTRlwSWrKgEtSUwZckpoy4JLUlAGXpKYMuCQ1ZcAlqSkDLklNGXBJasqAS1JTBlySmjLgktSUAZekpgy4JDVlwCWpqYkBT3JXkkNJHjvmutcn2ZvkyeH0vI0dU5K01jSPwO8Grllz3a3Avqq6GNg3XJYkLdDEgFfVg8Dza67eBdwznL8HePec55IkTbDefeAXVNVBgOH0/ONtmGR3ktUkq4cPH17n3UmS1trwJzGrak9V7aiqHSsrKxt9d5J02lhvwJ9NshVgOD00v5EkSdNYb8DvB24Yzt8AfGY+40iSpjXNywg/CXweuDTJM0luAm4HdiZ5Etg5XJYkLdCZkzaoquuOc9PVc55FknQSPBJTkpoy4JLUlAGXpKYMuCQ1ZcAlqSkDLklNGXBJasqAS1JTBlySmjLgktSUAZekpgy4JDVlwCWpKQMuSU0ZcElqyoBLUlMGXJKaMuCS1JQBl6SmDLgkNWXAJakpAy5JTRlwSWrKgEtSUwZckpoy4JLUlAGXpKYMuCQ1ZcAlqSkDLklNGXBJasqAS1JTBlySmjLgktSUAZekpgy4JDV15iyfnOTrwEvAK8CRqtoxj6EkSZPNFPDBO6rquTl8HUnSSXAXiiQ1NWvAC/j7JA8n2T1ugyS7k6wmWT18+PCMdydJOmrWgL+9qi4HrgXem+Sn125QVXuqakdV7VhZWZnx7iRJR80U8Kr61nB6CLgPuGIeQ0mSJlt3wJN8X5Jzjp4H3gk8Nq/BJEknNsurUC4A7kty9Ov8eVX97VymkiRNtO6AV9XTwFvnOIsk6ST4MkJJasqAS1JTBlySmprHofTaIHfsPbDsESSdwnwELklNGXBJasqAS1JTBlySmjLgktSUAZekpgy4JDVlwCWpKQMuSU0ZcElqyoBLUlMGXJKaMuCS1JQBl6SmDLgkNWXAJakpAy5JTfk/8kzB/xlH2hyW+W/5lp2XzP1r+ghckpoy4JLUlAGXpKYMuCQ1ZcAlqSkDLklNGXBJasqAS1JTbQ7k8WAaafPw3/N8+Ahckpoy4JLUlAGXpKYMuCQ1ZcAlqamZAp7kmiRPJHkqya3zGkqSNNm6A57kDOCPgGuB7cB1SbbPazBJ0onN8gj8CuCpqnq6qv4X+Atg13zGkiRNMsuBPG8E/vWYy88AP752oyS7gd3DxZeTPDHDfZ4qtgDPLXuIJTkt1/6B0ckW3nnpabf2wWn5fR/MZe0fmO3Tf2jclbMEPGOuq1ddUbUH2DPD/ZxykqxW1Y5lz7EMrt21n25O5bXPsgvlGeDCYy7/IPCt2caRJE1rloB/Ebg4yQ8neS3wy8D98xlLkjTJunehVNWRJO8D/g44A7irqh6f22Sntk21S+gkufbTk2s/BaXqVbutJUkNeCSmJDVlwCWpKQM+RpILk3w2yf4kjye5ecw235/kr5J8edjmxmXMuhGmXP95Se5L8miSf0ryI8uYdd6SfM+wnqPf198es02S/OHwFhKPJrl8GbPO25Rrf3OSzyf5nyQfWsacG2HKtf/K8P1+NMnnkrx1GbN+l6ryY80HsBW4fDh/DnAA2L5mm9uA3xvOrwDPA69d9uwLXP/vAx8ezr8Z2Lfsuee09gBnD+fPAh4CrlyzzbuAvxm2vRJ4aNlzL3Dt5wM/Bvwu8KFlz7zgtf8kcN5w/tpT4fvuI/AxqupgVX1pOP8SsJ/RkafftRlwTpIAZzMK+JGFDrpBplz/dmDfsM1XgW1JLljooBugRl4eLp41fKx9pn8X8PFh2y8A5ybZusg5N8I0a6+qQ1X1ReD/Fj3fRppy7Z+rqheGi19gdOzLUhnwCZJsA97G6DfysT4KvIXRwUtfAW6uqm8vdLgFOMH6vwz8/LDNFYwO9V36D/Q8JDkjySPAIWBvVa1d+7i3kVj7C66lKda+aZ3k2m9i9FfYUhnwE0hyNvBp4P1V9eKam38WeAR4A3AZ8NEkr1vwiBtqwvpvB84bfuB/A/hnNs9fIK9U1WWMfiFdMWb//lRvI9HRFGvftKZde5J3MAr4by5yvnEM+HEkOYtRvD5RVfeO2eRG4N7hT6+ngK8x2he8KUxaf1W9WFU3Dj/w1zN6HuBrCx5zQ1XVfwAPANesuWnTv43ECda+6Z1o7Ul+FPgYsKuq/n3Bo72KAR9j2K99J7C/qj5ynM2+CVw9bH8BcCnw9GIm3FjTrD/JucNbKAD8KvDgmEfp7SRZSXLucP57gZ8Bvrpms/uB64dXo1wJ/GdVHVzwqHM35do3pWnWnuQi4F7gPVV1YPFTvppHYo6R5KeAf2C0b/vofu3bgIsAqupPkrwBuJvRKzYC3F5Vf7b4aedvyvX/BPBx4BXgX4CbjnmCp63hEdY9jN4e4jXAp6rqd5L8Onxn7WH0HMg1wH8BN1bV6rJmnpcp1/4DwCrwOkY/Gy8zeoVS61/eU679Y8AvAN8YPu1ILfldCg24JDXlLhRJasqAS1JTBlySmjLgktSUAZekpgy4JDVlwCWpqf8HIs2lwUKgdC8AAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "index = np.array([_[\"index\"] for _ in results])\n", "plt.hist(index, bins=10, alpha=0.5)\n", "plt.axvline(x=model_simu.parameters[\"index\"].value, color=\"red\")\n", "print(f\"index: {index.mean()} += {index.std()}\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Exercises\n", "\n", "* Change the observation time to something longer or shorter. Does the observation and spectrum results change as you expected?\n", "* Change the spectral model, e.g. add a cutoff at 5 TeV, or put a steep-spectrum source with spectral index of 4.0\n", "* Simulate spectra with the spectral model we just defined. How much observation duration do you need to get back the injected parameters?" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.7.0" }, "nbsphinx": { "orphan": true } }, "nbformat": 4, "nbformat_minor": 4 }