.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/tutorials/ctapipe_handson.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code. .. rst-class:: sphx-glr-example-title .. _sphx_glr_auto_examples_tutorials_ctapipe_handson.py: Getting Started with ctapipe ============================ This hands-on was presented at the Paris CTAO Consoritum meeting (K. Kosack) .. GENERATED FROM PYTHON SOURCE LINES 12-15 Part 1: load and loop over data ------------------------------- .. GENERATED FROM PYTHON SOURCE LINES 15-31 .. code-block:: Python :lineno-start: 16 import glob import numpy as np import pandas as pd from ipywidgets import interact from matplotlib import pyplot as plt from ctapipe import utils from ctapipe.calib import CameraCalibrator from ctapipe.image import hillas_parameters, tailcuts_clean from ctapipe.io import EventSource, HDF5TableWriter from ctapipe.visualization import CameraDisplay # %matplotlib inline .. GENERATED FROM PYTHON SOURCE LINES 32-34 .. code-block:: Python :lineno-start: 32 path = utils.get_dataset_path("gamma_prod5.simtel.zst") .. GENERATED FROM PYTHON SOURCE LINES 35-40 .. code-block:: Python :lineno-start: 35 source = EventSource(path, max_events=5) for event in source: print(event.count, event.index.event_id, event.simulation.shower.energy) .. rst-class:: sphx-glr-script-out .. code-block:: none 0 4009 0.07287972420454025 TeV 1 5101 1.7454713582992554 TeV 2 5103 1.7454713582992554 TeV 3 5104 1.7454713582992554 TeV 4 5105 1.7454713582992554 TeV .. GENERATED FROM PYTHON SOURCE LINES 41-43 .. code-block:: Python :lineno-start: 41 event .. rst-class:: sphx-glr-script-out .. code-block:: none ctapipe.containers.ArrayEventContainer: index.*: event indexing information with default None r0.*: Raw Data with default None r1.*: R1 Calibrated Data with default None dl0.*: DL0 Data Volume Reduced Data with default None dl1.*: DL1 Calibrated image with default None dl2.*: DL2 reconstruction info with default None simulation.*: Simulated Event Information with default None with type trigger.*: central trigger information with default None count: number of events processed with default 0 monitoring.*: container for event-wise monitoring data (MON) with default None muon.*: Container for muon analysis results with default None .. GENERATED FROM PYTHON SOURCE LINES 44-46 .. code-block:: Python :lineno-start: 44 event.r1 .. rst-class:: sphx-glr-script-out .. code-block:: none ctapipe.containers.R1Container: tel[*]: map of tel_id to R1CameraContainer with default None .. GENERATED FROM PYTHON SOURCE LINES 47-50 .. code-block:: Python :lineno-start: 47 for event in EventSource(path, max_events=5): print(event.count, event.r1.tel.keys()) .. rst-class:: sphx-glr-script-out .. code-block:: none 0 dict_keys([9, 14, 104]) 1 dict_keys([26, 27, 29, 47, 53, 59, 69, 73, 121, 122, 124, 148, 162, 166]) 2 dict_keys([82, 92, 175, 179]) 3 dict_keys([17, 26, 27, 29, 43, 47, 53, 57, 69, 112, 121, 122, 124, 125, 128, 140, 144, 148, 162, 166]) 4 dict_keys([3, 4, 8, 9, 14, 15, 23, 24, 25, 29, 32, 36, 42, 46, 52, 56, 103, 104, 109, 110, 118, 119, 120, 126, 127, 129, 130, 137, 139, 143, 161]) .. GENERATED FROM PYTHON SOURCE LINES 51-53 .. code-block:: Python :lineno-start: 51 event.r0.tel[3] .. rst-class:: sphx-glr-script-out .. code-block:: none ctapipe.containers.R0CameraContainer: waveform: numpy array containing ADC samples: (n_channels, n_pixels, n_samples) with default None .. GENERATED FROM PYTHON SOURCE LINES 54-56 .. code-block:: Python :lineno-start: 54 r0tel = event.r0.tel[3] .. GENERATED FROM PYTHON SOURCE LINES 57-59 .. code-block:: Python :lineno-start: 57 r0tel.waveform .. rst-class:: sphx-glr-script-out .. code-block:: none array([[[432, 413, 395, ..., 382, 387, 382], [419, 420, 428, ..., 392, 396, 405], [381, 392, 413, ..., 375, 387, 387], ..., [391, 384, 385, ..., 379, 395, 404], [370, 376, 388, ..., 391, 383, 383], [393, 411, 413, ..., 406, 392, 396]], [[413, 396, 398, ..., 397, 401, 395], [392, 407, 403, ..., 403, 388, 399], [405, 403, 393, ..., 393, 400, 409], ..., [399, 406, 407, ..., 391, 405, 397], [401, 399, 396, ..., 395, 400, 388], [399, 393, 403, ..., 395, 407, 401]]], shape=(2, 1855, 40), dtype=uint16) .. GENERATED FROM PYTHON SOURCE LINES 60-63 .. code-block:: Python :lineno-start: 60 r0tel.waveform.shape .. rst-class:: sphx-glr-script-out .. code-block:: none (2, 1855, 40) .. GENERATED FROM PYTHON SOURCE LINES 64-67 note that this is (:math:`N_{channels}`, :math:`N_{pixels}`, :math:`N_{samples}`) .. GENERATED FROM PYTHON SOURCE LINES 67-70 .. code-block:: Python :lineno-start: 68 plt.pcolormesh(r0tel.waveform[0]) .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_001.png :alt: ctapipe handson :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 71-74 .. code-block:: Python :lineno-start: 71 brightest_pixel = np.argmax(r0tel.waveform[0].sum(axis=1)) print(f"pixel {brightest_pixel} has sum {r0tel.waveform[0,1535].sum()}") .. rst-class:: sphx-glr-script-out .. code-block:: none pixel 1831 has sum 15824 .. GENERATED FROM PYTHON SOURCE LINES 75-80 .. code-block:: Python :lineno-start: 75 plt.plot(r0tel.waveform[0, brightest_pixel], label="channel 0 (high-gain)") plt.plot(r0tel.waveform[1, brightest_pixel], label="channel 1 (low-gain)") plt.legend() .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_002.png :alt: ctapipe handson :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_002.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 81-86 .. code-block:: Python :lineno-start: 81 @interact def view_waveform(chan=0, pix_id=brightest_pixel): plt.plot(r0tel.waveform[chan, pix_id]) .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_003.png :alt: ctapipe handson :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_003.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none interactive(children=(IntSlider(value=0, description='chan', max=1), IntSlider(value=1831, description='pix_id', max=5493, min=-1831), Output()), _dom_classes=('widget-interact',)) .. GENERATED FROM PYTHON SOURCE LINES 87-89 try making this compare 2 waveforms .. GENERATED FROM PYTHON SOURCE LINES 92-101 Part 2: Explore the instrument description ------------------------------------------ This is all well and good, but we don’t really know what camera or telescope this is… how do we get instrumental description info? Currently this is returned *inside* the event (it will soon change to be separate in next version or so) .. GENERATED FROM PYTHON SOURCE LINES 101-104 .. code-block:: Python :lineno-start: 102 subarray = source.subarray .. GENERATED FROM PYTHON SOURCE LINES 105-107 .. code-block:: Python :lineno-start: 105 subarray .. rst-class:: sphx-glr-script-out .. code-block:: none SubarrayDescription(name='MonteCarloArray', n_tels=180) .. GENERATED FROM PYTHON SOURCE LINES 108-110 .. code-block:: Python :lineno-start: 108 subarray.peek() .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_004.png :alt: MonteCarloArray :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_004.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 111-113 .. code-block:: Python :lineno-start: 111 subarray.to_table() .. raw:: html
Table length=180
tel_idnametypepos_xpos_ypos_zcamera_nameoptics_namecamera_indexoptics_indextel_description
mmm
uint16str5str3float32float32float32str9str5int64int64str17
1LSTLST-20.643-64.81734.3LSTCamLST21LST_LST_LSTCam
2LSTLST79.993996-0.76829.4LSTCamLST21LST_LST_LSTCam
3LSTLST-19.39665.231.0LSTCamLST21LST_LST_LSTCam
4LSTLST-120.0331.15133.1LSTCamLST21LST_LST_LSTCam
5MSTMST-0.017-0.00124.35FlashCamMST12MST_MST_FlashCam
6MSTMST-1.468-151.22131.0FlashCamMST12MST_MST_FlashCam
7MSTMST-3.1379998-325.24539.0FlashCamMST12MST_MST_FlashCam
8MSTMST1.4339999151.2225.0FlashCamMST12MST_MST_FlashCam
9MSTMST3.1039999325.24323.5FlashCamMST12MST_MST_FlashCam
.................................
172ASTRISST260.0-920.065.0CHECASTRI00SST_ASTRI_CHEC
173ASTRISST-500.0815.015.0CHECASTRI00SST_ASTRI_CHEC
174ASTRISST-500.0-815.075.0CHECASTRI00SST_ASTRI_CHEC
175ASTRISST500.0815.045.0CHECASTRI00SST_ASTRI_CHEC
176ASTRISST500.0-815.053.0CHECASTRI00SST_ASTRI_CHEC
177ASTRISST-810.0655.012.0CHECASTRI00SST_ASTRI_CHEC
178ASTRISST-810.0-655.068.0CHECASTRI00SST_ASTRI_CHEC
179ASTRISST810.0655.020.0CHECASTRI00SST_ASTRI_CHEC
180ASTRISST810.0-655.041.0CHECASTRI00SST_ASTRI_CHEC


.. GENERATED FROM PYTHON SOURCE LINES 114-116 .. code-block:: Python :lineno-start: 114 subarray.tel[2] .. rst-class:: sphx-glr-script-out .. code-block:: none TelescopeDescription(type='LST', optics_name='LST', camera_name='LSTCam') .. GENERATED FROM PYTHON SOURCE LINES 117-119 .. code-block:: Python :lineno-start: 117 subarray.tel[2].camera .. rst-class:: sphx-glr-script-out .. code-block:: none CameraDescription(name=LSTCam, geometry=LSTCam, readout=LSTCam) .. GENERATED FROM PYTHON SOURCE LINES 120-122 .. code-block:: Python :lineno-start: 120 subarray.tel[2].optics .. rst-class:: sphx-glr-script-out .. code-block:: none OpticsDescription(name=LST, size_type=LST, reflector_shape=PARABOLIC, equivalent_focal_length=28.00 m, effective_focal_length=29.31 m, n_mirrors=1, mirror_area=386.73 m2) .. GENERATED FROM PYTHON SOURCE LINES 123-125 .. code-block:: Python :lineno-start: 123 tel = subarray.tel[2] .. GENERATED FROM PYTHON SOURCE LINES 126-128 .. code-block:: Python :lineno-start: 126 tel.camera .. rst-class:: sphx-glr-script-out .. code-block:: none CameraDescription(name=LSTCam, geometry=LSTCam, readout=LSTCam) .. GENERATED FROM PYTHON SOURCE LINES 129-131 .. code-block:: Python :lineno-start: 129 tel.optics .. rst-class:: sphx-glr-script-out .. code-block:: none OpticsDescription(name=LST, size_type=LST, reflector_shape=PARABOLIC, equivalent_focal_length=28.00 m, effective_focal_length=29.31 m, n_mirrors=1, mirror_area=386.73 m2) .. GENERATED FROM PYTHON SOURCE LINES 132-134 .. code-block:: Python :lineno-start: 132 tel.camera.geometry.pix_x .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 135-137 .. code-block:: Python :lineno-start: 135 tel.camera.geometry.to_table() .. raw:: html
Table length=1855
pix_idpix_xpix_ypix_area
mmm2
int64float64float64float64
00.00.00.002079326892271638
1-0.037796701200384926-0.0327324290367683250.002079326892271638
2-0.047245474963293450.016366662083976240.002079326892271638
3-0.0094487737629085230.0490990911207445650.002079326892271638
40.0377967012003849260.0327324290367683250.002079326892271638
50.04724547496329345-0.016366662083976240.002079326892271638
60.009448773762908523-0.0490990911207445650.002079326892271638
7-0.07559330728842001-0.06546483976983580.002079326892271638
8-0.08504208105132853-0.0163657486490912420.002079326892271638
............
1846-0.67088033276595160.96561882545319930.002079326892271638
18470.453564836008831741.08017546211512050.002079326892271638
18480.500810244496671.0638087511987850.002079326892271638
18491.16224173440920890.147289060475023960.002079326892271638
18501.17169050746812880.098189973012446710.002079326892271638
18510.708676884320607-0.93288632847675060.002079326892271638
18520.6708803327659516-0.96561882545319930.002079326892271638
1853-0.45356483600883174-1.08017546211512050.002079326892271638
1854-0.50081024449667-1.0638087511987850.002079326892271638


.. GENERATED FROM PYTHON SOURCE LINES 138-141 .. code-block:: Python :lineno-start: 138 tel.optics.mirror_area .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 142-144 .. code-block:: Python :lineno-start: 142 disp = CameraDisplay(tel.camera.geometry) .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_005.png :alt: LSTCam :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_005.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 145-151 .. code-block:: Python :lineno-start: 145 disp = CameraDisplay(tel.camera.geometry) disp.image = r0tel.waveform[ 0, :, 10 ] # display channel 0, sample 0 (try others like 10) .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_006.png :alt: LSTCam :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_006.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 152-155 \*\* aside: \*\* show demo using a CameraDisplay in interactive mode in ipython rather than notebook .. GENERATED FROM PYTHON SOURCE LINES 158-161 Part 3: Apply some calibration and trace integration ---------------------------------------------------- .. GENERATED FROM PYTHON SOURCE LINES 161-165 .. code-block:: Python :lineno-start: 163 calib = CameraCalibrator(subarray=subarray) .. GENERATED FROM PYTHON SOURCE LINES 166-170 .. code-block:: Python :lineno-start: 166 for event in EventSource(path, max_events=5): calib(event) # fills in r1, dl0, and dl1 print(event.dl1.tel.keys()) .. rst-class:: sphx-glr-script-out .. code-block:: none dict_keys([9, 14, 104]) dict_keys([26, 27, 29, 47, 53, 59, 69, 73, 121, 122, 124, 148, 162, 166]) dict_keys([82, 92, 175, 179]) dict_keys([17, 26, 27, 29, 43, 47, 53, 57, 69, 112, 121, 122, 124, 125, 128, 140, 144, 148, 162, 166]) dict_keys([3, 4, 8, 9, 14, 15, 23, 24, 25, 29, 32, 36, 42, 46, 52, 56, 103, 104, 109, 110, 118, 119, 120, 126, 127, 129, 130, 137, 139, 143, 161]) .. GENERATED FROM PYTHON SOURCE LINES 171-173 .. code-block:: Python :lineno-start: 171 event.dl1.tel[3] .. rst-class:: sphx-glr-script-out .. code-block:: none ctapipe.containers.DL1CameraContainer: image: Numpy array of camera image, after waveform extraction.Shape: (n_pixel) if n_channels is 1 or data is gain selectedelse: (n_channels, n_pixel) with default None peak_time: Numpy array containing position of the peak of the pulse as determined by the extractor.Shape: (n_pixel) if n_channels is 1 or data is gain selectedelse: (n_channels, n_pixel) with default None image_mask: Boolean numpy array where True means the pixel has passed cleaning. Shape: (n_pixel, ) with default None as a 1-D array with dtype bool is_valid: True if image extraction succeeded, False if failed or in the case of TwoPass methods, that the first pass only was returned. with default False parameters: Image parameters with default None with type .. GENERATED FROM PYTHON SOURCE LINES 174-176 .. code-block:: Python :lineno-start: 174 dl1tel = event.dl1.tel[3] .. GENERATED FROM PYTHON SOURCE LINES 177-179 .. code-block:: Python :lineno-start: 177 dl1tel.image.shape # note this will be gain-selected in next version, so will be just 1D array of 1855 .. rst-class:: sphx-glr-script-out .. code-block:: none (1855,) .. GENERATED FROM PYTHON SOURCE LINES 180-182 .. code-block:: Python :lineno-start: 180 dl1tel.peak_time .. rst-class:: sphx-glr-script-out .. code-block:: none array([27.381046, 17.746258, 18.80205 , ..., 11.64868 , 31.961287, 29.717525], shape=(1855,), dtype=float32) .. GENERATED FROM PYTHON SOURCE LINES 183-185 .. code-block:: Python :lineno-start: 183 CameraDisplay(tel.camera.geometry, image=dl1tel.image) .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_007.png :alt: LSTCam :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_007.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 186-189 .. code-block:: Python :lineno-start: 186 CameraDisplay(tel.camera.geometry, image=dl1tel.peak_time) .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_008.png :alt: LSTCam :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_008.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 190-192 Now for Hillas Parameters .. GENERATED FROM PYTHON SOURCE LINES 192-198 .. code-block:: Python :lineno-start: 194 image = dl1tel.image mask = tailcuts_clean(tel.camera.geometry, image, picture_thresh=10, boundary_thresh=5) mask .. rst-class:: sphx-glr-script-out .. code-block:: none array([False, False, False, ..., False, False, False], shape=(1855,)) .. GENERATED FROM PYTHON SOURCE LINES 199-201 .. code-block:: Python :lineno-start: 199 CameraDisplay(tel.camera.geometry, image=mask) .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_009.png :alt: LSTCam :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_009.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 202-205 .. code-block:: Python :lineno-start: 202 cleaned = image.copy() cleaned[~mask] = 0 .. GENERATED FROM PYTHON SOURCE LINES 206-212 .. code-block:: Python :lineno-start: 206 disp = CameraDisplay(tel.camera.geometry, image=cleaned) disp.cmap = plt.cm.coolwarm disp.add_colorbar() plt.xlim(0.5, 1.0) plt.ylim(-1.0, 0.0) .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_010.png :alt: LSTCam :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_010.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none (-1.0, 0.0) .. GENERATED FROM PYTHON SOURCE LINES 213-216 .. code-block:: Python :lineno-start: 213 params = hillas_parameters(tel.camera.geometry, cleaned) print(params) .. rst-class:: sphx-glr-script-out .. code-block:: none {'intensity': np.float64(32.78351974487305), 'kurtosis': np.float64(1.2333006676347782), 'length': , 'length_uncertainty': , 'phi': , 'psi': , 'psi_uncertainty': , 'r': , 'skewness': np.float64(-0.271615381552047), 'transverse_cog_uncertainty': , 'width': , 'width_uncertainty': , 'x': , 'y': } .. GENERATED FROM PYTHON SOURCE LINES 217-225 .. code-block:: Python :lineno-start: 217 disp = CameraDisplay(tel.camera.geometry, image=cleaned) disp.cmap = plt.cm.coolwarm disp.add_colorbar() plt.xlim(0.5, 1.0) plt.ylim(-1.0, 0.0) disp.overlay_moments(params, color="white", lw=2) .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_011.png :alt: LSTCam :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_011.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none Ignoring fixed x limits to fulfill fixed data aspect with adjustable data limits. .. GENERATED FROM PYTHON SOURCE LINES 226-236 Part 4: Let’s put it all together: ---------------------------------- - loop over events, selecting only telescopes of the same type (e.g. LST:LSTCam) - for each event, apply calibration/trace integration - calculate Hillas parameters - write out all hillas parameters to a file that can be loaded with Pandas .. GENERATED FROM PYTHON SOURCE LINES 239-241 first let’s select only those telescopes with LST:LSTCam .. GENERATED FROM PYTHON SOURCE LINES 241-244 .. code-block:: Python :lineno-start: 242 subarray.telescope_types .. rst-class:: sphx-glr-script-out .. code-block:: none (TelescopeDescription(type='SST', optics_name='ASTRI', camera_name='CHEC'), TelescopeDescription(type='LST', optics_name='LST', camera_name='LSTCam'), TelescopeDescription(type='MST', optics_name='MST', camera_name='NectarCam'), TelescopeDescription(type='MST', optics_name='MST', camera_name='FlashCam')) .. GENERATED FROM PYTHON SOURCE LINES 245-248 .. code-block:: Python :lineno-start: 245 subarray.get_tel_ids_for_type("LST_LST_LSTCam") .. rst-class:: sphx-glr-script-out .. code-block:: none (1, 2, 3, 4) .. GENERATED FROM PYTHON SOURCE LINES 249-251 Now let’s write out program .. GENERATED FROM PYTHON SOURCE LINES 251-255 .. code-block:: Python :lineno-start: 252 data = utils.get_dataset_path("gamma_prod5.simtel.zst") source = EventSource(data) # remove the max_events limit to get more stats .. GENERATED FROM PYTHON SOURCE LINES 256-266 .. code-block:: Python :lineno-start: 256 for event in source: calib(event) for tel_id, tel_data in event.dl1.tel.items(): tel = source.subarray.tel[tel_id] mask = tailcuts_clean(tel.camera.geometry, tel_data.image) if np.count_nonzero(mask) > 0: params = hillas_parameters(tel.camera.geometry[mask], tel_data.image[mask]) .. GENERATED FROM PYTHON SOURCE LINES 267-279 .. code-block:: Python :lineno-start: 267 with HDF5TableWriter(filename="hillas.h5", group_name="dl1", overwrite=True) as writer: source = EventSource(data, allowed_tels=[1, 2, 3, 4], max_events=10) for event in source: calib(event) for tel_id, tel_data in event.dl1.tel.items(): tel = source.subarray.tel[tel_id] mask = tailcuts_clean(tel.camera.geometry, tel_data.image) params = hillas_parameters(tel.camera.geometry[mask], tel_data.image[mask]) writer.write("hillas", params) .. GENERATED FROM PYTHON SOURCE LINES 280-283 We can now load in the file we created and plot it ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ .. GENERATED FROM PYTHON SOURCE LINES 283-286 .. code-block:: Python :lineno-start: 283 glob.glob("*.h5") .. rst-class:: sphx-glr-script-out .. code-block:: none ['hillas.h5'] .. GENERATED FROM PYTHON SOURCE LINES 287-290 .. code-block:: Python :lineno-start: 287 hillas = pd.read_hdf("hillas.h5", key="/dl1/hillas") hillas .. raw:: html
intensity skewness kurtosis x y r phi length length_uncertainty width width_uncertainty psi psi_uncertainty transverse_cog_uncertainty
0 47.098688 -0.530399 1.881643 0.695206 -0.847657 1.096282 -50.643096 0.032743 0.002240 0.026627 0.002372 -51.309393 138.997232 0.006644
1 126.314263 0.020218 1.742869 0.464446 -0.999899 1.102500 -65.085472 0.055256 0.002119 0.025090 0.001144 -67.386572 7.356117 0.004057


.. GENERATED FROM PYTHON SOURCE LINES 291-294 .. code-block:: Python :lineno-start: 291 _ = hillas.hist(figsize=(8, 8)) .. image-sg:: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_012.png :alt: intensity, skewness, kurtosis, x, y, r, phi, length, length_uncertainty, width, width_uncertainty, psi, psi_uncertainty, transverse_cog_uncertainty :srcset: /auto_examples/tutorials/images/sphx_glr_ctapipe_handson_012.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 295-298 If you do this yourself, chose a larger file to loop over more events to get better statistics .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 8.326 seconds) .. _sphx_glr_download_auto_examples_tutorials_ctapipe_handson.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: ctapipe_handson.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: ctapipe_handson.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: ctapipe_handson.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_