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