Jpp  16.0.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JAAnetToolkit.hh
Go to the documentation of this file.
1 #ifndef __JAANET__JAANETTOOLKIT__
2 #define __JAANET__JAANETTOOLKIT__
3 
4 #include <cstdlib>
5 #include <algorithm>
6 
14 
15 #include "JLang/JException.hh"
18 #include "JGeometry3D/JAxis3D.hh"
19 #include "JGeometry3D/JTrack3D.hh"
20 #include "JGeometry3D/JTrack3E.hh"
21 #include "JGeometry3D/JVertex3D.hh"
23 #include "JPhysics/JTimeRange.hh"
24 #include "JAAnet/JParticleTypes.hh"
25 
26 
27 /**
28  * \file
29  *
30  * Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
31  * \author mdejong
32  */
33 namespace JAANET {}
34 namespace JPP { using namespace JAANET; }
35 
36 /**
37  * Extensions to Evt data format.
38  */
39 namespace JAANET {
40 
50 
51  /**
52  * Index of rescaled weight.
53  */
54  static const int RESCALED_WEIGHT_INDEX = 3;
55 
56  /**
57  * AAnet shower fit reconstruction type.
58  */
59  static const int AASHOWER_RECONSTRUCTION_TYPE = 101;
60 
61  /**
62  * Enumeration of hit types based on km3 codes.
63  */
64  enum JHitType_t {
65  HIT_TYPE_MUON_DIRECT = +5, //!< Direct light from muon
66  HIT_TYPE_MUON_SCATTERED = -5, //!< Scattered light from muon
67  HIT_TYPE_DELTARAYS_DIRECT = +4, //!< Direct light from delta-rays
68  HIT_TYPE_DELTARAYS_SCATTERED = -4, //!< Scattered light from delta-rays
69  HIT_TYPE_BREMSSTRAHLUNG_DIRECT = +3, //!< Direct light from Bremsstrahlung
70  HIT_TYPE_BREMSSTRAHLUNG_SCATTERED = -3, //!< Scattered light from Bremsstrahlung
71  HIT_TYPE_SHOWER_DIRECT = +99, //!< Direct light from primary shower
72  HIT_TYPE_SHOWER_SCATTERED = -99, //!< Scattered light from primary shower
73  HIT_TYPE_NOISE = -1, //!< Random noise
74  HIT_TYPE_UNKNOWN = 0 //!< Unknown source
75  };
76 
77 
78  /**
79  * Get true time of hit.
80  *
81  * \param hit hit
82  * \return true time of hit
83  */
84  inline double getTime(const Hit& hit)
85  {
86  return hit.t;
87  }
88 
89 
90  /**
91  * Get true charge of hit.
92  *
93  * \param hit hit
94  * \return true charge of hit
95  */
96  inline double getNPE(const Hit& hit)
97  {
98  return hit.a;
99  }
100 
101 
102  /**
103  * Verify hit origin.
104  *
105  * \param hit hit
106  * \return true if noise; else false
107  */
108  inline bool is_noise(const Hit& hit)
109  {
110  return hit.type == HIT_TYPE_NOISE;
111  }
112 
113 
114  /**
115  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
116  * Note that the global event time in not taken into account.
117  *
118  * \param event event
119  * \return time range
120  */
121  inline JTimeRange getTimeRange(const Evt& event)
122  {
124 
125  for (std::vector<Hit>::const_iterator hit = event.mc_hits.begin(); hit != event.mc_hits.end(); ++hit) {
126  if (!is_noise(*hit)) {
127  time_range.include(getTime(*hit));
128  }
129  }
130 
131  return time_range;
132  }
133 
134 
135  /**
136  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
137  * The resulting time range is larger than or equal to the given time window.\n
138  * Note that the global event time in not taken into account.
139  *
140  * \param event event
141  * \param T_ns time window [ns]
142  * \return time range
143  */
144  inline JTimeRange getTimeRange(const Evt& event, const JTimeRange& T_ns)
145  {
146  JTimeRange time_range = getTimeRange(event);
147 
148  if (!time_range.is_valid()) {
149  time_range = T_ns;
150  }
151 
152  if (time_range.getLength() < T_ns.getLength()) {
153  time_range.setLength(T_ns.getLength());
154  }
155 
156  return time_range;
157  }
158 
159 
160  /**
161  * Get position.
162  *
163  * \param pos position
164  * \return position
165  */
166  inline JPosition3D getPosition(const Vec& pos)
167  {
168  return JPosition3D(pos.x, pos.y, pos.z);
169  }
170 
171 
172  /**
173  * Get position.
174  *
175  * \param pos position
176  * \return position
177  */
178  inline Vec getPosition(const JPosition3D& pos)
179  {
180  return Vec(pos.getX(), pos.getY(), pos.getZ());
181  }
182 
183 
184  /**
185  * Get position.
186  *
187  * \param track track
188  * \return position
189  */
190  inline JPosition3D getPosition(const Trk& track)
191  {
192  return getPosition(track.pos);
193  }
194 
195 
196  /**
197  * Get direction.
198  *
199  * \param dir direction
200  * \return direction
201  */
202  inline JDirection3D getDirection(const Vec& dir)
203  {
204  return JDirection3D(dir.x, dir.y, dir.z);
205  }
206 
207 
208  /**
209  * Get direction.
210  *
211  * \param dir direction
212  * \return direction
213  */
214  inline Vec getDirection(const JDirection3D& dir)
215  {
216  return Vec(dir.getDX(), dir.getDY(), dir.getDZ());
217  }
218 
219 
220  /**
221  * Get direction.
222  *
223  * \param track track
224  * \return direction
225  */
226  inline JDirection3D getDirection(const Trk& track)
227  {
228  return getDirection(track.dir);
229  }
230 
231 
232  /**
233  * Get axis.
234  *
235  * \param track track
236  * \return axis
237  */
238  inline JAxis3D getAxis(const Trk& track)
239  {
240  return JAxis3D(getPosition(track), getDirection(track));
241  }
242 
243 
244  /**
245  * Get track.
246  *
247  * \param track track
248  * \return track
249  */
250  inline JTrack3E getTrack(const Trk& track)
251  {
252  return JTrack3E(JTrack3D(getAxis(track), track.t), track.E);
253  }
254 
255 
256  /**
257  * Get vertex.
258  *
259  * \param track track
260  * \return vertex
261  */
262  inline JVertex3D getVertex(const Trk& track)
263  {
264  return JVertex3D(getPosition(track), track.t);
265  }
266 
267 
268  /**
269  * Get transformation.
270  *
271  * \param track track
272  * \return transformation
273  */
275  {
276  return JTransformation3D(getPosition(track), getDirection(track));
277  }
278 
279 
280  /**
281  * Get track information.
282  *
283  * \param track track
284  * \param index index
285  * \param value default value
286  * \return actual value
287  */
288  inline double getW(const Trk& track, const size_t index, const double value)
289  {
290  if (index < track.fitinf.size())
291  return track.fitinf[index];
292  else
293  return value;
294  }
295 
296 
297  /**
298  * Test whether given track is a photon or neutral pion.
299  *
300  * \param track track
301  * \return true if photon or neutral pion; else false
302  */
303  inline bool is_photon (const Trk& track) { return ( track.type == TRACK_TYPE_PHOTON ||
304  abs(track.type) == TRACK_TYPE_NEUTRAL_PION); }
305 
306  /**
307  * Test whether given track is a neutrino.
308  *
309  * \param track track
310  * \return true if neutrino; else false
311  */
312  inline bool is_neutrino(const Trk& track) { return (abs(track.type) == TRACK_TYPE_NUE ||
313  abs(track.type) == TRACK_TYPE_NUMU ||
314  abs(track.type) == TRACK_TYPE_NUTAU); }
315 
316  /**
317  * Test whether given track is a (anti-)electron.
318  *
319  * \param track track
320  * \return true if (anti-)electron; else false
321  */
322  inline bool is_electron(const Trk& track) { return (abs(track.type) == TRACK_TYPE_ELECTRON); }
323 
324  /**
325  * Test whether given track is a (anti-)muon.
326  *
327  * \param track track
328  * \return true if (anti-)muon; else false
329  */
330  inline bool is_muon (const Trk& track) { return (abs(track.type) == TRACK_TYPE_MUON); }
331 
332  /**
333  * Test whether given track is a (anti-)tau.
334  *
335  * \param track track
336  * \return true if (anti-)tau; else false
337  */
338  inline bool is_tau (const Trk& track) { return (abs(track.type) == TRACK_TYPE_TAU); }
339 
340  /**
341  * Test whether given track is a charged pion.
342  *
343  * \param track track
344  * \return true if charged pion; else false
345  */
346  inline bool is_pion (const Trk& track) { return (abs(track.type) == TRACK_TYPE_CHARGED_PION_PLUS); }
347 
348  /**
349  * Test whether given track is a (anti-)proton.
350  *
351  * \param track track
352  * \return true if (anti-)proton; else false
353  */
354  inline bool is_proton (const Trk& track) { return (abs(track.type) == TRACK_TYPE_PROTON); }
355 
356  /**
357  * Test whether given track is a (anti-)neutron.
358  *
359  * \param track track
360  * \return true if (anti-)neutron; else false
361  */
362  inline bool is_neutron (const Trk& track) { return (abs(track.type) == TRACK_TYPE_NEUTRON); }
363 
364  /**
365  * Test whether given track is a lepton
366  *
367  * \param track track
368  * \return true if lepton; else fails
369  */
370  inline bool is_lepton (const Trk& track) { return (is_neutrino(track) ||
371  is_electron(track) ||
372  is_muon (track) ||
373  is_tau (track)); }
374 
375  /**
376  * Test whether given track is a charged lepton
377  *
378  * \param track track
379  * \return true if lepton; else fails
380  */
381  inline bool is_charged_lepton (const Trk& track) { return (is_electron(track) ||
382  is_muon (track) ||
383  is_tau (track)); }
384 
385  /**
386  * Test whether given track is a hadron
387  *
388  * \param track track
389  * \return true if hadron; else fails
390  */
391  inline bool is_hadron (const Trk& track) { return (abs(track.type) != TRACK_TYPE_PHOTON &&
392  !is_lepton(track)); }
393 
394  /**
395  * Test whether given track is a meson (is hadron and third digit of PDG code is zero)
396  *
397  * \param track track
398  * \return true if hadron; else fails
399  */
400  inline bool is_meson (const Trk& track) { return (is_hadron(track) &&
401  ((int)(track.type / 1000)) % 10 == 0); }
402 
403  /**
404  * Test whether given track is a baryon (is hadron and third digit of PDG code is not zero)
405  *
406  * \param track track
407  * \return true if hadron; else fails
408  */
409  inline bool is_baryon (const Trk& track) { return (is_hadron(track) &&
410  ((int)(track.type / 1000)) % 10 != 0); }
411 
412  /**
413  * Test whether given track corresponds to given particle type.
414  *
415  * \param track track
416  * \param type particle type
417  * \return true if matches the corresponding particle; else false
418  */
419  inline bool has_particleID (const Trk& track, const int type) { return (track.type == type); }
420 
421  /**
422  * Test whether given event has an incoming neutrino.
423  *
424  * \param evt event
425  * \return true if neutrino; else false
426  */
427  inline bool has_neutrino(const Evt& evt)
428  {
429  return !evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0]);
430  }
431 
432  /**
433  * Get incoming neutrino.
434  *
435  * \param evt event
436  * \return neutrino
437  */
438  inline const Trk& get_neutrino(const Evt& evt)
439  {
440  if (has_neutrino(evt))
441  return evt.mc_trks[0];
442  else
443  THROW(JIndexOutOfRange, "JAANET::get_neutrino(): Event " << evt.id << " has no neutrino.");
444  }
445 
446  /**
447  * Count the number of electrons in a given event
448  *
449  * \param evt event
450  * \return number of electrons
451  */
452  inline int count_electrons(const Evt& evt)
453  {
454  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
455  }
456 
457  /**
458  * Test whether given event has an electron.
459  *
460  * \param evt event
461  * \return true if event has electron; else false
462  */
463  inline bool has_electron(const Evt& evt)
464  {
465  return count_electrons(evt) != 0;
466  }
467 
468  /**
469  * Get first electron from the event tracklist
470  *
471  * \param evt event
472  * \return electron
473  */
474  inline const Trk& get_electron(const Evt& evt)
475  {
476  if (count_electrons(evt) > 0)
477  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
478  else
479  THROW(JIndexOutOfRange, "This event has no electron.");
480  }
481 
482  /**
483  * Test whether given hit was produced by an electron
484  *
485  * Warning: This method will only work with the output of KM3Sim.
486  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
487  *
488  * \param hit hit
489  * \return true if hit was produced by an electron; else false
490  */
491  inline bool from_electron(const Hit& hit)
492  {
493  return hit.type == GEANT4_TYPE_ELECTRON || hit.type == GEANT4_TYPE_ANTIELECTRON;
494  }
495 
496  /**
497  * Count the number of muons in a given event
498  *
499  * \param evt event
500  * \return number of muons
501  */
502  inline int count_muons(const Evt& evt)
503  {
504  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
505  }
506 
507  /**
508  * Test whether given event has a muon.
509  *
510  * \param evt event
511  * \return true if event has muon; else false
512  */
513  inline bool has_muon(const Evt& evt)
514  {
515  return count_muons(evt) != 0;
516  }
517 
518  /**
519  * Get first muon from the event tracklist
520  *
521  * \param evt event
522  * \return muon
523  */
524  inline const Trk& get_muon(const Evt& evt)
525  {
526  if (count_muons(evt) > 0)
527  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
528  else
529  THROW(JIndexOutOfRange, "This event has no muon.");
530  }
531 
532  /**
533  * Test whether given hit was produced by a muon
534  *
535  * Warning: This method will only work with the output of KM3Sim.
536  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
537  *
538  * \param hit hit
539  * \return true if hit was produced by a muon; else false
540  */
541  inline bool from_muon(const Hit& hit)
542  {
543  return hit.type == GEANT4_TYPE_MUON || hit.type == GEANT4_TYPE_ANTIMUON;
544  }
545 
546  /**
547  * Count the number of taus in a given event
548  *
549  * \param evt event
550  * \return number of taus
551  */
552  inline int count_taus(const Evt& evt)
553  {
554  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
555  }
556 
557  /**
558  * Test whether given event has a tau.
559  *
560  * \param evt event
561  * \return true if event has tau; else false
562  */
563  inline bool has_tau(const Evt& evt)
564  {
565  return count_taus(evt) != 0;
566  }
567 
568  /**
569  * Get first tau from the event tracklist
570  *
571  * \param evt event
572  * \return tau
573  */
574  inline const Trk& get_tau(const Evt& evt)
575  {
576  if (count_taus(evt) > 0)
577  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
578  else
579  THROW(JIndexOutOfRange, "This event has no tau.");
580  }
581 
582  /**
583  * Test whether given hit was produced by a tau
584  *
585  * Warning: This method will only work with the output of KM3Sim.
586  * For JSirene or KM3, you should check if track.id == hit.origin instead.
587  *
588  * \param hit hit
589  * \return true if hit was produced by a tau; else false
590  */
591  inline bool from_tau(const Hit& hit)
592  {
593  return hit.type == GEANT4_TYPE_TAU || hit.type == GEANT4_TYPE_ANTITAU;
594  }
595 
596  /**
597  * Count the number of hadrons in a given event (not including neutral pions)
598  *
599  * \param evt event
600  * \return number of hadrons
601  */
602  inline int count_hadrons(const Evt& evt)
603  {
604  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_hadron);
605  }
606 
607  /**
608  * Test whether given hit was produced by a hadronic shower
609  *
610  * Warning: This method will only work with the output of KM3Sim.
611  * For JSirene or KM3, you should check if track.id == hit.origin instead.
612  *
613  * \param hit hit
614  * \return true if hit was produced by a hadron; else false
615  */
616  inline bool from_hadron(const Hit& hit)
617  {
618  return (!from_electron(hit) && !from_muon(hit) && !from_tau(hit));
619  }
620 
621  /**
622  * Auxiliary function to get average direction and total energy of a hadronic cascade
623  *
624  * \param evt event
625  * \return cascade (default if there are no hadrons)
626  */
627  inline Trk get_hadronic_cascade(const Evt& evt)
628  {
629  Trk cascade;
630 
631  for (std::vector<Trk>::const_iterator track = evt.mc_trks.begin(); track != evt.mc_trks.end(); ++track) {
632  if (is_hadron(*track)) {
633  cascade.E += track->E;
634  cascade.dir += track->dir * track->E;
635  cascade.t = track->t;
636  cascade.pos = track->pos;
637  }
638  }
639 
640  cascade.dir.normalize();
641 
642  return cascade;
643  }
644 
645  /**
646  * Add time offset to mc event;
647  * according to current definition, the absolute time of the event is defined by the track "t" attribute;
648  * this could change in the future if the global attribute mc_t is assigned to this purpose.
649  *
650  * \param evt event
651  * \param tOff time offset [ns]
652  */
653  inline void add_time_offset(Evt& evt, const double tOff)
654  {
655  for (std::vector<Trk>::iterator p = evt.mc_trks.begin(); p != evt.mc_trks.end(); p++) {
656  p->t += tOff;
657  }
658  }
659 }
660 
661 #endif
bool from_tau(const Hit &hit)
Test whether given hit was produced by a tau.
static const int RESCALED_WEIGHT_INDEX
Index of rescaled weight.
JVertex3D getVertex(const Trk &track)
Get vertex.
static const JRange< double, std::less< double > > DEFAULT_RANGE
Default range.
Definition: JRange.hh:555
bool is_electron(const Trk &track)
Test whether given track is a (anti-)electron.
int count_electrons(const Evt &evt)
Count the number of electrons in a given event.
JHitType_t
Enumeration of hit types based on km3 codes.
Exceptions.
Scattered light from primary shower.
bool is_proton(const Trk &track)
Test whether given track is a (anti-)proton.
const Trk & get_muon(const Evt &evt)
Get first muon from the event tracklist.
Data structure for direction in three dimensions.
Definition: JDirection3D.hh:33
JTrack3E getTrack(const Trk &track)
Get track.
Direct light from delta-rays.
double t
track time [ns] (when the particle is at pos )
Definition: Trk.hh:19
double z
Definition: Vec.hh:14
Scattered light from muon.
JTransformation3D getTransformation(const Trk &track)
Get transformation.
void add_time_offset(Evt &evt, const double tOff)
Add time offset to mc event; according to current definition, the absolute time of the event is defin...
bool has_neutrino(const Evt &evt)
Test whether given event has an incoming neutrino.
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
bool has_electron(const Evt &evt)
Test whether given event has an electron.
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [ns]).
int count_muons(const Evt &evt)
Count the number of muons in a given event.
Auxiliary methods for selection of reconstructed tracks.
Vec dir
track direction
Definition: Trk.hh:18
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:696
bool has_particleID(const Trk &track, const int type)
Test whether given track corresponds to given particle type.
bool is_charged_lepton(const Trk &track)
Test whether given track is a charged lepton.
Scattered light from Bremsstrahlung.
Direct light from Bremsstrahlung.
double getTime(const Hit &hit)
Get true time of hit.
bool is_baryon(const Trk &track)
Test whether given track is a baryon (is hadron and third digit of PDG code is not zero) ...
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:20
double y
Definition: Vec.hh:14
3D track with energy.
Definition: JTrack3E.hh:30
Direct light from primary shower.
static const int AASHOWER_RECONSTRUCTION_TYPE
AAnet shower fit reconstruction type.
bool is_noise(const Hit &hit)
Verify hit origin.
JTimeRange getTimeRange(const Evt &event)
Get time range (i.e. time between earliest and latest hit) of Monte Carlo event.
bool is_neutrino(const Trk &track)
Test whether given track is a neutrino.
Axis object.
Definition: JAxis3D.hh:38
bool is_lepton(const Trk &track)
Test whether given track is a lepton.
double getW(const Trk &track, const size_t index, const double value)
Get track information.
double a
hit amplitude (in p.e.)
Definition: Hit.hh:24
double x
Definition: Vec.hh:14
The Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition: Vec.hh:12
JDirection3D getDirection(const Vec &dir)
Get direction.
Trk get_hadronic_cascade(const Evt &evt)
Auxiliary function to get average direction and total energy of a hadronic cascade.
const Trk & get_electron(const Evt &evt)
Get first electron from the event tracklist.
double getDY() const
Get y direction.
Definition: JVersor3D.hh:106
double getDX() const
Get x direction.
Definition: JVersor3D.hh:95
JAxis3D getAxis(const Trk &track)
Get axis.
JPosition3D getPosition(const Vec &pos)
Get position.
std::vector< double > fitinf
place to store additional fit info, for jgandalf, see JFitParameters.hh
Definition: Trk.hh:31
double getY() const
Get y position.
Definition: JVector3D.hh:104
Vec & normalize()
Normalise this vector.
Definition: Vec.hh:159
bool is_photon(const Trk &track)
Test whether given track is a photon or neutral pion.
bool from_hadron(const Hit &hit)
Test whether given hit was produced by a hadronic shower.
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
bool is_pion(const Trk &track)
Test whether given track is a charged pion.
Vec pos
postion of the track at time t [m]
Definition: Trk.hh:17
int count_taus(const Evt &evt)
Count the number of taus in a given event.
Definition: Hit.hh:8
const Trk & get_tau(const Evt &evt)
Get first tau from the event tracklist.
Direct light from muon.
bool has_muon(const Evt &evt)
Test whether given event has a muon.
std::vector< Hit > mc_hits
MC: list of MC truth hits.
Definition: Evt.hh:45
double t
hit time (from tdc+calibration or MC truth)
Definition: Hit.hh:23
Scattered light from delta-rays.
Definition of particle types.
double getX() const
Get x position.
Definition: JVector3D.hh:94
bool has_tau(const Evt &evt)
Test whether given event has a tau.
int id
offline event identifier
Definition: Evt.hh:21
double getNPE(const Hit &hit)
Get true charge of hit.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
bool is_neutron(const Trk &track)
Test whether given track is a (anti-)neutron.
Exception for accessing an index in a collection that is outside of its range.
Definition: JException.hh:90
const Trk & get_neutrino(const Evt &evt)
Get incoming neutrino.
bool is_hadron(const Trk &track)
Test whether given track is a hadron.
bool is_tau(const Trk &track)
Test whether given track is a (anti-)tau.
bool from_electron(const Hit &hit)
Test whether given hit was produced by an electron.
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition: Trk.hh:14
int type
particle type or parametrisation used for hit (mc only)
Definition: Hit.hh:28
double getZ() const
Get z position.
Definition: JVector3D.hh:115
double getDZ() const
Get z direction.
Definition: JVersor3D.hh:117
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:46
bool is_meson(const Trk &track)
Test whether given track is a meson (is hadron and third digit of PDG code is zero) ...
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:19
int count_hadrons(const Evt &evt)
Count the number of hadrons in a given event (not including neutral pions)
bool from_muon(const Hit &hit)
Test whether given hit was produced by a muon.