Jpp  15.0.1-rc.1-highqe
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.getLength() < T_ns.getLength()) {
149  time_range.setLength(T_ns.getLength());
150  }
151 
152  return time_range;
153  }
154 
155 
156  /**
157  * Get position.
158  *
159  * \param pos position
160  * \return position
161  */
162  inline JPosition3D getPosition(const Vec& pos)
163  {
164  return JPosition3D(pos.x, pos.y, pos.z);
165  }
166 
167 
168  /**
169  * Get position.
170  *
171  * \param pos position
172  * \return position
173  */
174  inline Vec getPosition(const JPosition3D& pos)
175  {
176  return Vec(pos.getX(), pos.getY(), pos.getZ());
177  }
178 
179 
180  /**
181  * Get position.
182  *
183  * \param track track
184  * \return position
185  */
186  inline JPosition3D getPosition(const Trk& track)
187  {
188  return getPosition(track.pos);
189  }
190 
191 
192  /**
193  * Get direction.
194  *
195  * \param dir direction
196  * \return direction
197  */
198  inline JDirection3D getDirection(const Vec& dir)
199  {
200  return JDirection3D(dir.x, dir.y, dir.z);
201  }
202 
203 
204  /**
205  * Get direction.
206  *
207  * \param dir direction
208  * \return direction
209  */
210  inline Vec getDirection(const JDirection3D& dir)
211  {
212  return Vec(dir.getDX(), dir.getDY(), dir.getDZ());
213  }
214 
215 
216  /**
217  * Get direction.
218  *
219  * \param track track
220  * \return direction
221  */
222  inline JDirection3D getDirection(const Trk& track)
223  {
224  return getDirection(track.dir);
225  }
226 
227 
228  /**
229  * Get axis.
230  *
231  * \param track track
232  * \return axis
233  */
234  inline JAxis3D getAxis(const Trk& track)
235  {
236  return JAxis3D(getPosition(track), getDirection(track));
237  }
238 
239 
240  /**
241  * Get track.
242  *
243  * \param track track
244  * \return track
245  */
246  inline JTrack3E getTrack(const Trk& track)
247  {
248  return JTrack3E(JTrack3D(getAxis(track), track.t), track.E);
249  }
250 
251 
252  /**
253  * Get vertex.
254  *
255  * \param track track
256  * \return vertex
257  */
258  inline JVertex3D getVertex(const Trk& track)
259  {
260  return JVertex3D(getPosition(track), track.t);
261  }
262 
263 
264  /**
265  * Get transformation.
266  *
267  * \param track track
268  * \return transformation
269  */
271  {
272  return JTransformation3D(getPosition(track), getDirection(track));
273  }
274 
275 
276  /**
277  * Get track information.
278  *
279  * \param track track
280  * \param index index
281  * \param value default value
282  * \return actual value
283  */
284  inline double getW(const Trk& track, const size_t index, const double value)
285  {
286  if (index < track.fitinf.size())
287  return track.fitinf[index];
288  else
289  return value;
290  }
291 
292 
293  /**
294  * Test whether given track is a photon or neutral pion.
295  *
296  * \param track track
297  * \return true if photon or neutral pion; else false
298  */
299  inline bool is_photon (const Trk& track) { return ( track.type == TRACK_TYPE_PHOTON ||
300  abs(track.type) == TRACK_TYPE_NEUTRAL_PION); }
301 
302  /**
303  * Test whether given track is a neutrino.
304  *
305  * \param track track
306  * \return true if neutrino; else false
307  */
308  inline bool is_neutrino(const Trk& track) { return (abs(track.type) == TRACK_TYPE_NUE ||
309  abs(track.type) == TRACK_TYPE_NUMU ||
310  abs(track.type) == TRACK_TYPE_NUTAU); }
311 
312  /**
313  * Test whether given track is a (anti-)electron.
314  *
315  * \param track track
316  * \return true if (anti-)electron; else false
317  */
318  inline bool is_electron(const Trk& track) { return (abs(track.type) == TRACK_TYPE_ELECTRON); }
319 
320  /**
321  * Test whether given track is a (anti-)muon.
322  *
323  * \param track track
324  * \return true if (anti-)muon; else false
325  */
326  inline bool is_muon (const Trk& track) { return (abs(track.type) == TRACK_TYPE_MUON); }
327 
328  /**
329  * Test whether given track is a (anti-)tau.
330  *
331  * \param track track
332  * \return true if (anti-)tau; else false
333  */
334  inline bool is_tau (const Trk& track) { return (abs(track.type) == TRACK_TYPE_TAU); }
335 
336  /**
337  * Test whether given track is a charged pion.
338  *
339  * \param track track
340  * \return true if charged pion; else false
341  */
342  inline bool is_pion (const Trk& track) { return (abs(track.type) == TRACK_TYPE_CHARGED_PION_PLUS); }
343 
344  /**
345  * Test whether given track is a (anti-)proton.
346  *
347  * \param track track
348  * \return true if (anti-)proton; else false
349  */
350  inline bool is_proton (const Trk& track) { return (abs(track.type) == TRACK_TYPE_PROTON); }
351 
352  /**
353  * Test whether given track is a (anti-)neutron.
354  *
355  * \param track track
356  * \return true if (anti-)neutron; else false
357  */
358  inline bool is_neutron (const Trk& track) { return (abs(track.type) == TRACK_TYPE_NEUTRON); }
359 
360  /**
361  * Test whether given track is a lepton
362  *
363  * \param track track
364  * \return true if lepton; else fails
365  */
366  inline bool is_lepton (const Trk& track) { return (is_neutrino(track) ||
367  is_electron(track) ||
368  is_muon (track) ||
369  is_tau (track)); }
370 
371  /**
372  * Test whether given track is a charged lepton
373  *
374  * \param track track
375  * \return true if lepton; else fails
376  */
377  inline bool is_charged_lepton (const Trk& track) { return (is_electron(track) ||
378  is_muon (track) ||
379  is_tau (track)); }
380 
381  /**
382  * Test whether given track is a hadron
383  *
384  * \param track track
385  * \return true if hadron; else fails
386  */
387  inline bool is_hadron (const Trk& track) { return (abs(track.type) != TRACK_TYPE_PHOTON &&
388  !is_lepton(track)); }
389 
390  /**
391  * Test whether given track is a meson (is hadron and third digit of PDG code is zero)
392  *
393  * \param track track
394  * \return true if hadron; else fails
395  */
396  inline bool is_meson (const Trk& track) { return (is_hadron(track) &&
397  ((int)(track.type / 1000)) % 10 == 0); }
398 
399  /**
400  * Test whether given track is a baryon (is hadron and third digit of PDG code is not zero)
401  *
402  * \param track track
403  * \return true if hadron; else fails
404  */
405  inline bool is_baryon (const Trk& track) { return (is_hadron(track) &&
406  ((int)(track.type / 1000)) % 10 != 0); }
407 
408  /**
409  * Test whether given track corresponds to given particle type.
410  *
411  * \param track track
412  * \param type particle type
413  * \return true if matches the corresponding particle; else false
414  */
415  inline bool has_particleID (const Trk& track, const int type) { return (track.type == type); }
416 
417  /**
418  * Test whether given event has an incoming neutrino.
419  *
420  * \param evt event
421  * \return true if neutrino; else false
422  */
423  inline bool has_neutrino(const Evt& evt)
424  {
425  return !evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0]);
426  }
427 
428  /**
429  * Get incoming neutrino.
430  *
431  * \param evt event
432  * \return neutrino
433  */
434  inline const Trk& get_neutrino(const Evt& evt)
435  {
436  if (has_neutrino(evt))
437  return evt.mc_trks[0];
438  else
439  THROW(JIndexOutOfRange, "This event has no neutrino.");
440  }
441 
442  /**
443  * Count the number of electrons in a given event
444  *
445  * \param evt event
446  * \return number of electrons
447  */
448  inline int count_electrons(const Evt& evt)
449  {
450  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
451  }
452 
453  /**
454  * Test whether given event has an electron.
455  *
456  * \param evt event
457  * \return true if event has electron; else false
458  */
459  inline bool has_electron(const Evt& evt)
460  {
461  return count_electrons(evt) != 0;
462  }
463 
464  /**
465  * Get first electron from the event tracklist
466  *
467  * \param evt event
468  * \return electron
469  */
470  inline const Trk& get_electron(const Evt& evt)
471  {
472  if (count_electrons(evt) > 0)
473  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
474  else
475  THROW(JIndexOutOfRange, "This event has no electron.");
476  }
477 
478  /**
479  * Test whether given hit was produced by an electron
480  *
481  * Warning: This method will only work with the output of KM3Sim.
482  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
483  *
484  * \param hit hit
485  * \return true if hit was produced by an electron; else false
486  */
487  inline bool from_electron(const Hit& hit)
488  {
489  return hit.type == GEANT4_TYPE_ELECTRON || hit.type == GEANT4_TYPE_ANTIELECTRON;
490  }
491 
492  /**
493  * Count the number of muons in a given event
494  *
495  * \param evt event
496  * \return number of muons
497  */
498  inline int count_muons(const Evt& evt)
499  {
500  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
501  }
502 
503  /**
504  * Test whether given event has a muon.
505  *
506  * \param evt event
507  * \return true if event has muon; else false
508  */
509  inline bool has_muon(const Evt& evt)
510  {
511  return count_muons(evt) != 0;
512  }
513 
514  /**
515  * Get first muon from the event tracklist
516  *
517  * \param evt event
518  * \return muon
519  */
520  inline const Trk& get_muon(const Evt& evt)
521  {
522  if (count_muons(evt) > 0)
523  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
524  else
525  THROW(JIndexOutOfRange, "This event has no muon.");
526  }
527 
528  /**
529  * Test whether given hit was produced by a muon
530  *
531  * Warning: This method will only work with the output of KM3Sim.
532  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
533  *
534  * \param hit hit
535  * \return true if hit was produced by a muon; else false
536  */
537  inline bool from_muon(const Hit& hit)
538  {
539  return hit.type == GEANT4_TYPE_MUON || hit.type == GEANT4_TYPE_ANTIMUON;
540  }
541 
542  /**
543  * Count the number of taus in a given event
544  *
545  * \param evt event
546  * \return number of taus
547  */
548  inline int count_taus(const Evt& evt)
549  {
550  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
551  }
552 
553  /**
554  * Test whether given event has a tau.
555  *
556  * \param evt event
557  * \return true if event has tau; else false
558  */
559  inline bool has_tau(const Evt& evt)
560  {
561  return count_taus(evt) != 0;
562  }
563 
564  /**
565  * Get first tau from the event tracklist
566  *
567  * \param evt event
568  * \return tau
569  */
570  inline const Trk& get_tau(const Evt& evt)
571  {
572  if (count_taus(evt) > 0)
573  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
574  else
575  THROW(JIndexOutOfRange, "This event has no tau.");
576  }
577 
578  /**
579  * Test whether given hit was produced by a tau
580  *
581  * Warning: This method will only work with the output of KM3Sim.
582  * For JSirene or KM3, you should check if track.id == hit.origin instead.
583  *
584  * \param hit hit
585  * \return true if hit was produced by a tau; else false
586  */
587  inline bool from_tau(const Hit& hit)
588  {
589  return hit.type == GEANT4_TYPE_TAU || hit.type == GEANT4_TYPE_ANTITAU;
590  }
591 
592  /**
593  * Count the number of hadrons in a given event (not including neutral pions)
594  *
595  * \param evt event
596  * \return number of hadrons
597  */
598  inline int count_hadrons(const Evt& evt)
599  {
600  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_hadron);
601  }
602 
603  /**
604  * Test whether given hit was produced by a hadronic shower
605  *
606  * Warning: This method will only work with the output of KM3Sim.
607  * For JSirene or KM3, you should check if track.id == hit.origin instead.
608  *
609  * \param hit hit
610  * \return true if hit was produced by a hadron; else false
611  */
612  inline bool from_hadron(const Hit& hit)
613  {
614  return (!from_electron(hit) && !from_muon(hit) && !from_tau(hit));
615  }
616 
617  /**
618  * Auxiliary function to get average direction and total energy of a hadronic cascade
619  *
620  * \param evt event
621  * \return cascade (default if there are no hadrons)
622  */
623  inline Trk get_hadronic_cascade(const Evt& evt)
624  {
625  Trk cascade;
626 
627  for (std::vector<Trk>::const_iterator track = evt.mc_trks.begin(); track != evt.mc_trks.end(); ++track) {
628  if (is_hadron(*track)) {
629  cascade.E += track->E;
630  cascade.dir += track->dir * track->E;
631  cascade.t = track->t;
632  cascade.pos = track->pos;
633  }
634  }
635 
636  cascade.dir.normalize();
637 
638  return cascade;
639  }
640 
641  /**
642  * Add time offset to mc event;
643  * according to current definition, the absolute time of the event is defined by the track "t" attribute;
644  * this could change in the future if the global attribute mc_t is assigned to this purpose.
645  *
646  * \param evt event
647  * \param tOff time offset [ns]
648  */
649  inline void add_time_offset(Evt& evt, const double tOff)
650  {
651  for (std::vector<Trk>::iterator p = evt.mc_trks.begin(); p != evt.mc_trks.end(); p++) {
652  p->t += tOff;
653  }
654  }
655 }
656 
657 #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:568
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: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: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.
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.