__init__.py 59 KB
Newer Older
1
2
3
4
5
6
# Event source for MAGIC calibrated data files.
# Requires uproot package (https://github.com/scikit-hep/uproot).

import glob
import re

7
import scipy
8
9
10
11
12
13
import numpy as np
import scipy.interpolate

from astropy import units as u
from astropy.time import Time
from ctapipe.io.eventsource import EventSource
14
from ctapipe.io.containers import DataContainer, TelescopePointingContainer
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
from ctapipe.instrument import TelescopeDescription, SubarrayDescription, OpticsDescription, CameraGeometry

__all__ = ['MAGICEventSource']


class MAGICEventSource(EventSource):
    """
    EventSource for MAGIC calibrated data.

    This class operates with the MAGIC data run-wise. This means that the files
    corresponding to the same data run are loaded and processed together.
    """
    _count = 0

    def __init__(self, config=None, tool=None, **kwargs):
        """
        Constructor

        Parameters
        ----------
        config: traitlets.loader.Config
            Configuration specified by config file or cmdline arguments.
            Used to set traitlet values.
            Set to None if no configuration to pass.
        tool: ctapipe.core.Tool
            Tool executable that is calling this component.
            Passes the correct logger to the component.
            Set to None if no Tool to pass.
        kwargs: dict
            Additional parameters to be passed.
            NOTE: The file mask of the data to read can be passed with
            the 'input_url' parameter.
        """

        try:
            import uproot
        except ImportError:
            msg = "The `uproot` python module is required to access the MAGIC data"
            self.log.error(msg)
            raise

56
57
        self.file_list = glob.glob(kwargs['input_url'])
        self.file_list.sort()
58
59
60
61
62

        # EventSource can not handle file wild cards as input_url
        # To overcome this we substitute the input_url with first file matching
        # the specified file mask.
        del kwargs['input_url']
63
        super().__init__(input_url=self.file_list[0], **kwargs)
64
65

        # Retrieving the list of run numbers corresponding to the data files
66
        run_numbers = list(map(self._get_run_number, self.file_list))
67
68
69
70
71
72
73
74
75
76
77
78
79
80
        self.run_numbers = np.unique(run_numbers)

        # # Setting up the current run with the first run present in the data
        # self.current_run = self._set_active_run(run_number=0)
        self.current_run = None

        # MAGIC telescope positions in m wrt. to the center of CTA simulations
        self.magic_tel_positions = {
            1: [-27.24, -146.66, 50.00] * u.m,
            2: [-96.44, -96.77, 51.00] * u.m
        }
        # MAGIC telescope description
        optics = OpticsDescription.from_name('MAGIC')
        geom = CameraGeometry.from_name('MAGICCam')
81
        self.magic_tel_description = TelescopeDescription(name='MAGIC', tel_type='MAGIC', optics=optics, camera=geom)
82
        self.magic_tel_descriptions = {1: self.magic_tel_description, 2: self.magic_tel_description}
83
        self._subarray_info = SubarrayDescription('MAGIC', self.magic_tel_positions, self.magic_tel_descriptions)
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

    @staticmethod
    def is_compatible(file_mask):
        """
        This method checks if the specified file mask corresponds
        to MAGIC data files. The result will be True only if all
        the files are of ROOT format and contain an 'Events' tree.

        Parameters
        ----------
        file_mask: str
            A file mask to check

        Returns
        -------
        bool:
            True if the masked files are MAGIC data runs, False otherwise.

        """

        is_magic_root_file = True

        file_list = glob.glob(file_mask)

        for file_path in file_list:
            try:
                import uproot

                try:
                    with uproot.open(file_path) as input_data:
                        if 'Events' not in input_data:
                            is_magic_root_file = False
                except ValueError:
                    # uproot raises ValueError if the file is not a ROOT file
                    is_magic_root_file = False
                    pass

            except ImportError:
Ievgen Vovk's avatar
Ievgen Vovk committed
122
                if re.match(r'.+_m\d_.+root', file_path.lower()) is None:
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
                    is_magic_root_file = False

        return is_magic_root_file

    @staticmethod
    def _get_run_number(file_name):
        """
        This internal method extracts the run number from
        the specified file name.

        Parameters
        ----------
        file_name: str
            A file name to process.

        Returns
        -------
        int:
            A run number of the file.
        """

144
        mask = r".*\d+_M\d+_(\d+)\.\d+_.*"
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
        parsed_info = re.findall(mask, file_name)

        try:
            run_number = int(parsed_info[0])
        except IndexError:
            raise IndexError('Can not identify the run number of the file {:s}'.format(file_name))

        return run_number

    def _set_active_run(self, run_number):
        """
        This internal method sets the run that will be used for data loading.

        Parameters
        ----------
        run_number: int
            The run number to use.

        Returns
        -------

        """

        input_path = '/'.join(self.input_url.split('/')[:-1])
        this_run_mask = input_path + '/*{:d}*root'.format(run_number)

        run = dict()
        run['number'] = run_number
        run['read_events'] = 0
174
        run['data'] = MarsDataRun(run_file_mask=this_run_mask, filter_list=self.file_list)
175
176
177

        return run

178
179
180
181
    @property
    def subarray(self):
        return self._subarray_info

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
256
257
258
259
260
261
262
263
264
265
266
267
268
    def _generator(self):
        """
        The default event generator. Return the stereo event
        generator instance.

        Returns
        -------

        """

        return self._stereo_event_generator()

    def _stereo_event_generator(self):
        """
        Stereo event generator. Yields DataContainer instances, filled
        with the read event data.

        Returns
        -------

        """

        counter = 0

        # Data container - is initialized once, and data is replaced within it after each yield
        data = DataContainer()
        data.meta['origin'] = "MAGIC"
        data.meta['input_url'] = self.input_url
        data.meta['is_simulation'] = False

        # Telescopes with data:
        tels_in_file = ["m1", "m2"]
        tels_with_data = {1, 2}

        # Loop over the available data runs
        for run_number in self.run_numbers:

            # Removing the previously read data run from memory
            if self.current_run is not None:
                if 'data' in self.current_run:
                    del self.current_run['data']

            # Setting the new active run
            self.current_run = self._set_active_run(run_number)

            # Loop over the events
            for event_i in range(self.current_run['data'].n_stereo_events):
                # Event and run ids
                event_order_number = self.current_run['data'].stereo_ids[event_i][0]
                event_id = self.current_run['data'].event_data['M1']['stereo_event_number'][event_order_number]
                obs_id = self.current_run['number']

                # Reading event data
                event_data = self.current_run['data'].get_stereo_event_data(event_i)

                # Event counter
                data.count = counter

                # Setting up the R0 container
                data.r0.obs_id = obs_id
                data.r0.event_id = event_id
                data.r0.tel.clear()

                # Setting up the R1 container
                data.r1.obs_id = obs_id
                data.r1.event_id = event_id
                data.r1.tel.clear()

                # Setting up the DL0 container
                data.dl0.obs_id = obs_id
                data.dl0.event_id = event_id
                data.dl0.tel.clear()

                # Filling the DL1 container with the event data
                for tel_i, tel_id in enumerate(tels_in_file):
                    # Creating the telescope pointing container
                    pointing = TelescopePointingContainer()
                    pointing.azimuth = np.deg2rad(event_data['pointing_az']) * u.rad
                    pointing.altitude = np.deg2rad(90 - event_data['pointing_zd']) * u.rad
                    pointing.ra = np.deg2rad(event_data['pointing_ra']) * u.rad
                    pointing.dec = np.deg2rad(event_data['pointing_dec']) * u.rad

                    # Adding the pointing container to the event data
                    data.pointing[tel_i + 1] = pointing

                    # Adding event charge and peak positions per pixel
                    data.dl1.tel[tel_i + 1].image = event_data['{:s}_image'.format(tel_id)]
Ievgen Vovk's avatar
Ievgen Vovk committed
269
                    data.dl1.tel[tel_i + 1].pulse_time = event_data['{:s}_pulse_time'.format(tel_id)]
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
                    # data.dl1.tel[i_tel + 1].badpixels = np.array(
                    #     file['dl1/tel' + str(i_tel + 1) + '/badpixels'], dtype=np.bool)

                # Adding the event arrival time
                time_tmp = Time(event_data['mjd'], scale='utc', format='mjd')
                data.trig.gps_time = Time(time_tmp, format='unix', scale='utc', precision=9)

                # Setting the telescopes with data
                data.r0.tels_with_data = tels_with_data
                data.r1.tels_with_data = tels_with_data
                data.dl0.tels_with_data = tels_with_data
                data.trig.tels_with_trigger = tels_with_data

                yield data
                counter += 1

        return

    def _mono_event_generator(self, telescope):
        """
        Mono event generator. Yields DataContainer instances, filled
        with the read event data.

        Parameters
        ----------
        telescope: str
            The telescope for which to return events. Can be either "M1" or "M2".

        Returns
        -------

        """

        counter = 0
        telescope = telescope.upper()

        # Data container - is initialized once, and data is replaced within it after each yield
        data = DataContainer()
        data.meta['origin'] = "MAGIC"
        data.meta['input_url'] = self.input_url
        data.meta['is_simulation'] = False

        # Telescopes with data:
        tels_in_file = ["M1", "M2"]

        if telescope not in tels_in_file:
            raise ValueError("Specified telescope {:s} is not in the allowed list {}".format(telescope, tels_in_file))

        tel_i = tels_in_file.index(telescope)
        tels_with_data = {tel_i + 1, }

        # Loop over the available data runs
        for run_number in self.run_numbers:

            # Removing the previously read data run from memory
            if self.current_run is not None:
                if 'data' in self.current_run:
                    del self.current_run['data']

            # Setting the new active run
            self.current_run = self._set_active_run(run_number)

            if telescope == 'M1':
                n_events = self.current_run['data'].n_mono_events_m1
            else:
                n_events = self.current_run['data'].n_mono_events_m2

            # Loop over the events
            for event_i in range(n_events):
                # Event and run ids
                event_order_number = self.current_run['data'].mono_ids[telescope][event_i]
                event_id = self.current_run['data'].event_data[telescope]['stereo_event_number'][event_order_number]
                obs_id = self.current_run['number']

                # Reading event data
                event_data = self.current_run['data'].get_mono_event_data(event_i, telescope=telescope)

                # Event counter
                data.count = counter

                # Setting up the R0 container
                data.r0.obs_id = obs_id
                data.r0.event_id = event_id
                data.r0.tel.clear()

                # Setting up the R1 container
                data.r1.obs_id = obs_id
                data.r1.event_id = event_id
                data.r1.tel.clear()

                # Setting up the DL0 container
                data.dl0.obs_id = obs_id
                data.dl0.event_id = event_id
                data.dl0.tel.clear()

                # Creating the telescope pointing container
                pointing = TelescopePointingContainer()
                pointing.azimuth = np.deg2rad(event_data['pointing_az']) * u.rad
                pointing.altitude = np.deg2rad(90 - event_data['pointing_zd']) * u.rad
                pointing.ra = np.deg2rad(event_data['pointing_ra']) * u.rad
                pointing.dec = np.deg2rad(event_data['pointing_dec']) * u.rad

                # Adding the pointing container to the event data
                data.pointing[tel_i + 1] = pointing

                # Adding event charge and peak positions per pixel
                data.dl1.tel[tel_i + 1].image = event_data['image']
Ievgen Vovk's avatar
Ievgen Vovk committed
377
                data.dl1.tel[tel_i + 1].pulse_time = event_data['pulse_time']
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
                # data.dl1.tel[tel_i + 1].badpixels = np.array(
                #     file['dl1/tel' + str(i_tel + 1) + '/badpixels'], dtype=np.bool)

                # Adding the event arrival time
                time_tmp = Time(event_data['mjd'], scale='utc', format='mjd')
                data.trig.gps_time = Time(time_tmp, format='unix', scale='utc', precision=9)

                # Setting the telescopes with data
                data.r0.tels_with_data = tels_with_data
                data.r1.tels_with_data = tels_with_data
                data.dl0.tels_with_data = tels_with_data
                data.trig.tels_with_trigger = tels_with_data

                yield data
                counter += 1

        return

    def _pedestal_event_generator(self, telescope):
        """
        Pedestal event generator. Yields DataContainer instances, filled
        with the read event data.

        Parameters
        ----------
        telescope: str
            The telescope for which to return events. Can be either "M1" or "M2".

        Returns
        -------

        """

        counter = 0
        telescope = telescope.upper()

        # Data container - is initialized once, and data is replaced within it after each yield
        data = DataContainer()
        data.meta['origin'] = "MAGIC"
        data.meta['input_url'] = self.input_url
        data.meta['is_simulation'] = False

        # Telescopes with data:
        tels_in_file = ["M1", "M2"]

        if telescope not in tels_in_file:
            raise ValueError("Specified telescope {:s} is not in the allowed list {}".format(telescope, tels_in_file))

        tel_i = tels_in_file.index(telescope)
        tels_with_data = {tel_i + 1, }

        # Loop over the available data runs
        for run_number in self.run_numbers:

            # Removing the previously read data run from memory
            if self.current_run is not None:
                if 'data' in self.current_run:
                    del self.current_run['data']

            # Setting the new active run
            self.current_run = self._set_active_run(run_number)

            if telescope == 'M1':
                n_events = self.current_run['data'].n_pedestal_events_m1
            else:
                n_events = self.current_run['data'].n_pedestal_events_m2

            # Loop over the events
            for event_i in range(n_events):
                # Event and run ids
                event_order_number = self.current_run['data'].pedestal_ids[telescope][event_i]
                event_id = self.current_run['data'].event_data[telescope]['stereo_event_number'][event_order_number]
                obs_id = self.current_run['number']

                # Reading event data
                event_data = self.current_run['data'].get_pedestal_event_data(event_i, telescope=telescope)

                # Event counter
                data.count = counter

                # Setting up the R0 container
                data.r0.obs_id = obs_id
                data.r0.event_id = event_id
                data.r0.tel.clear()

                # Setting up the R1 container
                data.r1.obs_id = obs_id
                data.r1.event_id = event_id
                data.r1.tel.clear()

                # Setting up the DL0 container
                data.dl0.obs_id = obs_id
                data.dl0.event_id = event_id
                data.dl0.tel.clear()

                # Creating the telescope pointing container
                pointing = TelescopePointingContainer()
                pointing.azimuth = np.deg2rad(event_data['pointing_az']) * u.rad
                pointing.altitude = np.deg2rad(90 - event_data['pointing_zd']) * u.rad
                pointing.ra = np.deg2rad(event_data['pointing_ra']) * u.rad
                pointing.dec = np.deg2rad(event_data['pointing_dec']) * u.rad

                # Adding the pointing container to the event data
                data.pointing[tel_i + 1] = pointing

                # Adding event charge and peak positions per pixel
                data.dl1.tel[tel_i + 1].image = event_data['image']
Ievgen Vovk's avatar
Ievgen Vovk committed
485
                data.dl1.tel[tel_i + 1].pulse_time = event_data['pulse_time']
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
                # data.dl1.tel[tel_i + 1].badpixels = np.array(
                #     file['dl1/tel' + str(i_tel + 1) + '/badpixels'], dtype=np.bool)

                # Adding the event arrival time
                time_tmp = Time(event_data['mjd'], scale='utc', format='mjd')
                data.trig.gps_time = Time(time_tmp, format='unix', scale='utc', precision=9)

                # Setting the telescopes with data
                data.r0.tels_with_data = tels_with_data
                data.r1.tels_with_data = tels_with_data
                data.dl0.tels_with_data = tels_with_data
                data.trig.tels_with_trigger = tels_with_data

                yield data
                counter += 1

        return


505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
class MAGICEventSourceMC(EventSource):
    """
    EventSource for MAGIC calibrated MCs.
    """
    _count = 0

    def __init__(self, config=None, tool=None, **kwargs):
        """
        Constructor

        Parameters
        ----------
        config: traitlets.loader.Config
            Configuration specified by config file or cmdline arguments.
            Used to set traitlet values.
            Set to None if no configuration to pass.
        tool: ctapipe.core.Tool
            Tool executable that is calling this component.
            Passes the correct logger to the component.
            Set to None if no Tool to pass.
        kwargs: dict
            Additional parameters to be passed.
            NOTE: The file mask of the data to read can be passed with
            the 'input_url' parameter.
        """

        try:
            import uproot
        except ImportError:
            msg = "The `uproot` python module is required to access the MAGIC data"
            self.log.error(msg)
            raise
537
538
539
        
        if len(glob.glob(kwargs['input_url'])) > 1:
            raise ImportError('MC data can no be loaded with wildcards. Please load them run by run')  
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554

        self.file_name = kwargs['input_url']

        super().__init__(**kwargs)

        self.mc_file = MarsMCFile(self.file_name)

        # MAGIC telescope positions in m wrt. to the center of CTA simulations
        self.magic_tel_positions = {
            1: [-27.24, -146.66, 50.00] * u.m,
            2: [-96.44, -96.77, 51.00] * u.m
        }
        # MAGIC telescope description
        optics = OpticsDescription.from_name('MAGIC')
        geom = CameraGeometry.from_name('MAGICCam')
555
        self.magic_tel_description = TelescopeDescription(name='MAGIC', tel_type='MAGIC', optics=optics, camera=geom)
556
        self.magic_tel_descriptions = {1: self.magic_tel_description, 2: self.magic_tel_description}
557
        self._subarray_info = SubarrayDescription('MAGIC', self.magic_tel_positions, self.magic_tel_descriptions)
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597

    @staticmethod
    def is_compatible(file_name):
        """
        This method checks if the specified file name corresponds
        to a MAGIC data file. The result will be True only if all
        the files are of ROOT format and contain an 'Events' tree.

        Parameters
        ----------
        file_name: str
            A file name to check

        Returns
        -------
        bool:
            True if the given file is MAGIC data file, False otherwise.

        """

        is_magic_root_file = True

        try:
            import uproot

            try:
                with uproot.open(file_name) as input_data:
                    if 'Events' not in input_data:
                        is_magic_root_file = False
            except ValueError:
                # uproot raises ValueError if the file is not a ROOT file
                is_magic_root_file = False
                pass

        except ImportError:
            if re.match(r'.+_m\d_.+root', file_name.lower()) is None:
                is_magic_root_file = False

        return is_magic_root_file

598
599
600
601
    @property
    def subarray(self):
        return self._subarray_info

602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
    def _generator(self):
        """
        The default event generator. Return the stereo event
        generator instance.

        Returns
        -------

        """

        return self._mono_event_generator()

    def _mono_event_generator(self):
        """
        Mono event generator. Yields DataContainer instances, filled
        with the read event data.

        Returns
        -------
        ctapipe.io.containers.DataContainer

        """

        counter = 0

        # Data container - is initialized once, and data is replaced within it after each yield
        data = DataContainer()
        data.meta['origin'] = "MAGIC MC"
        data.meta['input_url'] = self.input_url
631
        data.meta['is_simulation'] = True
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699

        tels_with_data = {self.mc_file.telescope, }

        # Loop over the events
        for event_i in range(self.mc_file.n_mono_events):
            # Event and run ids
            event_order_number = self.mc_file.mono_ids[event_i]
            event_id = self.mc_file.event_data['daq_event_number'][event_order_number]
            obs_id = self.mc_file.run_number

            # Reading event data
            event_data = self.mc_file.get_mono_event_data(event_i)

            # Event counter
            data.count = counter

            # Setting up the R0 container
            data.r0.obs_id = obs_id
            data.r0.event_id = event_id
            data.r0.tel.clear()

            # Setting up the R1 container
            data.r1.obs_id = obs_id
            data.r1.event_id = event_id
            data.r1.tel.clear()

            # Setting up the DL0 container
            data.dl0.obs_id = obs_id
            data.dl0.event_id = event_id
            data.dl0.tel.clear()

            # Setting up the MC container
            data.mc.tel.clear()

            # Creating the telescope pointing container
            pointing = TelescopePointingContainer()
            pointing.azimuth = np.deg2rad(event_data['pointing_az']) * u.rad
            pointing.altitude = np.deg2rad(90 - event_data['pointing_zd']) * u.rad

            # Adding the pointing container to the event data
            data.pointing[self.mc_file.telescope] = pointing

            # Adding event charge and peak positions per pixel
            data.dl1.tel[self.mc_file.telescope].image = event_data['image']
            data.dl1.tel[self.mc_file.telescope].pulse_time = event_data['pulse_time']

            # Setting the telescopes with data
            data.r0.tels_with_data = tels_with_data
            data.r1.tels_with_data = tels_with_data
            data.dl0.tels_with_data = tels_with_data
            data.trig.tels_with_trigger = tels_with_data

            # mc = data.mc.tel[self.mc_file.telescope]
            # mc.dc_to_pe = array_event['laser_calibrations'][tel_id]['calib']
            # mc.pedestal = array_event['camera_monitorings'][tel_id]['pedestal']
            # mc.reference_pulse_shape = pixel_settings['ref_shape'].astype('float64')
            # mc.meta['refstep'] = float(pixel_settings['ref_step'])
            # mc.time_slice = float(pixel_settings['time_slice'])
            # mc.photo_electron_image = (
            #     array_event
            #         .get('photoelectrons', {})
            #         .get(tel_index, {})
            #         .get('photoelectrons', np.zeros(n_pixel, dtype='float32'))
            # )

            data.mc.energy = event_data['true_energy'] * u.GeV
            data.mc.alt = (90 - event_data['true_zd']) * u.deg
            data.mc.az = event_data['true_az'] * u.deg
700
701
702
703
            data.mc.shower_primary_id = 1 - event_data['true_shower_primary_id']
            data.mc.h_first_int = event_data['true_h_first_int'] * u.cm
            data.mc.core_x = event_data['true_core_x'] * u.cm
            data.mc.core_y = event_data['true_core_y'] * u.cm
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726

            yield data
            counter += 1

        return

    def _pedestal_event_generator(self):
        """
        Pedestal event generator. Yields DataContainer instances, filled
        with the read event data.

        Returns
        -------
        ctapipe.io.containers.DataContainer

        """

        counter = 0

        # Data container - is initialized once, and data is replaced within it after each yield
        data = DataContainer()
        data.meta['origin'] = "MAGIC MC"
        data.meta['input_url'] = self.input_url
727
        data.meta['is_simulation'] = True
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795

        tels_with_data = {self.mc_file.telescope, }

        # Loop over the events
        for event_i in range(self.mc_file.n_pedestal_events):
            # Event and run ids
            event_order_number = self.mc_file.pedestal_ids[event_i]
            event_id = self.mc_file.event_data['daq_event_number'][event_order_number]
            obs_id = self.mc_file.run_number

            # Reading event data
            event_data = self.mc_file.get_mono_event_data(event_i)

            # Event counter
            data.count = counter

            # Setting up the R0 container
            data.r0.obs_id = obs_id
            data.r0.event_id = event_id
            data.r0.tel.clear()

            # Setting up the R1 container
            data.r1.obs_id = obs_id
            data.r1.event_id = event_id
            data.r1.tel.clear()

            # Setting up the DL0 container
            data.dl0.obs_id = obs_id
            data.dl0.event_id = event_id
            data.dl0.tel.clear()

            # Setting up the MC container
            data.mc.tel.clear()

            # Creating the telescope pointing container
            pointing = TelescopePointingContainer()
            pointing.azimuth = np.deg2rad(event_data['pointing_az']) * u.rad
            pointing.altitude = np.deg2rad(90 - event_data['pointing_zd']) * u.rad

            # Adding the pointing container to the event data
            data.pointing[self.mc_file.telescope] = pointing

            # Adding event charge and peak positions per pixel
            data.dl1.tel[self.mc_file.telescope].image = event_data['image']
            data.dl1.tel[self.mc_file.telescope].pulse_time = event_data['pulse_time']

            # Setting the telescopes with data
            data.r0.tels_with_data = tels_with_data
            data.r1.tels_with_data = tels_with_data
            data.dl0.tels_with_data = tels_with_data
            data.trig.tels_with_trigger = tels_with_data

            # mc = data.mc.tel[self.mc_file.telescope]
            # mc.dc_to_pe = array_event['laser_calibrations'][tel_id]['calib']
            # mc.pedestal = array_event['camera_monitorings'][tel_id]['pedestal']
            # mc.reference_pulse_shape = pixel_settings['ref_shape'].astype('float64')
            # mc.meta['refstep'] = float(pixel_settings['ref_step'])
            # mc.time_slice = float(pixel_settings['time_slice'])
            # mc.photo_electron_image = (
            #     array_event
            #         .get('photoelectrons', {})
            #         .get(tel_index, {})
            #         .get('photoelectrons', np.zeros(n_pixel, dtype='float32'))
            # )

            data.mc.energy = event_data['true_energy'] * u.GeV
            data.mc.alt = (90 - event_data['true_zd']) * u.deg
            data.mc.az = event_data['true_az'] * u.deg
796
797
798
799
            data.mc.shower_primary_id = 1 - event_data['true_shower_primary_id']
            data.mc.h_first_int = event_data['true_h_first_int'] * u.m
            data.mc.core_x = event_data['true_core_x'] * u.cm
            data.mc.core_y = event_data['true_core_y'] * u.cm
800
801
802
803
804
805
806

            yield data
            counter += 1

        return


807
808
809
810
811
class MarsDataRun:
    """
    This class implements reading of the event data from a single MAGIC data run.
    """

812
    def __init__(self, run_file_mask, filter_list=None):
813
814
815
816
817
818
819
820
821
        """
        Constructor of the class. Defines the run to use and the camera pixel arrangement.

        Parameters
        ----------
        run_file_mask: str
            A path mask for files belonging to the run. Must correspond to a single run
            or an exception will be raised. Must correspond to calibrated ("sorcerer"-level)
            data.
822
823
824
        filter_list: list, optional
            A list of files, to which the run_file_mask should be applied. If None, all the
            files satisfying run_file_mask will be used. Defaults to None.
825
826
827
828
829
830
        """

        self.run_file_mask = run_file_mask

        # Preparing the lists of M1/2 data files
        file_list = glob.glob(run_file_mask)
831
832
833
834
835

        # Filtering out extra files if necessary
        if filter_list is not None:
            file_list = list(set(file_list) & set(filter_list))

836
837
        self.m1_file_list = list(filter(lambda name: '_M1_' in name, file_list))
        self.m2_file_list = list(filter(lambda name: '_M2_' in name, file_list))
838
        self.m1_file_list.sort()
839
840
841
842
        self.m2_file_list.sort()

        # Retrieving the list of run numbers corresponding to the data files
        run_numbers = list(map(self._get_run_number, file_list))
843
        run_numbers = scipy.unique(run_numbers)
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860

        # Checking if a single run is going to be read
        if len(run_numbers) > 1:
            raise ValueError("Run mask corresponds to more than one run: {}".format(run_numbers))

        # Reading the event data
        self.event_data = dict()
        self.event_data['M1'] = self.load_events(self.m1_file_list)
        self.event_data['M2'] = self.load_events(self.m2_file_list)

        # Detecting pedestal events
        self.pedestal_ids = self._find_pedestal_events()
        # Detecting stereo events
        self.stereo_ids = self._find_stereo_events()
        # Detecting mono events
        self.mono_ids = self._find_mono_events()

861
862
        self.n_camera_pixels = 1039

863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
    @property
    def n_events_m1(self):
        return len(self.event_data['M1']['MJD'])

    @property
    def n_events_m2(self):
        return len(self.event_data['M2']['MJD'])

    @property
    def n_stereo_events(self):
        return len(self.stereo_ids)

    @property
    def n_mono_events_m1(self):
        return len(self.mono_ids['M1'])

    @property
    def n_mono_events_m2(self):
        return len(self.mono_ids['M2'])

    @property
    def n_pedestal_events_m1(self):
        return len(self.pedestal_ids['M1'])

    @property
    def n_pedestal_events_m2(self):
        return len(self.pedestal_ids['M2'])

    @staticmethod
    def _get_run_number(file_name):
        """
        This internal method extracts the run number from
        a specified file name.

        Parameters
        ----------
        file_name: str
            A file name to process.

        Returns
        -------
        int:
            A run number of the file.
        """

908
        mask = r".*\d+_M\d+_(\d+)\.\d+_.*"
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
        parsed_info = re.findall(mask, file_name)

        run_number = int(parsed_info[0])

        return run_number

    @staticmethod
    def load_events(file_list):
        """
        This method loads events from the pre-defiled file and returns them as a dictionary.

        Parameters
        ----------
        file_name: str
            Name of the MAGIC calibrated file to use.

        Returns
        -------
        dict:
            A dictionary with the even properties: charge / arrival time data, trigger, direction etc.
        """

        try:
            import uproot
        except ImportError:
            msg = "The `uproot` python module is required to access the MAGIC data"
            raise ImportError(msg)

        event_data = dict()

        event_data['charge'] = []
        event_data['arrival_time'] = []
941
942
943
944
945
946
947
        event_data['trigger_pattern'] = scipy.array([])
        event_data['stereo_event_number'] = scipy.array([])
        event_data['pointing_zd'] = scipy.array([])
        event_data['pointing_az'] = scipy.array([])
        event_data['pointing_ra'] = scipy.array([])
        event_data['pointing_dec'] = scipy.array([])
        event_data['MJD'] = scipy.array([])
948
949
950

        event_data['file_edges'] = [0]

951
952
953
954
955
956
957
958
959
960
        array_list = ['MCerPhotEvt.fPixels.fPhot', 'MArrivalTime.fData',
                      'MTriggerPattern.fPrescaled',
                      'MRawEvtHeader.fStereoEvtNumber', 'MRawEvtHeader.fDAQEvtNumber',
                      'MTime.fMjd', 'MTime.fTime.fMilliSec', 'MTime.fNanoSec'
                      ]

        pointing_array_list = ['MPointingPos.fZd', 'MPointingPos.fAz', 'MPointingPos.fRa', 'MPointingPos.fDec']
        drive_array_list = ['MReportDrive.fMjd', 'MReportDrive.fCurrentZd', 'MReportDrive.fCurrentAz',
                            'MReportDrive.fRa', 'MReportDrive.fDec']

961
962
963
964
        for file_name in file_list:

            input_file = uproot.open(file_name)

965
966
            events = input_file['Events'].arrays(array_list)

967
            # Reading the info common to MC and real data
968
969
970
971
            charge = events[b'MCerPhotEvt.fPixels.fPhot']
            arrival_time = events[b'MArrivalTime.fData']
            trigger_pattern = events[b'MTriggerPattern.fPrescaled']
            stereo_event_number = events[b'MRawEvtHeader.fStereoEvtNumber']
972
973

            # Computing the event arrival time
974
975
976
            mjd = events[b'MTime.fMjd']
            millisec = events[b'MTime.fTime.fMilliSec']
            nanosec = events[b'MTime.fNanoSec']
977
978
979
980
981

            mjd = mjd + (millisec / 1e3 + nanosec / 1e9) / 86400.0

            if 'MPointingPos.' in input_file['Events']:
                # Retrieving the telescope pointing direction
982
983
984
985
986
987
                pointing = input_file['Events'].arrays(pointing_array_list)

                pointing_zd = pointing[b'MPointingPos.fZd']
                pointing_az = pointing[b'MPointingPos.fAz']
                pointing_ra = pointing[b'MPointingPos.fRa']
                pointing_dec = pointing[b'MPointingPos.fDec']
988
989
            else:
                # Getting the telescope drive info
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
                drive = input_file['Drive'].arrays(drive_array_list)

                drive_mjd = drive[b'MReportDrive.fMjd']
                drive_zd = drive[b'MReportDrive.fCurrentZd']
                drive_az = drive[b'MReportDrive.fCurrentAz']
                drive_ra = drive[b'MReportDrive.fRa']
                drive_dec = drive[b'MReportDrive.fDec']

                # Finding only non-repeating drive entries
                # Repeating entries lead to failure in pointing interpolation
                non_repeating = scipy.diff(drive_mjd) > 0
                non_repeating = scipy.concatenate((non_repeating, [True]))

                # Filtering out the repeating ones
                drive_mjd = drive_mjd[non_repeating]
                drive_zd = drive_zd[non_repeating]
                drive_az = drive_az[non_repeating]
                drive_ra = drive_ra[non_repeating]
                drive_dec = drive_dec[non_repeating]

                if len(drive_zd) > 2:
                    # If there are enough drive data, creating azimuth and zenith angles interpolators
                    drive_zd_pointing_interpolator = scipy.interpolate.interp1d(drive_mjd, drive_zd, fill_value="extrapolate")
                    drive_az_pointing_interpolator = scipy.interpolate.interp1d(drive_mjd, drive_az, fill_value="extrapolate")

                    # Creating azimuth and zenith angles interpolators
                    drive_ra_pointing_interpolator = scipy.interpolate.interp1d(drive_mjd, drive_ra, fill_value="extrapolate")
                    drive_dec_pointing_interpolator = scipy.interpolate.interp1d(drive_mjd, drive_dec, fill_value="extrapolate")

                    # Interpolating the drive pointing to the event time stamps
                    pointing_zd = drive_zd_pointing_interpolator(mjd)
                    pointing_az = drive_az_pointing_interpolator(mjd)
                    pointing_ra = drive_ra_pointing_interpolator(mjd)
                    pointing_dec = drive_dec_pointing_interpolator(mjd)

                else:
                    # Not enough data to interpolate the pointing direction.
                    pointing_zd = scipy.repeat(-1, len(mjd))
                    pointing_az = scipy.repeat(-1, len(mjd))
                    pointing_ra = scipy.repeat(-1, len(mjd))
                    pointing_dec = scipy.repeat(-1, len(mjd))
1031
1032
1033

            event_data['charge'].append(charge)
            event_data['arrival_time'].append(arrival_time)
1034
1035
1036
1037
1038
1039
1040
1041
            event_data['trigger_pattern'] = scipy.concatenate((event_data['trigger_pattern'], trigger_pattern))
            event_data['stereo_event_number'] = scipy.concatenate((event_data['stereo_event_number'], stereo_event_number))
            event_data['pointing_zd'] = scipy.concatenate((event_data['pointing_zd'], pointing_zd))
            event_data['pointing_az'] = scipy.concatenate((event_data['pointing_az'], pointing_az))
            event_data['pointing_ra'] = scipy.concatenate((event_data['pointing_ra'], pointing_ra))
            event_data['pointing_dec'] = scipy.concatenate((event_data['pointing_dec'], pointing_dec))

            event_data['MJD'] = scipy.concatenate((event_data['MJD'], mjd))
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243

            event_data['file_edges'].append(len(event_data['trigger_pattern']))

        return event_data

    def _find_pedestal_events(self):
        """
        This internal method identifies the IDs (order numbers) of the
        pedestal events in the run.

        Returns
        -------
        dict:
            A dictionary of pedestal event IDs in M1/2 separately.
        """

        pedestal_ids = dict()

        pedestal_trigger_pattern = 8

        for telescope in self.event_data:
            ped_triggers = np.where(self.event_data[telescope]['trigger_pattern'] == pedestal_trigger_pattern)
            pedestal_ids[telescope] = ped_triggers[0]

        return pedestal_ids

    def _find_stereo_events(self):
        """
        This internal methods identifies stereo events in the run.

        Returns
        -------
        list:
            A list of pairs (M1_id, M2_id) corresponding to stereo events in the run.
        """

        data_trigger_pattern = 128

        m2_data_condition = (self.event_data['M2']['trigger_pattern'] == data_trigger_pattern)

        stereo_ids = []
        n_m1_events = len(self.event_data['M1']['stereo_event_number'])

        for m1_id in range(0, n_m1_events):
            if self.event_data['M1']['trigger_pattern'][m1_id] == data_trigger_pattern:
                m2_stereo_condition = (self.event_data['M2']['stereo_event_number'] ==
                                       self.event_data['M1']['stereo_event_number'][m1_id])

                m12_match = np.where(m2_data_condition & m2_stereo_condition)

                if len(m12_match[0]) > 0:
                    stereo_pair = (m1_id, m12_match[0][0])
                    stereo_ids.append(stereo_pair)

        return stereo_ids

    def _find_mono_events(self):
        """
        This internal method identifies the IDs (order numbers) of the
        pedestal events in the run.

        Returns
        -------
        dict:
            A dictionary of pedestal event IDs in M1/2 separately.
        """

        mono_ids = dict()
        mono_ids['M1'] = []
        mono_ids['M2'] = []

        data_trigger_pattern = 128

        m1_data_condition = self.event_data['M1']['trigger_pattern'] == data_trigger_pattern
        m2_data_condition = self.event_data['M2']['trigger_pattern'] == data_trigger_pattern

        n_m1_events = len(self.event_data['M1']['stereo_event_number'])
        n_m2_events = len(self.event_data['M2']['stereo_event_number'])

        for m1_id in range(0, n_m1_events):
            if m1_data_condition[m1_id]:
                m2_stereo_condition = (self.event_data['M2']['stereo_event_number'] ==
                                       self.event_data['M1']['stereo_event_number'][m1_id])

                m12_match = np.where(m2_data_condition & m2_stereo_condition)

                if len(m12_match[0]) == 0:
                    mono_ids['M1'].append(m1_id)

        for m2_id in range(0, n_m2_events):
            if m2_data_condition[m2_id]:
                m1_stereo_condition = (self.event_data['M1']['stereo_event_number'] ==
                                       self.event_data['M2']['stereo_event_number'][m2_id])

                m12_match = np.where(m1_data_condition & m1_stereo_condition)

                if len(m12_match[0]) == 0:
                    mono_ids['M2'].append(m2_id)

        return mono_ids

    def _get_pedestal_file_num(self, pedestal_event_num, telescope):
        """
        This internal method identifies the M1/2 file number of the
        given pedestal event in M1/2 file lists, corresponding to this run.

        Parameters
        ----------
        pedestal_event_num: int
            Order number of the event in the list of pedestal events
            of the specified telescope, corresponding to this run.
        telescope: str
            The name of the telescope to which this event corresponds.
            May be "M1" or "M2".

        Returns
        -------
        file_num:
            Order number of the corresponding file in the M1 or M2 file list.
        """

        event_id = self.pedestal_ids[telescope][pedestal_event_num]
        file_num = np.digitize([event_id], self.event_data[telescope]['file_edges'])
        file_num = file_num[0] - 1

        return file_num

    def _get_stereo_file_num(self, stereo_event_num):
        """
        This internal method identifies the M1/2 file numbers of the
        given stereo event in M1/2 file lists, corresponding to this run.

        Parameters
        ----------
        stereo_event_num: int
            Order number of the event in the list of stereo events corresponding
            to this run.

        Returns
        -------
        m1_file_num:
            Order number of the corresponding file in the M1 file list.
        m2_file_num:
            Order number of the corresponding file in the M2 file list.
        """

        m1_id = self.stereo_ids[stereo_event_num][0]
        m2_id = self.stereo_ids[stereo_event_num][1]
        m1_file_num = np.digitize([m1_id], self.event_data['M1']['file_edges'])
        m2_file_num = np.digitize([m2_id], self.event_data['M2']['file_edges'])

        m1_file_num = m1_file_num[0] - 1
        m2_file_num = m2_file_num[0] - 1

        return m1_file_num, m2_file_num

    def _get_mono_file_num(self, mono_event_num, telescope):
        """
        This internal method identifies the M1/2 file number of the
        given mono event in M1/2 file lists, corresponding to this run.

        Parameters
        ----------
        mono_event_num: int
            Order number of the event in the list of stereo events corresponding
            to this run.
        telescope: str
            The name of the telescope to which this event corresponds.
            May be "M1" or "M2".

        Returns
        -------
        file_num:
            Order number of the corresponding file in the M1 or M2 file list.
        """

        event_id = self.mono_ids[telescope][mono_event_num]
        file_num = np.digitize([event_id], self.event_data[telescope]['file_edges'])
        file_num = file_num[0] - 1

        return file_num

    def get_pedestal_event_data(self, pedestal_event_num, telescope):
        """
        This method read the photon content and arrival time (per pixel)
        for the specified pedestal event. Also returned is the event telescope pointing
        data.

        Parameters
        ----------
        pedestal_event_num: int
            Order number of the event in the list of pedestal events for the
            given telescope, corresponding to this run.
        telescope: str
            The name of the telescope to which this event corresponds.
            May be "M1" or "M2".

        Returns
        -------
        dict:
            The output has the following structure:
            'image' - photon_content in requested telescope
Ievgen Vovk's avatar
Ievgen Vovk committed
1244
            'pulse_time' - arrival_times in requested telescope
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
            'pointing_az' - pointing azimuth
            'pointing_zd' - pointing zenith angle
            'pointing_ra' - pointing right ascension
            'pointing_dec' - pointing declination
            'mjd' - event arrival time
        """

        file_num = self._get_pedestal_file_num(pedestal_event_num, telescope)
        event_id = self.pedestal_ids[telescope][pedestal_event_num]

        id_in_file = event_id - self.event_data[telescope]['file_edges'][file_num]

        photon_content = self.event_data[telescope]['charge'][file_num][id_in_file][:self.n_camera_pixels]
        arrival_times = self.event_data[telescope]['arrival_time'][file_num][id_in_file][:self.n_camera_pixels]

        event_data = dict()
        event_data['image'] = photon_content
Ievgen Vovk's avatar
Ievgen Vovk committed
1262
        event_data['pulse_time'] = arrival_times
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
        event_data['pointing_az'] = self.event_data[telescope]['pointing_az'][event_id]
        event_data['pointing_zd'] = self.event_data[telescope]['pointing_zd'][event_id]
        event_data['pointing_ra'] = self.event_data[telescope]['pointing_ra'][event_id]
        event_data['pointing_dec'] = self.event_data[telescope]['pointing_dec'][event_id]
        event_data['mjd'] = self.event_data[telescope]['MJD'][event_id]

        return event_data

    def get_stereo_event_data(self, stereo_event_num):
        """
        This method read the photon content and arrival time (per pixel)
        for the specified stereo event. Also returned is the event telescope pointing
        data.

        Parameters
        ----------
        stereo_event_num: int
            Order number of the event in the list of stereo events corresponding
            to this run.

        Returns
        -------
        dict:
            The output has the following structure:
            'm1_image' - M1 photon_content
Ievgen Vovk's avatar
Ievgen Vovk committed
1288
            'm1_pulse_time' - M1 arrival_times
1289
            'm2_image' - M2 photon_content
Ievgen Vovk's avatar
Ievgen Vovk committed
1290
            'm2_pulse_time' - M2 arrival_times
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
            'pointing_az' - pointing azimuth
            'pointing_zd' - pointing zenith angle
            'pointing_ra' - pointing right ascension
            'pointing_dec' - pointing declination
            'mjd' - event arrival time
        """

        m1_file_num, m2_file_num = self._get_stereo_file_num(stereo_event_num)
        m1_id = self.stereo_ids[stereo_event_num][0]
        m2_id = self.stereo_ids[stereo_event_num][1]

        m1_id_in_file = m1_id - self.event_data['M1']['file_edges'][m1_file_num]
        m2_id_in_file = m2_id - self.event_data['M2']['file_edges'][m2_file_num]

        m1_photon_content = self.event_data['M1']['charge'][m1_file_num][m1_id_in_file][:self.n_camera_pixels]
        m1_arrival_times = self.event_data['M1']['arrival_time'][m1_file_num][m1_id_in_file][:self.n_camera_pixels]

        m2_photon_content = self.event_data['M2']['charge'][m2_file_num][m2_id_in_file][:self.n_camera_pixels]
        m2_arrival_times = self.event_data['M2']['arrival_time'][m2_file_num][m2_id_in_file][:self.n_camera_pixels]

        event_data = dict()
        event_data['m1_image'] = m1_photon_content
Ievgen Vovk's avatar
Ievgen Vovk committed
1313
        event_data['m1_pulse_time'] = m1_arrival_times
1314
        event_data['m2_image'] = m2_photon_content
Ievgen Vovk's avatar
Ievgen Vovk committed
1315
        event_data['m2_pulse_time'] = m2_arrival_times
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
        event_data['pointing_az'] = self.event_data['M1']['pointing_az'][m1_id]
        event_data['pointing_zd'] = self.event_data['M1']['pointing_zd'][m1_id]
        event_data['pointing_ra'] = self.event_data['M1']['pointing_ra'][m1_id]
        event_data['pointing_dec'] = self.event_data['M1']['pointing_dec'][m1_id]
        event_data['mjd'] = self.event_data['M1']['MJD'][m1_id]

        return event_data

    def get_mono_event_data(self, mono_event_num, telescope):
        """
        This method read the photon content and arrival time (per pixel)
        for the specified mono event. Also returned is the event telescope pointing
        data.

        Parameters
        ----------
        mono_event_num: int
            Order number of the event in the list of mono events for the
            given telescope, corresponding to this run.
        telescope: str
            The name of the telescope to which this event corresponds.
            May be "M1" or "M2".

        Returns
        -------
        dict:
            The output has the following structure:
            'image' - photon_content in requested telescope
Ievgen Vovk's avatar
Ievgen Vovk committed
1344
            'pulse_time' - arrival_times in requested telescope
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
            'pointing_az' - pointing azimuth
            'pointing_zd' - pointing zenith angle
            'pointing_ra' - pointing right ascension
            'pointing_dec' - pointing declination
            'mjd' - event arrival time
        """

        file_num = self._get_mono_file_num(mono_event_num, telescope)
        event_id = self.mono_ids[telescope][mono_event_num]

        id_in_file = event_id - self.event_data[telescope]['file_edges'][file_num]

        photon_content = self.event_data[telescope]['charge'][file_num][id_in_file][:self.n_camera_pixels]
        arrival_times = self.event_data[telescope]['arrival_time'][file_num][id_in_file][:self.n_camera_pixels]

        event_data = dict()
        event_data['image'] = photon_content
Ievgen Vovk's avatar
Ievgen Vovk committed
1362
        event_data['pulse_time'] = arrival_times
1363
1364
1365
1366
1367
1368
1369
        event_data['pointing_az'] = self.event_data[telescope]['pointing_az'][event_id]
        event_data['pointing_zd'] = self.event_data[telescope]['pointing_zd'][event_id]
        event_data['pointing_ra'] = self.event_data[telescope]['pointing_ra'][event_id]
        event_data['pointing_dec'] = self.event_data[telescope]['pointing_dec'][event_id]
        event_data['mjd'] = self.event_data[telescope]['MJD'][event_id]

        return event_data
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390


class MarsMCFile:
    """
    This class implements reading of the event data from a MAGIC MC data file.
    It actually corresponds to 2 MC files - for M1 and M2 telescopes, which are
    treated as one as the stored events are recorded in the stereo mode.
    """

    def __init__(self, mc_file):
        """
        Constructor of the class. Defines the run to use and the camera pixel arrangement.

        Parameters
        ----------
        mc_file: str
            MC file path.
        """

        self.mc_file = mc_file
        self.telescope = int(re.findall(r'.*_M(\d)_.*', mc_file)[0])
1391
1392
1393
1394
1395
1396
1397
        try:
            self.run_number = int(re.findall(r'.*_M\d_za\d+to\d+_\d_(\d+)_Y_.*', mc_file)[0])
        except:
            self.run_number = int(re.findall(r'.*_M\d_\d_(\d+)_.*', mc_file)[0])
            msg = "WARNING: MAGIC MC file format is unusual."
            print(msg)
        
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
        # Reading the event data
        self.event_data = self.load_events(self.mc_file)

        # Detecting pedestal events
        self.pedestal_ids = self._find_pedestal_events()
        # Detecting mono events
        self.mono_ids = self._find_mono_events()

        self.n_camera_pixels = 1039

    @property
    def n_events(self):
        return len(self.event_data['stereo_event_number'])

    @property
    def n_mono_events(self):
        return len(self.mono_ids)

    @property
    def n_pedestal_events(self):
        return len(self.pedestal_ids)

    @staticmethod
    def load_events(file_name):
        """
        This method loads events from the pre-defiled file and returns them as a dictionary.

        Parameters
        ----------
        file_name: str
            Name of the MAGIC calibrated file to use.

        Returns
        -------
        dict:
            A dictionary with the even properties: charge / arrival time data, trigger, direction etc.
        """

        try:
            import uproot
        except ImportError:
            msg = "The `uproot` python module is required to access the MAGIC data"
            raise ImportError(msg)

        event_data = dict()

        array_list = ['MCerPhotEvt.fPixels.fPhot', 'MArrivalTime.fData',
                      'MTriggerPattern.fPrescaled',
                      'MRawEvtHeader.fStereoEvtNumber', 'MRawEvtHeader.fDAQEvtNumber',
                      'MPointingPos.fZd', 'MPointingPos.fAz', 'MPointingPos.fRa', 'MPointingPos.fDec',
1448
1449
                      'MMcEvt.fEnergy', 'MMcEvt.fTheta', 'MMcEvt.fPhi', 'MMcEvt.fPartId',
                      'MMcEvt.fZFirstInteraction', 'MMcEvt.fCoreX', 'MMcEvt.fCoreY', 
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
                      ]

        names_mapping = {
            b'MCerPhotEvt.fPixels.fPhot': 'charge',
            b'MArrivalTime.fData': 'arrival_time',
            b'MTriggerPattern.fPrescaled': 'trigger_pattern',
            b'MRawEvtHeader.fStereoEvtNumber': 'stereo_event_number',
            b'MRawEvtHeader.fDAQEvtNumber': 'daq_event_number',
            b'MPointingPos.fZd': 'pointing_zd',
            b'MPointingPos.fAz': 'pointing_az',
            b'MPointingPos.fRa': 'pointing_ra',
            b'MPointingPos.fDec': 'pointing_dec',
            b'MMcEvt.fEnergy': 'true_energy',
            b'MMcEvt.fTheta': 'true_zd',
1464
1465
1466
1467
1468
            b'MMcEvt.fPhi': 'true_az',
            b'MMcEvt.fPartId': 'true_shower_primary_id',
            b'MMcEvt.fZFirstInteraction': 'true_h_first_int',
            b'MMcEvt.fCoreX': 'true_core_x',
            b'MMcEvt.fCoreY': 'true_core_y'
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
        }

        with uproot.open(file_name) as input_file:
            if 'Events' in input_file:
                data = input_file['Events'].arrays(array_list)

                # Mapping the read structure to the alternative names
                for key in data:
                    name = names_mapping[key]
                    event_data[name] = data[key]

                # Post processing
                event_data['true_zd'] = scipy.degrees(event_data['true_zd'])
                event_data['true_az'] = scipy.degrees(event_data['true_az'])
                # Transformation from Monte Carlo to usual azimuth
                event_data['true_az'] = -1 * (event_data['true_az'] - 180 + 7)

            else:
                # The file is likely corrupted, so return empty arrays
                for key in names_mapping:
                    name = names_mapping[key]
                    event_data[name] = scipy.zeros(0)

        return event_data

    def _find_pedestal_events(self):
        """
        This internal method identifies the IDs (order numbers) of the
        pedestal events in the run.

        Returns
        -------
        dict:
            A dictionary of pedestal event IDs in M1/2 separately.
        """

        ped_triggers = scipy.where(self.event_data['trigger_pattern'] == 8)
        pedestal_ids = ped_triggers[0]

        return pedestal_ids

    def _find_mono_events(self):
        """
        This internal method identifies the IDs (order numbers) of the
        pedestal events in the run.

        Returns
        -------
        dict:
            A dictionary of pedestal event IDs.
        """

        mono_triggers = scipy.where(self.event_data['trigger_pattern'] == 1)
        mono_ids = mono_triggers[0]

        return mono_ids

    def get_pedestal_event_data(self, pedestal_event_num):
        """
        This method read the photon content and arrival time (per pixel)
        for the specified pedestal event. Also returned is the event telescope pointing
        data.

        Parameters
        ----------
        pedestal_event_num: int
            Order number of the event in the list of pedestal events for the
            given telescope, corresponding to this run.

        Returns
        -------
        event_data: dict
            A dictionary containing the event image, arrival time map, true energy and Az/Zd
            as well as Az/Zd of the current telescope pointing.
        """

        event_id = self.pedestal_ids[pedestal_event_num]

        photon_content = self.event_data['charge'][event_id][:self.n_camera_pixels]
        arrival_times = self.event_data['arrival_time'][event_id][:self.n_camera_pixels]

        pointing_zd = self.event_data['pointing_zd'][event_id]
        pointing_az = self.event_data['pointing_az'][event_id]

        event_true_energy = 0.0
        event_true_zd = 0.0
        event_true_az = 0.0
1556
1557
1558
1559
        event_true_shower_primary_id = 0.0
        event_true_h_first_int = 0.0
        event_true_core_x = 0.0
        event_true_core_y = 0.0
1560
1561
1562
1563
1564
1565
1566
1567
1568

        event_data = dict()
        event_data['image'] = photon_content
        event_data['pulse_time'] = arrival_times
        event_data['pointing_az'] = pointing_az
        event_data['pointing_zd'] = pointing_zd
        event_data['true_energy'] = event_true_energy
        event_data['true_az'] = event_true_az
        event_data['true_zd'] = event_true_zd
1569
1570
1571
1572
        event_data['true_shower_primary_id'] = event_true_shower_primary_id
        event_data['true_h_first_int'] = event_true_h_first_int
        event_data['true_core_x'] = event_true_core_x
        event_data['true_core_y'] = event_true_core_y
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602

        return event_data

    def get_mono_event_data(self, mono_event_num):
        """
        This method read the photon content and arrival time (per pixel)
        for the specified mono event. Also returned is the event telescope pointing
        data.

        Parameters
        ----------
        mono_event_num: int
            Order number of the event in the list of mono events for the
            given telescope, corresponding to this run.

        Returns
        -------
        event_data: dict
            A dictionary containing the event image, arrival time map, true energy and Az/Zd
            as well as Az/Zd of the current telescope pointing.
        """

        event_id = self.mono_ids[mono_event_num]

        photon_content = self.event_data['charge'][event_id][:self.n_camera_pixels]
        arrival_times = self.event_data['arrival_time'][event_id][:self.n_camera_pixels]

        pointing_zd = self.event_data['pointing_zd'][event_id]
        pointing_az = self.event_data['pointing_az'][event_id]

1603
1604
1605
        event_true_energy = self.event_data['true_energy'][event_id]
        event_true_zd = self.event_data['true_zd'][event_id]
        event_true_az = self.event_data['true_az'][event_id]
1606
1607
1608
1609
        event_true_shower_primary_id = self.event_data['true_shower_primary_id'][event_id]
        event_true_h_first_int = self.event_data['true_h_first_int'][event_id]
        event_true_core_x = self.event_data['true_core_x'][event_id]
        event_true_core_y = self.event_data['true_core_y'][event_id]
1610
1611
1612
1613
1614
1615
1616
1617
1618

        event_data = dict()
        event_data['image'] = photon_content
        event_data['pulse_time'] = arrival_times
        event_data['pointing_az'] = pointing_az
        event_data['pointing_zd'] = pointing_zd
        event_data['true_energy'] = event_true_energy
        event_data['true_az'] = event_true_az
        event_data['true_zd'] = event_true_zd
1619
1620
1621
1622
        event_data['true_shower_primary_id'] = event_true_shower_primary_id
        event_data['true_h_first_int'] = event_true_h_first_int
        event_data['true_core_x'] = event_true_core_x
        event_data['true_core_y'] = event_true_core_y
1623
1624

        return event_data