Jpp  18.1.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 
15 
16 #include "JLang/JException.hh"
17 #include "JLang/JPredicate.hh"
18 
21 #include "JGeometry3D/JAxis3D.hh"
22 #include "JGeometry3D/JTrack3D.hh"
23 #include "JGeometry3D/JTrack3E.hh"
24 #include "JGeometry3D/JVertex3D.hh"
26 
27 #include "JPhysics/JConstants.hh"
28 #include "JPhysics/JTimeRange.hh"
29 #include "JPhysics/JGeane.hh"
30 
31 #include "JAAnet/JParticleTypes.hh"
32 #include "JAAnet/JPDB.hh"
33 
34 
35 /**
36  * \file
37  *
38  * Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
39  * \author mdejong
40  */
41 namespace JAANET {}
42 namespace JPP { using namespace JAANET; }
43 
44 /**
45  * Extensions to Evt data format.
46  */
47 namespace JAANET {
48 
52 
60 
64 
65  /**
66  * AAnet shower fit reconstruction type.
67  */
68  static const int AASHOWER_RECONSTRUCTION_TYPE = 101;
69 
70  /**
71  * Enumeration of hit types based on km3 codes.
72  */
73  enum JHitType_t {
74  HIT_TYPE_MUON_DIRECT = +5, //!< Direct light from muon
75  HIT_TYPE_MUON_SCATTERED = -5, //!< Scattered light from muon
76  HIT_TYPE_DELTARAYS_DIRECT = +4, //!< Direct light from delta-rays
77  HIT_TYPE_DELTARAYS_SCATTERED = -4, //!< Scattered light from delta-rays
78  HIT_TYPE_BREMSSTRAHLUNG_DIRECT = +3, //!< Direct light from Bremsstrahlung
79  HIT_TYPE_BREMSSTRAHLUNG_SCATTERED = -3, //!< Scattered light from Bremsstrahlung
80  HIT_TYPE_SHOWER_DIRECT = +99, //!< Direct light from primary shower
81  HIT_TYPE_SHOWER_SCATTERED = -99, //!< Scattered light from primary shower
82  HIT_TYPE_NOISE = -1, //!< Random noise
83  HIT_TYPE_UNKNOWN = 0 //!< Unknown source
84  };
85 
86 
87  /**
88  * Get true time of hit.
89  *
90  * \param hit hit
91  * \return true time of hit
92  */
93  inline double getTime(const Hit& hit)
94  {
95  return hit.t;
96  }
97 
98 
99  /**
100  * Get true charge of hit.
101  *
102  * \param hit hit
103  * \return true charge of hit
104  */
105  inline double getNPE(const Hit& hit)
106  {
107  return hit.a;
108  }
109 
110 
111  /**
112  * Verify hit origin.
113  *
114  * \param hit hit
115  * \return true if noise; else false
116  */
117  inline bool is_noise(const Hit& hit)
118  {
119  return hit.type == HIT_TYPE_NOISE;
120  }
121 
122 
123  /**
124  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
125  * Note that the global event time in not taken into account.
126  *
127  * \param event event
128  * \return time range
129  */
130  inline JTimeRange getTimeRange(const Evt& event)
131  {
133 
134  for (std::vector<Hit>::const_iterator hit = event.mc_hits.begin(); hit != event.mc_hits.end(); ++hit) {
135  if (!is_noise(*hit)) {
136  time_range.include(getTime(*hit));
137  }
138  }
139 
140  return time_range;
141  }
142 
143 
144  /**
145  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
146  * The resulting time range is larger than or equal to the given time window.\n
147  * Note that the global event time in not taken into account.
148  *
149  * \param event event
150  * \param T_ns time window [ns]
151  * \return time range
152  */
153  inline JTimeRange getTimeRange(const Evt& event, const JTimeRange& T_ns)
154  {
155  JTimeRange time_range = getTimeRange(event);
156 
157  if (!time_range.is_valid()) {
158  time_range = T_ns;
159  }
160 
161  if (time_range.getLength() < T_ns.getLength()) {
162  time_range.setLength(T_ns.getLength());
163  }
164 
165  return time_range;
166  }
167 
168 
169  /**
170  * Get position.
171  *
172  * \param pos position
173  * \return position
174  */
175  inline JPosition3D getPosition(const Vec& pos)
176  {
177  return JPosition3D(pos.x, pos.y, pos.z);
178  }
179 
180 
181  /**
182  * Get position.
183  *
184  * \param pos position
185  * \return position
186  */
187  inline Vec getPosition(const JPosition3D& pos)
188  {
189  return Vec(pos.getX(), pos.getY(), pos.getZ());
190  }
191 
192 
193  /**
194  * Get position.
195  *
196  * \param track track
197  * \return position
198  */
199  inline JPosition3D getPosition(const Trk& track)
200  {
201  return getPosition(track.pos);
202  }
203 
204 
205  /**
206  * Get direction.
207  *
208  * \param dir direction
209  * \return direction
210  */
211  inline JDirection3D getDirection(const Vec& dir)
212  {
213  return JDirection3D(dir.x, dir.y, dir.z);
214  }
215 
216 
217  /**
218  * Get direction.
219  *
220  * \param dir direction
221  * \return direction
222  */
223  inline Vec getDirection(const JDirection3D& dir)
224  {
225  return Vec(dir.getDX(), dir.getDY(), dir.getDZ());
226  }
227 
228 
229  /**
230  * Get direction.
231  *
232  * \param track track
233  * \return direction
234  */
235  inline JDirection3D getDirection(const Trk& track)
236  {
237  return getDirection(track.dir);
238  }
239 
240 
241  /**
242  * Get axis.
243  *
244  * \param track track
245  * \return axis
246  */
247  inline JAxis3D getAxis(const Trk& track)
248  {
249  return JAxis3D(getPosition(track), getDirection(track));
250  }
251 
252 
253  /**
254  * Get track.
255  *
256  * \param track track
257  * \return track
258  */
259  inline JTrack3E getTrack(const Trk& track)
260  {
261  return JTrack3E(JTrack3D(getAxis(track), track.t), track.E);
262  }
263 
264 
265  /**
266  * Get vertex.
267  *
268  * \param track track
269  * \return vertex
270  */
271  inline JVertex3D getVertex(const Trk& track)
272  {
273  return JVertex3D(getPosition(track), track.t);
274  }
275 
276 
277  /**
278  * Get transformation.
279  *
280  * \param track track
281  * \return transformation
282  */
284  {
285  return JTransformation3D(getPosition(track), getDirection(track));
286  }
287 
288 
289  /**
290  * Get track information.
291  *
292  * \param track track
293  * \param index index
294  * \param value default value
295  * \return actual value
296  */
297  inline double getW(const Trk& track, const size_t index, const double value)
298  {
299  if (index < track.fitinf.size())
300  return track.fitinf[index];
301  else
302  return value;
303  }
304 
305 
306  /**
307  * Test whether given track is a photon or neutral pion.
308  *
309  * \param track track
310  * \return true if photon or neutral pion; else false
311  */
312  inline bool is_photon (const Trk& track) { return ( track.type == TRACK_TYPE_PHOTON ||
313  abs(track.type) == TRACK_TYPE_NEUTRAL_PION); }
314 
315  /**
316  * Test whether given track is a neutrino.
317  *
318  * \param track track
319  * \return true if neutrino; else false
320  */
321  inline bool is_neutrino(const Trk& track) { return (abs(track.type) == TRACK_TYPE_NUE ||
322  abs(track.type) == TRACK_TYPE_NUMU ||
323  abs(track.type) == TRACK_TYPE_NUTAU); }
324 
325  /**
326  * Test whether given track is a (anti-)electron.
327  *
328  * \param track track
329  * \return true if (anti-)electron; else false
330  */
331  inline bool is_electron(const Trk& track) { return (abs(track.type) == TRACK_TYPE_ELECTRON); }
332 
333  /**
334  * Test whether given track is a (anti-)muon.
335  *
336  * \param track track
337  * \return true if (anti-)muon; else false
338  */
339  inline bool is_muon (const Trk& track) { return (abs(track.type) == TRACK_TYPE_MUON); }
340 
341  /**
342  * Test whether given track is a (anti-)tau.
343  *
344  * \param track track
345  * \return true if (anti-)tau; else false
346  */
347  inline bool is_tau (const Trk& track) { return (abs(track.type) == TRACK_TYPE_TAU); }
348 
349  /**
350  * Test whether given track is a charged pion.
351  *
352  * \param track track
353  * \return true if charged pion; else false
354  */
355  inline bool is_pion (const Trk& track) { return (abs(track.type) == TRACK_TYPE_CHARGED_PION_PLUS); }
356 
357  /**
358  * Test whether given track is a (anti-)proton.
359  *
360  * \param track track
361  * \return true if (anti-)proton; else false
362  */
363  inline bool is_proton (const Trk& track) { return (abs(track.type) == TRACK_TYPE_PROTON); }
364 
365  /**
366  * Test whether given track is a (anti-)neutron.
367  *
368  * \param track track
369  * \return true if (anti-)neutron; else false
370  */
371  inline bool is_neutron (const Trk& track) { return (abs(track.type) == TRACK_TYPE_NEUTRON); }
372 
373  /**
374  * Test whether given track is a lepton
375  *
376  * \param track track
377  * \return true if lepton; else fails
378  */
379  inline bool is_lepton (const Trk& track) { return (is_neutrino(track) ||
380  is_electron(track) ||
381  is_muon (track) ||
382  is_tau (track)); }
383 
384  /**
385  * Test whether given track is a charged lepton
386  *
387  * \param track track
388  * \return true if lepton; else fails
389  */
390  inline bool is_charged_lepton (const Trk& track) { return (is_electron(track) ||
391  is_muon (track) ||
392  is_tau (track)); }
393 
394  /**
395  * Test whether given track is a hadron
396  *
397  * \param track track
398  * \return true if hadron; else fails
399  */
400  inline bool is_hadron (const Trk& track) { return (abs(track.type) != TRACK_TYPE_PHOTON &&
401  !is_lepton(track)); }
402 
403  /**
404  * Test whether given track is a meson (is hadron and third digit of PDG code is zero)
405  *
406  * \param track track
407  * \return true if hadron; else fails
408  */
409  inline bool is_meson (const Trk& track) { return (is_hadron(track) &&
410  ((int)(track.type / 1000)) % 10 == 0); }
411 
412  /**
413  * Test whether given track is a baryon (is hadron and third digit of PDG code is not zero)
414  *
415  * \param track track
416  * \return true if hadron; else fails
417  */
418  inline bool is_baryon (const Trk& track) { return (is_hadron(track) &&
419  ((int)(track.type / 1000)) % 10 != 0); }
420 
421  /**
422  * Test whether given track corresponds to given particle type.
423  *
424  * \param track track
425  * \param type particle type
426  * \return true if matches the corresponding particle; else false
427  */
428  inline bool has_particleID (const Trk& track, const int type) { return (track.type == type); }
429 
430  /**
431  * Test whether given track corresponds to an initial state particle.
432  *
433  * \param track track
434  * \return true if track corresponds to an initial state particle; else false
435  */
436  inline bool is_initialstate(const Trk& track)
437  {
438  if (Evt::ROOT_IO_VERSION >= 14) {
439 
440  return (track.status == TRK_ST_PRIMARYNEUTRINO ||
441  track.status == TRK_ST_PRIMARYCOSMIC ||
442  track.status == TRK_ST_MUONBUNDLE ||
443  track.status == TRK_ST_ININUCLEI);
444 
445  } else if (Evt::ROOT_IO_VERSION >= 11) {
446 
447  return (track.mother_id == TRK_MOTHER_UNDEFINED ||
448  track.mother_id == TRK_MOTHER_NONE);
449 
450  } else {
451 
452  return is_neutrino(track) && track.id == 0;
453  }
454  }
455 
456  /**
457  * Test whether given track corresponds to a final state particle.
458  *
459  * \param track track
460  * \return true if track corresponds to a final state particle; else false
461  */
462  inline bool is_finalstate(const Trk& track)
463  {
464  if (Evt::ROOT_IO_VERSION >= 16) {
465  return track.status == TRK_ST_FINALSTATE;
466  } else {
467  return !is_initialstate(track);
468  }
469  }
470 
471  /**
472  * Test whether given event has an incoming neutrino.
473  *
474  * \param evt event
475  * \return true if neutrino; else false
476  */
477  inline bool has_neutrino(const Evt& evt)
478  {
479  if (Evt::ROOT_IO_VERSION >= 14) {
480 
481  std::vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
483  return i != evt.mc_trks.cend();
484 
485  } else {
486 
487  return !evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0]);
488  }
489  }
490 
491  /**
492  * Get incoming neutrino.
493  *
494  * \param evt event
495  * \return neutrino
496  */
497  inline const Trk& get_neutrino(const Evt& evt)
498  {
499  if (Evt::ROOT_IO_VERSION >= 14) {
500 
501  std::vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
503 
504  if (i != evt.mc_trks.cend()) { return *i; }
505 
506  } else if (!evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0])) {
507 
508  return evt.mc_trks[0];
509  }
510 
511  THROW(JIndexOutOfRange, "JAANET::get_neutrino(): Event " << evt.id << " has no neutrino.");
512  }
513 
514  /**
515  * Count the number of electrons in a given event
516  *
517  * \param evt event
518  * \return number of electrons
519  */
520  inline int count_electrons(const Evt& evt)
521  {
522  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
523  }
524 
525  /**
526  * Test whether given event has an electron.
527  *
528  * \param evt event
529  * \return true if event has electron; else false
530  */
531  inline bool has_electron(const Evt& evt)
532  {
533  return count_electrons(evt) != 0;
534  }
535 
536  /**
537  * Get first electron from the event tracklist
538  *
539  * \param evt event
540  * \return electron
541  */
542  inline const Trk& get_electron(const Evt& evt)
543  {
544  if (count_electrons(evt) > 0)
545  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
546  else
547  THROW(JIndexOutOfRange, "This event has no electron.");
548  }
549 
550  /**
551  * Test whether given hit was produced by an electron
552  *
553  * Warning: This method will only work with the output of KM3Sim.
554  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
555  *
556  * \param hit hit
557  * \return true if hit was produced by an electron; else false
558  */
559  inline bool from_electron(const Hit& hit)
560  {
561  return hit.type == GEANT4_TYPE_ELECTRON || hit.type == GEANT4_TYPE_ANTIELECTRON;
562  }
563 
564  /**
565  * Count the number of muons in a given event
566  *
567  * \param evt event
568  * \return number of muons
569  */
570  inline int count_muons(const Evt& evt)
571  {
572  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
573  }
574 
575  /**
576  * Test whether given event has a muon.
577  *
578  * \param evt event
579  * \return true if event has muon; else false
580  */
581  inline bool has_muon(const Evt& evt)
582  {
583  return count_muons(evt) != 0;
584  }
585 
586  /**
587  * Get first muon from the event tracklist
588  *
589  * \param evt event
590  * \return muon
591  */
592  inline const Trk& get_muon(const Evt& evt)
593  {
594  if (count_muons(evt) > 0)
595  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
596  else
597  THROW(JIndexOutOfRange, "This event has no muon.");
598  }
599 
600  /**
601  * Test whether given hit was produced by a muon
602  *
603  * Warning: This method will only work with the output of KM3Sim.
604  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
605  *
606  * \param hit hit
607  * \return true if hit was produced by a muon; else false
608  */
609  inline bool from_muon(const Hit& hit)
610  {
611  return hit.type == GEANT4_TYPE_MUON || hit.type == GEANT4_TYPE_ANTIMUON;
612  }
613 
614  /**
615  * Count the number of taus in a given event
616  *
617  * \param evt event
618  * \return number of taus
619  */
620  inline int count_taus(const Evt& evt)
621  {
622  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
623  }
624 
625  /**
626  * Test whether given event has a tau.
627  *
628  * \param evt event
629  * \return true if event has tau; else false
630  */
631  inline bool has_tau(const Evt& evt)
632  {
633  return count_taus(evt) != 0;
634  }
635 
636  /**
637  * Get first tau from the event tracklist
638  *
639  * \param evt event
640  * \return tau
641  */
642  inline const Trk& get_tau(const Evt& evt)
643  {
644  if (count_taus(evt) > 0)
645  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
646  else
647  THROW(JIndexOutOfRange, "This event has no tau.");
648  }
649 
650  /**
651  * Test whether given hit was produced by a tau
652  *
653  * Warning: This method will only work with the output of KM3Sim.
654  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
655  *
656  * \param hit hit
657  * \return true if hit was produced by a tau; else false
658  */
659  inline bool from_tau(const Hit& hit)
660  {
661  return hit.type == GEANT4_TYPE_TAU || hit.type == GEANT4_TYPE_ANTITAU;
662  }
663 
664  /**
665  * Count the number of hadrons in a given event (not including neutral pions)
666  *
667  * \param evt event
668  * \return number of hadrons
669  */
670  inline int count_hadrons(const Evt& evt)
671  {
672  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_hadron);
673  }
674 
675  /**
676  * Test whether given hit was produced by a hadronic shower
677  *
678  * Warning: This method will only work with the output of KM3Sim.
679  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
680  *
681  * \param hit hit
682  * \return true if hit was produced by a hadron; else false
683  */
684  inline bool from_hadron(const Hit& hit)
685  {
686  return (!from_electron(hit) && !from_muon(hit) && !from_tau(hit));
687  }
688 
689 
690  /**
691  * Add time offset to mc event;
692  * according to current definition, the absolute time of the event is defined by the track "t" attribute;
693  * this could change in the future if the global attribute mc_t is assigned to this purpose.
694  *
695  * \param evt event
696  * \param tOff time offset [ns]
697  */
698  inline void add_time_offset(Evt& evt, const double tOff)
699  {
700  for (std::vector<Trk>::iterator p = evt.mc_trks.begin(); p != evt.mc_trks.end(); p++) {
701  p->t += tOff;
702  }
703  }
704 
705 
706  /**
707  * Get kinetic energy of particle with given mass.
708  *
709  * \param E energy [GeV]
710  * \param m mass [GeV]
711  * \return energy [GeV]
712  */
713  inline double getKineticEnergy(const double E, const double m)
714  {
715  if (E > m) {
716  return sqrt((E - m) * (E + m));
717  } else {
718  return 0.0;
719  }
720  }
721 
722 
723  /**
724  * Get track kinetic energy.
725  *
726  * \param trk track
727  * \return kinetic energy [GeV]
728  */
729  inline double getKineticEnergy(const Trk& trk)
730  {
731  const double energy = trk.E;
732  const double mass = JPDB::getInstance().getPDG(trk.type).mass;
733 
734  return getKineticEnergy(energy, mass);
735  }
736 
737 
738  /**
739  * Get initial state energy of a neutrino interaction.\n
740  * This method includes baryon number conservation.
741  *
742  * \param evt event
743  * \return energy [GeV]
744  */
745  inline double getE0(const Evt& evt)
746  {
747  const Trk& neutrino = get_neutrino(evt);
748 
749  double E0 = neutrino.E;
750 
751  for (std::vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.end(); ++track) {
752 
753  if (!is_finalstate(*track)) { continue; }
754 
755  if (track->type == TRACK_TYPE_NEUTRON ||
756  track->type == TRACK_TYPE_SIGMA_MINUS ||
757  track->type == TRACK_TYPE_NEUTRAL_SIGMA) {
758 
759  E0 += MASS_NEUTRON;
760 
761  } else if (track->type == TRACK_TYPE_PROTON ||
762  track->type == TRACK_TYPE_LAMBDA ||
763  track->type == TRACK_TYPE_SIGMA_PLUS) {
764 
765  E0 += MASS_PROTON;
766 
767  } else if (track->type == TRACK_TYPE_ANTINEUTRON) {
768 
769  E0 -= MASS_NEUTRON;
770 
771  } else if (track->type == TRACK_TYPE_ANTIPROTON ||
772  track->type == TRACK_TYPE_ANTILAMBDA) {
773 
774  E0 -= MASS_PROTON;
775  }
776  }
777 
778  return E0;
779  }
780 
781 
782  /**
783  * Get final state energy of a neutrino interaction.\n
784  * This method includes muon energy loss.
785  *
786  * \param evt event
787  * \return energy [GeV]
788  */
789  inline double getE1(const Evt& evt)
790  {
791  double E1 = 0.0;
792 
793  const Trk& neutrino = get_neutrino(evt);
794 
795  for (std::vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
796 
797  if (!is_finalstate(*track)) { continue; }
798 
799  if (is_muon(*track)) {
800 
801  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
802  evt.mc_trks[track->mother_id] : neutrino );
803 
804  const double distance = (track->pos - mother.pos).len();
805 
806  E1 += JPHYSICS::gWater(track->E, -distance);
807 
808  } else {
809 
810  E1 += track->E;
811  }
812  }
813 
814  return E1;
815  }
816 
817 
818  /**
819  * Get momentum vector of the initial state of a neutrino interaction.\n
820  * This method assumes neutrino DIS on a stationary nucleus
821  *
822  * \param evt event
823  * \return energy [GeV]
824  */
825  inline Vec getP0(const Evt& evt)
826  {
827  const Trk& neutrino = get_neutrino(evt);
828 
829  return neutrino.E * neutrino.dir;
830  }
831 
832 
833  /**
834  * Get momentum vector of the final state of a neutrino interaction.\n
835  * This method includes muon energy losses.
836  *
837  * \param evt event
838  * \return final state momentum vector [GeV]
839  */
840  inline Vec getP1(const Evt& evt)
841  {
842  Vec P1(0,0,0);
843 
844  const Trk& neutrino = get_neutrino(evt);
845 
846  for (std::vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
847 
848  if (!is_finalstate(*track)) { continue; }
849 
850  double kineticEnergy = 0.0;
851 
852  if (is_muon(*track)) {
853 
854  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
855  evt.mc_trks[track->mother_id] : neutrino );
856 
857  const double distance = (track->pos - mother.pos).len();
858  const double energy = JPHYSICS::gWater(track->E, -distance);
859 
860  kineticEnergy = getKineticEnergy(energy, JPHYSICS::MASS_MUON);
861 
862  } else {
863 
864  kineticEnergy = getKineticEnergy(*track);
865  }
866 
867  P1 += kineticEnergy * track->dir;
868  }
869 
870  return P1;
871  }
872 }
873 
874 #endif
bool from_tau(const Hit &hit)
Test whether given hit was produced by a tau.
Vec getP1(const Evt &evt)
Get momentum vector of the final state of a neutrino interaction.
JVertex3D getVertex(const Trk &track)
Get vertex.
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.
JPredicate< JResult_t T::*, JComparison::eq > make_predicate(JResult_t T::*member, const JResult_t value)
Helper method to create predicate for data member.
Definition: JPredicate.hh:128
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
then usage $script< input file >[option[primary[working directory]]] nWhere option can be E
Definition: JMuonPostfit.sh:40
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.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
double getE0(const Evt &evt)
Get initial state energy of a neutrino interaction.
static const int TRK_MOTHER_NONE
mother id of a particle if it has no parent
Definition: trkmembers.hh:13
Energy loss of muon.
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:712
double getE1(const Evt &evt)
Get final state energy of a neutrino interaction.
static const double MASS_MUON
muon mass [GeV]
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.
static const JGeaneWater gWater
Function object for energy loss of muon in sea water.
Definition: JGeane.hh:381
const JParticle & getPDG(const int pdg) const
Get particle corresponding to given PDG code.
Definition: JPDB.hh:234
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.
static const JPDB & getInstance()
Get particle data book.
Definition: JPDB.hh:125
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.
static const double MASS_NEUTRON
neutron mass [GeV]
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.
bool is_finalstate(const Trk &track)
Test whether given track corresponds to a final state particle.
int mother_id
MC id of the parent particle.
Definition: Trk.hh:29
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
double getKineticEnergy(const double E, const double m)
Get kinetic energy of particle with given mass.
JDirection3D getDirection(const Vec &dir)
Get direction.
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, see km3net-dataformat/definitions/fitparameters.csv
Definition: Trk.hh:32
Physics constants.
static const int TRK_ST_PRIMARYNEUTRINO
initial state neutrino (&#39;neutrino&#39; tag in evt files from gseagen and genhen).
Definition: trkmembers.hh:16
static const int TRK_ST_MUONBUNDLE
initial state muon bundle (mupage)
Definition: trkmembers.hh:18
static const int TRK_ST_ININUCLEI
Initial state nuclei (gseagen)
Definition: trkmembers.hh:19
static const int TRK_ST_PRIMARYCOSMIC
initial state cosmic ray (&#39;track_primary&#39; tag in evt files from corant).
Definition: trkmembers.hh:17
double getY() const
Get y position.
Definition: JVector3D.hh:104
bool is_initialstate(const Trk &track)
Test whether given track corresponds to an initial state particle.
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.
static const int TRK_MOTHER_UNDEFINED
KM3NeT Data Definitions v3.1.0 https://git.km3net.de/common/km3net-dataformat.
Definition: trkmembers.hh:12
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
int id
track identifier
Definition: Trk.hh:16
int status
MC status code, see km3net-dataformat/definitions/trkmembers.csv for values.
Definition: Trk.hh:28
then for APP in event gandalf start energy
Definition: JMuonMCEvt.sh:44
static const int TRK_ST_FINALSTATE
for MC: the particle must be processed by detector simulation (&#39;track_in&#39; tag in evt files)...
Definition: trkmembers.hh:15
bool is_pion(const Trk &track)
Test whether given track is a charged pion.
Vec pos
postion [m] of the track at time t
Definition: Trk.hh:17
double mass
mass of particle [GeV]
Definition: JPDB.hh:92
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.
then P1
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:48
double t
hit time (from tdc+calibration or MC truth)
Definition: Hit.hh:23
Scattered light from delta-rays.
static const double MASS_PROTON
proton mass [GeV]
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:22
double getNPE(const Hit &hit)
Get true charge of hit.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
static int ROOT_IO_VERSION
Streamer version as obtained from ROOT file.
Definition: Evt.hh:274
bool is_neutron(const Trk &track)
Test whether given track is a (anti-)neutron.
Exception for accessing a value in a collection that is outside of its range.
Definition: JException.hh:178
static JRange< double, std::less< double > > DEFAULT_RANGE()
Default range.
Definition: JRange.hh:555
Exception for accessing an index in a collection that is outside of its range.
Definition: JException.hh:106
Vec getP0(const Evt &evt)
Get momentum vector of the initial state of a neutrino interaction.
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:49
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:20
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [s]).
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.