.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "tutorials/api/theta_square_plot.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end <sphx_glr_download_tutorials_api_theta_square_plot.py>` to download the full example code. or to run this example in your browser via Binder .. rst-class:: sphx-glr-example-title .. _sphx_glr_tutorials_api_theta_square_plot.py: Make a theta-square plot ======================== This tutorial explains how to make such a plot, that is the distribution of event counts as a function of the squared angular distance, to a test position. .. GENERATED FROM PYTHON SOURCE LINES 13-16 Setup ----- .. GENERATED FROM PYTHON SOURCE LINES 16-27 .. code-block:: Python # %matplotlib inline import matplotlib.pyplot as plt from astropy.coordinates import SkyCoord from astropy import units as u from gammapy.data import DataStore from gammapy.maps import MapAxis from gammapy.makers.utils import make_theta_squared_table from gammapy.visualization import plot_theta_squared_table .. GENERATED FROM PYTHON SOURCE LINES 28-31 Check setup ----------- .. GENERATED FROM PYTHON SOURCE LINES 31-37 .. code-block:: Python from gammapy.utils.check import check_tutorials_setup check_tutorials_setup() .. rst-class:: sphx-glr-script-out .. code-block:: none System: python_executable : /home/runner/work/gammapy-docs/gammapy-docs/gammapy/.tox/build_docs/bin/python python_version : 3.9.22 machine : x86_64 system : Linux Gammapy package: version : 2.0.dev1326+g6403c0061 path : /home/runner/work/gammapy-docs/gammapy-docs/gammapy/.tox/build_docs/lib/python3.9/site-packages/gammapy Other packages: numpy : 1.26.4 scipy : 1.13.1 astropy : 6.0.1 regions : 0.8 click : 8.1.8 yaml : 6.0.2 IPython : 8.18.1 jupyterlab : not installed matplotlib : 3.9.4 pandas : not installed healpy : 1.17.3 iminuit : 2.31.1 sherpa : 4.16.1 naima : 0.10.0 emcee : 3.1.6 corner : 2.2.3 ray : 2.46.0 Gammapy environment variables: GAMMAPY_DATA : /home/runner/work/gammapy-docs/gammapy-docs/gammapy-datasets/dev .. GENERATED FROM PYTHON SOURCE LINES 38-43 Get some data ------------- Some data taken on the Crab by H.E.S.S. are used. .. GENERATED FROM PYTHON SOURCE LINES 43-48 .. code-block:: Python data_store = DataStore.from_dir("$GAMMAPY_DATA/hess-dl3-dr1") observations = data_store.get_observations([23523, 23526]) .. GENERATED FROM PYTHON SOURCE LINES 49-52 Define a test position ----------------------- .. GENERATED FROM PYTHON SOURCE LINES 52-57 .. code-block:: Python position = SkyCoord.from_name("crab") print(position) .. rst-class:: sphx-glr-script-out .. code-block:: none <SkyCoord (ICRS): (ra, dec) in deg (83.6324, 22.0174)> .. GENERATED FROM PYTHON SOURCE LINES 58-61 Creation of the theta2 plot --------------------------- .. GENERATED FROM PYTHON SOURCE LINES 61-74 .. code-block:: Python theta2_axis = MapAxis.from_bounds(0, 0.2, nbin=20, interp="lin", unit="deg2") theta2_table = make_theta_squared_table( observations=observations, position=position, theta_squared_axis=theta2_axis, ) plt.figure(figsize=(10, 5)) plot_theta_squared_table(theta2_table) plt.show() .. image-sg:: /tutorials/api/images/sphx_glr_theta_square_plot_001.png :alt: theta square plot :srcset: /tutorials/api/images/sphx_glr_theta_square_plot_001.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 75-80 Making a theta2 plot for a given energy range --------------------------------------------- with the function `make_theta_squared_table`, one can also select a fixed energy range. .. GENERATED FROM PYTHON SOURCE LINES 80-92 .. code-block:: Python theta2_table_en = make_theta_squared_table( observations=observations, position=position, theta_squared_axis=theta2_axis, energy_edges=[1.2, 11] * u.TeV, ) plt.figure(figsize=(10, 5)) plot_theta_squared_table(theta2_table_en) plt.show() .. image-sg:: /tutorials/api/images/sphx_glr_theta_square_plot_002.png :alt: theta square plot :srcset: /tutorials/api/images/sphx_glr_theta_square_plot_002.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 93-104 Statistical significance of a detection --------------------------------------- To get the significance of a signal, the usual method consists of using the reflected background method (see the maker tutorial: :doc:`/user-guide/makers/reflected`) to compute the WStat statistics (see :ref:`wstat`, :doc:`/user-guide/stats/fit_statistics`). This is the well-known method of Li&Ma [LiMa1983]_ using ON and OFF regions. The following tutorials show how to get an excess significance: * :doc:`/tutorials/analysis-1d/spectral_analysis` * :doc:`/tutorials/analysis-1d/extended_source_spectral_analysis` .. _sphx_glr_download_tutorials_api_theta_square_plot.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: binder-badge .. image:: images/binder_badge_logo.svg :target: https://mybinder.org/v2/gh/gammapy/gammapy-webpage/main?urlpath=lab/tree/notebooks/dev/tutorials/api/theta_square_plot.ipynb :alt: Launch binder :width: 150 px .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: theta_square_plot.ipynb <theta_square_plot.ipynb>` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: theta_square_plot.py <theta_square_plot.py>` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: theta_square_plot.zip <theta_square_plot.zip>` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery <https://sphinx-gallery.github.io>`_