Jpp  master_rocky-37-gf0c5bc59d
the software that should make you happy
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"
29 #include "JPhysics/JTimeRange.hh"
30 #include "JPhysics/JGeane.hh"
31 
32 #include "JAAnet/JParticleTypes.hh"
33 #include "JAAnet/JPDB.hh"
34 
35 
36 /**
37  * \file
38  *
39  * Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
40  * \author mdejong
41  */
42 namespace JAANET {}
43 namespace JPP { using namespace JAANET; }
44 
45 /**
46  * Extensions to Evt data format.
47  */
48 namespace JAANET {
49 
53 
61 
66 
67  /**
68  * AAnet shower fit reconstruction type.
69  */
70  static const int AASHOWER_RECONSTRUCTION_TYPE = 101;
71 
72 
73  /**
74  * Enumeration of hit types based on km3 codes.
75  */
76  enum JHitType_t {
77  HIT_TYPE_MUON_DIRECT = +5, //!< Direct light from muon
78  HIT_TYPE_MUON_SCATTERED = -5, //!< Scattered light from muon
79  HIT_TYPE_DELTARAYS_DIRECT = +4, //!< Direct light from delta-rays
80  HIT_TYPE_DELTARAYS_SCATTERED = -4, //!< Scattered light from delta-rays
81  HIT_TYPE_BREMSSTRAHLUNG_DIRECT = +3, //!< Direct light from Bremsstrahlung
82  HIT_TYPE_BREMSSTRAHLUNG_SCATTERED = -3, //!< Scattered light from Bremsstrahlung
83  HIT_TYPE_SHOWER_DIRECT = +99, //!< Direct light from primary shower
84  HIT_TYPE_SHOWER_SCATTERED = -99, //!< Scattered light from primary shower
85  HIT_TYPE_NOISE = -1, //!< Random noise
86  HIT_TYPE_UNKNOWN = 0 //!< Unknown source
87  };
88 
89 
90  /**
91  * Get true time of hit.
92  *
93  * \param hit hit
94  * \return true time of hit
95  */
96  inline double getTime(const Hit& hit)
97  {
98  return hit.t;
99  }
100 
101 
102  /**
103  * Get true charge of hit.
104  *
105  * \param hit hit
106  * \return true charge of hit
107  */
108  inline double getNPE(const Hit& hit)
109  {
110  return hit.a;
111  }
112 
113 
114  /**
115  * Verify hit origin.
116  *
117  * \param hit hit
118  * \return true if noise; else false
119  */
120  inline bool is_noise(const Hit& hit)
121  {
122  return hit.type == HIT_TYPE_NOISE;
123  }
124 
125 
126  /**
127  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
128  * Note that the global event time in not taken into account.
129  *
130  * \param event event
131  * \return time range
132  */
133  inline JTimeRange getTimeRange(const Evt& event)
134  {
136 
137  for (std::vector<Hit>::const_iterator hit = event.mc_hits.begin(); hit != event.mc_hits.end(); ++hit) {
138  if (!is_noise(*hit)) {
139  time_range.include(getTime(*hit));
140  }
141  }
142 
143  return time_range;
144  }
145 
146 
147  /**
148  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
149  * The resulting time range is larger than or equal to the given time window.\n
150  * Note that the global event time in not taken into account.
151  *
152  * \param event event
153  * \param T_ns time window [ns]
154  * \return time range
155  */
156  inline JTimeRange getTimeRange(const Evt& event, const JTimeRange& T_ns)
157  {
158  JTimeRange time_range = getTimeRange(event);
159 
160  if (!time_range.is_valid()) {
161  time_range = T_ns;
162  }
163 
164  if (time_range.getLength() < T_ns.getLength()) {
165  time_range.setLength(T_ns.getLength());
166  }
167 
168  return time_range;
169  }
170 
171 
172  /**
173  * Get position.
174  *
175  * \param pos position
176  * \return position
177  */
178  inline JPosition3D getPosition(const Vec& pos)
179  {
180  return JPosition3D(pos.x, pos.y, pos.z);
181  }
182 
183 
184  /**
185  * Get position.
186  *
187  * \param pos position
188  * \return position
189  */
190  inline Vec getPosition(const JPosition3D& pos)
191  {
192  return Vec(pos.getX(), pos.getY(), pos.getZ());
193  }
194 
195 
196  /**
197  * Get position.
198  *
199  * \param track track
200  * \return position
201  */
202  inline JPosition3D getPosition(const Trk& track)
203  {
204  return getPosition(track.pos);
205  }
206 
207 
208  /**
209  * Get direction.
210  *
211  * \param dir direction
212  * \return direction
213  */
214  inline JDirection3D getDirection(const Vec& dir)
215  {
216  return JDirection3D(dir.x, dir.y, dir.z);
217  }
218 
219 
220  /**
221  * Get direction.
222  *
223  * \param dir direction
224  * \return direction
225  */
226  inline Vec getDirection(const JDirection3D& dir)
227  {
228  return Vec(dir.getDX(), dir.getDY(), dir.getDZ());
229  }
230 
231 
232  /**
233  * Get direction.
234  *
235  * \param track track
236  * \return direction
237  */
238  inline JDirection3D getDirection(const Trk& track)
239  {
240  return getDirection(track.dir);
241  }
242 
243 
244  /**
245  * Get axis.
246  *
247  * \param track track
248  * \return axis
249  */
250  inline JAxis3D getAxis(const Trk& track)
251  {
252  return JAxis3D(getPosition(track), getDirection(track));
253  }
254 
255 
256  /**
257  * Get track.
258  *
259  * \param track track
260  * \return track
261  */
262  inline JTrack3E getTrack(const Trk& track)
263  {
264  return JTrack3E(JTrack3D(getAxis(track), track.t), track.E);
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 (anti-)muon.
308  *
309  * \param track track
310  * \return true if (anti-)muon; else false
311  */
312  inline bool is_muonbundle(const Trk& track) { return (abs(track.type) == PDG_MUONBUNDLE); }
313 
314  /**
315  * Test whether given track is a neutrino.
316  *
317  * \param track track
318  * \return true if neutrino; else false
319  */
320  inline bool is_neutrino(const Trk& track) { return (abs(track.type) == TRACK_TYPE_NUE ||
321  abs(track.type) == TRACK_TYPE_NUMU ||
322  abs(track.type) == TRACK_TYPE_NUTAU); }
323 
324  /**
325  * Test whether given track is a (anti-)electron.
326  *
327  * \param track track
328  * \return true if (anti-)electron; else false
329  */
330  inline bool is_electron(const Trk& track) { return (abs(track.type) == TRACK_TYPE_ELECTRON); }
331 
332  /**
333  * Test whether given track is a (anti-)muon.
334  *
335  * \param track track
336  * \return true if (anti-)muon; else false
337  */
338  inline bool is_muon (const Trk& track) { return (abs(track.type) == TRACK_TYPE_MUON); }
339 
340  /**
341  * Test whether given track is a (anti-)tau.
342  *
343  * \param track track
344  * \return true if (anti-)tau; else false
345  */
346  inline bool is_tau (const Trk& track) { return (abs(track.type) == TRACK_TYPE_TAU); }
347 
348  /**
349  * Test whether given track is a charged pion.
350  *
351  * \param track track
352  * \return true if charged pion; else false
353  */
354  inline bool is_pion (const Trk& track) { return (abs(track.type) == TRACK_TYPE_CHARGED_PION_PLUS); }
355 
356  /**
357  * Test whether given track is a (anti-)proton.
358  *
359  * \param track track
360  * \return true if (anti-)proton; else false
361  */
362  inline bool is_proton (const Trk& track) { return (abs(track.type) == TRACK_TYPE_PROTON); }
363 
364  /**
365  * Test whether given track is a (anti-)neutron.
366  *
367  * \param track track
368  * \return true if (anti-)neutron; else false
369  */
370  inline bool is_neutron (const Trk& track) { return (abs(track.type) == TRACK_TYPE_NEUTRON); }
371 
372  /**
373  * Test whether given track is a lepton
374  *
375  * \param track track
376  * \return true if lepton; else fails
377  */
378  inline bool is_lepton (const Trk& track) { return (is_neutrino(track) ||
379  is_electron(track) ||
380  is_muon (track) ||
381  is_tau (track)); }
382 
383  /**
384  * Test whether given track is a charged lepton
385  *
386  * \param track track
387  * \return true if lepton; else fails
388  */
389  inline bool is_charged_lepton (const Trk& track) { return (is_electron(track) ||
390  is_muon (track) ||
391  is_tau (track)); }
392 
393  /**
394  * Test whether given track is a hadron
395  *
396  * \param track track
397  * \return true if hadron; else fails
398  */
399  inline bool is_hadron (const Trk& track) { return (abs(track.type) != TRACK_TYPE_PHOTON &&
400  !is_lepton(track)); }
401 
402  /**
403  * Test whether given track is a meson (is hadron and third digit of PDG code is zero)
404  *
405  * \param track track
406  * \return true if hadron; else fails
407  */
408  inline bool is_meson (const Trk& track) { return (is_hadron(track) &&
409  ((int)(track.type / 1000)) % 10 == 0); }
410 
411  /**
412  * Test whether given track is a baryon (is hadron and third digit of PDG code is not zero)
413  *
414  * \param track track
415  * \return true if hadron; else fails
416  */
417  inline bool is_baryon (const Trk& track) { return (is_hadron(track) &&
418  ((int)(track.type / 1000)) % 10 != 0); }
419 
420  /**
421  * Test whether given track is a nucleus (PDG code corresponds to a baryon or has 10 digits, starting with '10')
422  *
423  * \param track track
424  * \return true if hadron; else fails
425  */
426  inline bool is_nucleus (const Trk& track)
427  {
428  using namespace std;
429 
430  const string typestr = to_string(abs(track.type));
431 
432  return (is_baryon(track) || (typestr.size() == 10 && typestr.substr(0,2) == "10"));
433  }
434 
435  /**
436  * Test whether given track corresponds to given particle type.
437  *
438  * \param track track
439  * \param type particle type
440  * \return true if matches the corresponding particle; else false
441  */
442  inline bool has_particleID (const Trk& track, const int type) { return (track.type == type); }
443 
444  /**
445  * Test whether given track corresponds to an initial state particle.
446  *
447  * \param track track
448  * \return true if track corresponds to an initial state particle; else false
449  */
450  inline bool is_initialstate(const Trk& track)
451  {
452  if (Evt::ROOT_IO_VERSION >= 14) {
453 
454  return (track.status == TRK_ST_PRIMARYNEUTRINO ||
455  track.status == TRK_ST_PRIMARYCOSMIC ||
456  track.status == TRK_ST_MUONBUNDLE ||
457  track.status == TRK_ST_ININUCLEI);
458 
459  } else if (Evt::ROOT_IO_VERSION >= 11) {
460 
461  return (track.mother_id == TRK_MOTHER_UNDEFINED ||
462  track.mother_id == TRK_MOTHER_NONE);
463 
464  } else {
465 
466  return is_neutrino(track) && track.id == 0;
467  }
468  }
469 
470  /**
471  * Test whether given track corresponds to a final state particle.
472  *
473  * \param track track
474  * \return true if track corresponds to a final state particle; else false
475  */
476  inline bool is_finalstate(const Trk& track)
477  {
478  if (Evt::ROOT_IO_VERSION >= 15) {
479  return track.status == TRK_ST_FINALSTATE;
480  } else {
481  return !is_initialstate(track);
482  }
483  }
484 
485  /**
486  * Test whether given event has an incoming neutrino.
487  *
488  * \param evt event
489  * \return true if neutrino; else false
490  */
491  inline bool has_neutrino(const Evt& evt)
492  {
493  using namespace std;
494 
495  if (Evt::ROOT_IO_VERSION >= 14) {
496 
497  vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
499  return i != evt.mc_trks.cend();
500 
501  } else {
502 
503  return !evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0]);
504  }
505  }
506 
507  /**
508  * Get incoming neutrino.
509  *
510  * \param evt event
511  * \return neutrino
512  */
513  inline const Trk& get_neutrino(const Evt& evt)
514  {
515  using namespace std;
516 
517  if (Evt::ROOT_IO_VERSION >= 14) {
518 
519  vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
521 
522  if (i != evt.mc_trks.cend()) { return *i; }
523 
524  } else if (!evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0])) {
525 
526  return evt.mc_trks[0];
527  }
528 
529  THROW(JIndexOutOfRange, "JAANET::get_neutrino(): Event " << evt.id << " has no neutrino.");
530  }
531 
532  /**
533  * Test whether given event has well-defined target nucleus.
534  *
535  * \param evt event
536  * \return true if target nucleus is defined; else false
537  */
538  inline bool has_target_nucleus(const Evt& evt)
539  {
540  using namespace std;
541 
542  if (Evt::ROOT_IO_VERSION >= 14) {
543 
544  vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
546 
547  return i != evt.mc_trks.cend();
548 
549  } else {
550 
551  return false;
552  }
553  }
554 
555  /**
556  * Get target nucleus.
557  *
558  * \param evt event
559  * \return target nucleus
560  */
561  inline const Trk& get_target_nucleus(const Evt& evt)
562  {
563  using namespace std;
564  using namespace JPP;
565 
566  vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
568 
569  if (i != evt.mc_trks.cend()) {
570  return *i;
571  } else {
572  THROW(JIndexOutOfRange, "JAANET::get_target_nucleus(): Cannot retrieve target nucleus for event " << evt.id << ".");
573  }
574  }
575 
576 
577  /**
578  * Get primary.
579  *
580  * \param evt event
581  * \return primary track
582  */
583  inline const Trk& get_primary(const Evt& evt)
584  {
585  using namespace std;
586  using namespace JPP;
587 
588  for (vector<Trk>::const_iterator i = evt.mc_trks.cbegin(); i != evt.mc_trks.cend(); ++i) {
589  if (is_initialstate(*i) && i->status != TRK_ST_ININUCLEI) {
590  return *i;
591  }
592  }
593 
594  THROW(JIndexOutOfRange, "JAANET::get_primary(): Cannot retrieve primary track for event " << evt.id << ".");
595  }
596 
597 
598  /**
599  * Auxiliary function to retrieve the leading lepton of a neutrino interaction.
600  *
601  * \param event event
602  * \return leading lepton
603  */
604  inline const Trk& get_leading_lepton(const Evt& event)
605  {
606  using namespace std;
607  using namespace JPP;
608 
609  const Trk& primary = get_primary(event);
610 
611  if (is_neutrino(primary)) {
612 
613  for (vector<Trk>::const_iterator i = event.mc_trks.cbegin(); i != event.mc_trks.cend(); ++i) {
614  if (is_lepton(*i) && i->mother_id == primary.id) {
615  return *i;
616  }
617  }
618 
619  THROW(JValueOutOfRange, "get_leading_lepton(): Could not find leading lepton for given event");
620 
621  } else {
622  THROW(JValueOutOfRange, "get_leading_lepton(): Given event does not correspond to a neutrino interaction.");
623  }
624  }
625 
626 
627  /**
628  * Get vertex.
629  *
630  * \param track track
631  * \return vertex
632  */
633  inline JVertex3D getVertex(const Trk& track)
634  {
635  return JVertex3D(getPosition(track), track.t);
636  }
637 
638  /**
639  * Get event vertex.
640  *
641  * \param event event
642  * \return event vertex
643  */
644  inline JVertex3D getVertex(const Evt& event)
645  {
646  using namespace std;
647  using namespace JPP;
648 
649  if (has_neutrino(event)) {
650 
651  const Trk& neutrino = get_neutrino(event);
652 
653  return getVertex(neutrino);
654 
655  } else if (!event.mc_trks.empty()) {
656 
657  vector<Trk>::const_iterator i = find_if(event.mc_trks.begin(), event.mc_trks.end(), is_initialstate);
658 
659  if (i != event.mc_trks.end()) {
660  return getVertex(*i);
661  } else {
662  THROW(JValueOutOfRange, "getVertex(): No initial state track found.");
663  }
664 
665  } else if (!event.trks.empty()) { // For reconstructed DAQ events
666 
667  const Trk& track = get_best_reconstructed_track<JPP_RECONSTRUCTION_TYPE>(event);
668 
669  return getVertex(track);
670 
671  } else {
672 
673  THROW(JValueOutOfRange, "getVertex(): Could not find a valid event vertex.");
674  }
675  }
676 
677  /**
678  * Count the number of electrons in a given event
679  *
680  * \param evt event
681  * \return number of electrons
682  */
683  inline int count_electrons(const Evt& evt)
684  {
685  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
686  }
687 
688  /**
689  * Test whether given event has an electron.
690  *
691  * \param evt event
692  * \return true if event has electron; else false
693  */
694  inline bool has_electron(const Evt& evt)
695  {
696  return count_electrons(evt) != 0;
697  }
698 
699  /**
700  * Get first electron from the event tracklist
701  *
702  * \param evt event
703  * \return electron
704  */
705  inline const Trk& get_electron(const Evt& evt)
706  {
707  if (count_electrons(evt) > 0)
708  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
709  else
710  THROW(JIndexOutOfRange, "This event has no electron.");
711  }
712 
713  /**
714  * Test whether given hit was produced by an electron
715  *
716  * Warning: This method will only work with the output of KM3Sim.
717  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
718  *
719  * \param hit hit
720  * \return true if hit was produced by an electron; else false
721  */
722  inline bool from_electron(const Hit& hit)
723  {
724  return hit.type == GEANT4_TYPE_ELECTRON || hit.type == GEANT4_TYPE_ANTIELECTRON;
725  }
726 
727  /**
728  * Count the number of muons in a given event
729  *
730  * \param evt event
731  * \return number of muons
732  */
733  inline int count_muons(const Evt& evt)
734  {
735  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
736  }
737 
738  /**
739  * Test whether given event has a muon.
740  *
741  * \param evt event
742  * \return true if event has muon; else false
743  */
744  inline bool has_muon(const Evt& evt)
745  {
746  return count_muons(evt) != 0;
747  }
748 
749  /**
750  * Get first muon from the event tracklist
751  *
752  * \param evt event
753  * \return muon
754  */
755  inline const Trk& get_muon(const Evt& evt)
756  {
757  if (count_muons(evt) > 0)
758  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
759  else
760  THROW(JIndexOutOfRange, "This event has no muon.");
761  }
762 
763  /**
764  * Test whether given hit was produced by a muon
765  *
766  * Warning: This method will only work with the output of KM3Sim.
767  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
768  *
769  * \param hit hit
770  * \return true if hit was produced by a muon; else false
771  */
772  inline bool from_muon(const Hit& hit)
773  {
774  return hit.type == GEANT4_TYPE_MUON || hit.type == GEANT4_TYPE_ANTIMUON;
775  }
776 
777  /**
778  * Count the number of taus in a given event
779  *
780  * \param evt event
781  * \return number of taus
782  */
783  inline int count_taus(const Evt& evt)
784  {
785  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
786  }
787 
788  /**
789  * Test whether given event has a tau.
790  *
791  * \param evt event
792  * \return true if event has tau; else false
793  */
794  inline bool has_tau(const Evt& evt)
795  {
796  return count_taus(evt) != 0;
797  }
798 
799  /**
800  * Get first tau from the event tracklist
801  *
802  * \param evt event
803  * \return tau
804  */
805  inline const Trk& get_tau(const Evt& evt)
806  {
807  if (count_taus(evt) > 0)
808  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
809  else
810  THROW(JIndexOutOfRange, "This event has no tau.");
811  }
812 
813 
814  /**
815  * Test whether given event contains a tau-lepton decay into a muon.
816  *
817  * \param event event
818  * \return true if event contains a muonic tau-lepton decay; else false
819  */
820  inline bool has_muonic_taudecay(const Evt& event)
821  {
822  using namespace std;
823 
824  if (has_tau(event)) {
825 
826  const Trk& tau = get_tau(event);
827 
828  for (vector<Trk>::const_iterator i = event.mc_trks.cbegin(); i != event.mc_trks.cend(); ++i) {
829  if (is_muon(*i) && i->mother_id == tau.id) {
830  return true;
831  }
832  }
833  }
834 
835  return false;
836  }
837 
838  /**
839  * Test whether given event contains a leptonic tau-decay.
840  *
841  * \param event event
842  * \return true if event contains a leptonic tau-decay; else false
843  */
844  inline bool has_leptonic_taudecay(const Evt& event)
845  {
846  using namespace std;
847 
848  if (has_tau(event)) {
849 
850  const Trk& tau = get_tau(event);
851 
852  for (vector<Trk>::const_iterator i = event.mc_trks.cbegin(); i != event.mc_trks.cend(); ++i) {
853  if ((is_electron(*i) || is_muon(*i)) && i->mother_id == tau.id) {
854  return true;
855  }
856  }
857  }
858 
859  return false;
860  }
861 
862  /**
863  * Test whether given event contains a hadronic tau-decay.
864  *
865  * \param event event
866  * \return true if event contains a hadronic tau-decay; else false
867  */
868  inline bool has_hadronic_taudecay(const Evt& event)
869  {
870  return has_tau(event) && !has_leptonic_taudecay(event);
871  }
872 
873  /**
874  * Count the number of tau-lepton decay prongs in a given event.\n
875  * The number of prongs is defined as the number of charged tau-lepton decay products.
876  *
877  * \param event event
878  * \return number of tau-lepton decay prongs
879  */
880  inline int count_taudecay_prongs(const Evt& event)
881  {
882  using namespace std;
883 
884  int n = 0;
885 
886  if (has_tau(event)) {
887 
888  const Trk& tau = get_tau(event);
889 
890  for (vector<Trk>::const_iterator i = event.mc_trks.cbegin(); i != event.mc_trks.cend(); ++i) {
891  if (i->mother_id == tau.id && JPDB::getInstance().getPDG(i->type).charge != 0) {
892  ++n;
893  }
894  }
895  }
896 
897  return n;
898  }
899 
900  /**
901  * Test whether given hit was produced by a tau
902  *
903  * Warning: This method will only work with the output of KM3Sim.
904  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
905  *
906  * \param hit hit
907  * \return true if hit was produced by a tau; else false
908  */
909  inline bool from_tau(const Hit& hit)
910  {
911  return hit.type == GEANT4_TYPE_TAU || hit.type == GEANT4_TYPE_ANTITAU;
912  }
913 
914  /**
915  * Count the number of hadrons in a given event (not including neutral pions)
916  *
917  * \param evt event
918  * \return number of hadrons
919  */
920  inline int count_hadrons(const Evt& evt)
921  {
922  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_hadron);
923  }
924 
925  /**
926  * Test whether given hit was produced by a hadronic shower
927  *
928  * Warning: This method will only work with the output of KM3Sim.
929  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
930  *
931  * \param hit hit
932  * \return true if hit was produced by a hadron; else false
933  */
934  inline bool from_hadron(const Hit& hit)
935  {
936  return (!from_electron(hit) && !from_muon(hit) && !from_tau(hit));
937  }
938 
939 
940  /**
941  * Add time offset to mc event;
942  * according to current definition, the absolute time of the event is defined by the track "t" attribute;
943  * this could change in the future if the global attribute mc_t is assigned to this purpose.
944  *
945  * \param evt event
946  * \param tOff time offset [ns]
947  */
948  inline void add_time_offset(Evt& evt, const double tOff)
949  {
950  for (std::vector<Trk>::iterator p = evt.mc_trks.begin(); p != evt.mc_trks.end(); p++) {
951  p->t += tOff;
952  }
953  }
954 
955 
956  /**
957  * Get track kinetic energy.
958  *
959  * \param trk track
960  * \return kinetic energy [GeV]
961  */
962  inline double getKineticEnergy(const Trk& trk)
963  {
964  const double energy = trk.E;
965  const double mass = JPDB::getInstance().getPDG(trk.type).mass;
966 
967  return getKineticEnergy(energy, mass);
968  }
969 
970 
971  /**
972  * Get initial state energy of a neutrino interaction.\n
973  * This method includes baryon number conservation.
974  *
975  * \param evt event
976  * \return energy [GeV]
977  */
978  inline double getE0(const Evt& evt)
979  {
980  using namespace std;
981 
982  const Trk& neutrino = get_neutrino(evt);
983 
984  double E0 = neutrino.E;
985 
986  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.end(); ++track) {
987 
988  if (!is_finalstate(*track)) { continue; }
989 
990  if (track->type == TRACK_TYPE_NEUTRON ||
991  track->type == TRACK_TYPE_SIGMA_MINUS ||
992  track->type == TRACK_TYPE_NEUTRAL_SIGMA) {
993 
994  E0 += MASS_NEUTRON;
995 
996  } else if (track->type == TRACK_TYPE_PROTON ||
997  track->type == TRACK_TYPE_LAMBDA ||
998  track->type == TRACK_TYPE_SIGMA_PLUS) {
999 
1000  E0 += MASS_PROTON;
1001 
1002  } else if (track->type == TRACK_TYPE_ANTINEUTRON) {
1003 
1004  E0 -= MASS_NEUTRON;
1005 
1006  } else if (track->type == TRACK_TYPE_ANTIPROTON ||
1007  track->type == TRACK_TYPE_ANTILAMBDA) {
1008 
1009  E0 -= MASS_PROTON;
1010  }
1011  }
1012 
1013  return E0;
1014  }
1015 
1016 
1017  /**
1018  * Get final state energy of a neutrino interaction.\n
1019  * This method includes muon energy loss.
1020  *
1021  * \param evt event
1022  * \return energy [GeV]
1023  */
1024  inline double getE1(const Evt& evt)
1025  {
1026  using namespace std;
1027  using namespace JPP;
1028 
1029  double E1 = 0.0;
1030 
1031  const Trk& neutrino = get_neutrino(evt);
1032 
1033  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
1034 
1035  if (!is_finalstate(*track)) { continue; }
1036 
1037  if (is_muon(*track)) {
1038 
1039  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
1040  evt.mc_trks[track->mother_id] : neutrino );
1041 
1042  const double distance = (track->pos - mother.pos).len();
1043 
1044  E1 += gWater(track->E, -distance);
1045 
1046  } else {
1047 
1048  E1 += track->E;
1049  }
1050  }
1051 
1052  return E1;
1053  }
1054 
1055 
1056  /**
1057  * Get momentum vector of the initial state of a neutrino interaction.\n
1058  * This method assumes neutrino DIS on a stationary nucleus
1059  *
1060  * \param evt event
1061  * \return energy [GeV]
1062  */
1063  inline Vec getP0(const Evt& evt)
1064  {
1065  const Trk& neutrino = get_neutrino(evt);
1066 
1067  return neutrino.E * neutrino.dir;
1068  }
1069 
1070 
1071  /**
1072  * Get momentum vector of the final state of a neutrino interaction.\n
1073  * This method includes muon energy losses.
1074  *
1075  * \param evt event
1076  * \return final state momentum vector [GeV]
1077  */
1078  inline Vec getP1(const Evt& evt)
1079  {
1080  using namespace std;
1081  using namespace JPP;
1082 
1083  Vec P1(0,0,0);
1084 
1085  const Trk& neutrino = get_neutrino(evt);
1086 
1087  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
1088 
1089  if (!is_finalstate(*track)) { continue; }
1090 
1091  double kineticEnergy = 0.0;
1092 
1093  if (is_muon(*track)) {
1094 
1095  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
1096  evt.mc_trks[track->mother_id] : neutrino );
1097 
1098  const double distance = (track->pos - mother.pos).len();
1099  const double energy = gWater(track->E, -distance);
1100 
1101  kineticEnergy = getKineticEnergy(energy, MASS_MUON);
1102 
1103  } else {
1104 
1105  kineticEnergy = getKineticEnergy(*track);
1106  }
1107 
1108  P1 += kineticEnergy * track->dir;
1109  }
1110 
1111  return P1;
1112  }
1113 
1114 
1115  /**
1116  * Get final state invariant mass.
1117  *
1118  * \param event event
1119  * \return invariant mass [GeV]
1120  */
1121  inline double getInvariantMass(const Evt& event)
1122  {
1123  using namespace std;
1124  using namespace JPP;
1125 
1126  double Minv = 0.0;
1127 
1128  for (vector<Trk>::const_iterator track = event.mc_trks.cbegin(); track != event.mc_trks.cend(); ++track) {
1129 
1130  if (!is_finalstate(*track)) { continue; }
1131 
1132  const double mass = JPDB::getInstance().getPDG(track->type).mass;
1133 
1134  Minv += mass;
1135  }
1136 
1137  return Minv;
1138  }
1139 }
1140 
1141 #endif
Exceptions.
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:712
Energy loss of muon.
Definition of particle types.
Auxiliary methods for physics calculations.
Physics constants.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Axis object.
Definition: JAxis3D.hh:41
Data structure for direction in three dimensions.
Definition: JDirection3D.hh:35
Data structure for position in three dimensions.
Definition: JPosition3D.hh:38
3D track with energy.
Definition: JTrack3E.hh:32
double getY() const
Get y position.
Definition: JVector3D.hh:104
double getZ() const
Get z position.
Definition: JVector3D.hh:115
double getX() const
Get x position.
Definition: JVector3D.hh:94
double getDY() const
Get y direction.
Definition: JVersor3D.hh:106
double getDX() const
Get x direction.
Definition: JVersor3D.hh:95
double getDZ() const
Get z direction.
Definition: JVersor3D.hh:117
Exception for accessing an index in a collection that is outside of its range.
Definition: JException.hh:108
Exception for accessing a value in a collection that is outside of its range.
Definition: JException.hh:180
static JRange< double, std::less< double > > DEFAULT_RANGE()
Default range.
Definition: JRange.hh:555
Extensions to Evt data format.
JAxis3D getAxis(const Trk &track)
Get axis.
bool is_neutron(const Trk &track)
Test whether given track is a (anti-)neutron.
bool is_charged_lepton(const Trk &track)
Test whether given track is a charged lepton.
const Trk & get_primary(const Evt &evt)
Get primary.
bool is_photon(const Trk &track)
Test whether given track is a photon or neutral pion.
const Trk & get_neutrino(const Evt &evt)
Get incoming neutrino.
const Trk & get_muon(const Evt &evt)
Get first muon from the event tracklist.
int count_taudecay_prongs(const Evt &event)
Count the number of tau-lepton decay prongs in a given event.
int count_muons(const Evt &evt)
Count the number of muons in a given event.
bool from_tau(const Hit &hit)
Test whether given hit was produced by a tau.
bool has_target_nucleus(const Evt &evt)
Test whether given event has well-defined target nucleus.
JTrack3E getTrack(const Trk &track)
Get track.
const Trk & get_target_nucleus(const Evt &evt)
Get target nucleus.
bool is_electron(const Trk &track)
Test whether given track is a (anti-)electron.
bool has_tau(const Evt &evt)
Test whether given event has a tau.
bool is_nucleus(const Trk &track)
Test whether given track is a nucleus (PDG code corresponds to a baryon or has 10 digits,...
bool from_muon(const Hit &hit)
Test whether given hit was produced by a 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...
Vec getP0(const Evt &evt)
Get momentum vector of the initial state of a neutrino interaction.
double getInvariantMass(const Evt &event)
Get final state invariant mass.
double getKineticEnergy(const Trk &trk)
Get track kinetic energy.
bool from_hadron(const Hit &hit)
Test whether given hit was produced by a hadronic shower.
int count_taus(const Evt &evt)
Count the number of taus in a given event.
const Trk & get_electron(const Evt &evt)
Get first electron from the event tracklist.
JVertex3D getVertex(const Evt &event)
Get event vertex.
bool is_muonbundle(const Trk &track)
Test whether given track is a (anti-)muon.
JTransformation3D getTransformation(const Trk &track)
Get transformation.
bool has_muon(const Evt &evt)
Test whether given event has a muon.
const Trk & get_leading_lepton(const Evt &event)
Auxiliary function to retrieve the leading lepton of a neutrino interaction.
bool is_pion(const Trk &track)
Test whether given track is a charged pion.
JTimeRange getTimeRange(const Evt &event, const JTimeRange &T_ns)
Get time range (i.e. time between earliest and latest hit) of Monte Carlo event.
static const int AASHOWER_RECONSTRUCTION_TYPE
AAnet shower fit reconstruction type.
bool has_neutrino(const Evt &evt)
Test whether given event has an incoming neutrino.
bool is_finalstate(const Trk &track)
Test whether given track corresponds to a final state particle.
bool has_electron(const Evt &evt)
Test whether given event has an electron.
bool is_neutrino(const Trk &track)
Test whether given track is a neutrino.
JDirection3D getDirection(const Trk &track)
Get direction.
double getE1(const Evt &evt)
Get final state energy of a neutrino interaction.
double getTime(const Hit &hit)
Get true time of hit.
double getE0(const Evt &evt)
Get initial state energy of a neutrino interaction.
bool is_initialstate(const Trk &track)
Test whether given track corresponds to an initial state particle.
double getW(const Trk &track, const size_t index, const double value)
Get track information.
JPosition3D getPosition(const Trk &track)
Get position.
int count_electrons(const Evt &evt)
Count the number of electrons in a given event.
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
JHitType_t
Enumeration of hit types based on km3 codes.
@ HIT_TYPE_UNKNOWN
Unknown source.
@ HIT_TYPE_NOISE
Random noise.
@ HIT_TYPE_DELTARAYS_DIRECT
Direct light from delta-rays.
@ HIT_TYPE_MUON_DIRECT
Direct light from muon.
@ HIT_TYPE_DELTARAYS_SCATTERED
Scattered light from delta-rays.
@ HIT_TYPE_MUON_SCATTERED
Scattered light from muon.
@ HIT_TYPE_SHOWER_DIRECT
Direct light from primary shower.
@ HIT_TYPE_BREMSSTRAHLUNG_SCATTERED
Scattered light from Bremsstrahlung.
@ HIT_TYPE_SHOWER_SCATTERED
Scattered light from primary shower.
@ HIT_TYPE_BREMSSTRAHLUNG_DIRECT
Direct light from Bremsstrahlung.
bool is_noise(const Hit &hit)
Verify hit origin.
const Trk & get_tau(const Evt &evt)
Get first tau from the event tracklist.
bool is_lepton(const Trk &track)
Test whether given track is a lepton.
bool from_electron(const Hit &hit)
Test whether given hit was produced by an electron.
bool is_proton(const Trk &track)
Test whether given track is a (anti-)proton.
bool is_baryon(const Trk &track)
Test whether given track is a baryon (is hadron and third digit of PDG code is not zero)
@ GEANT4_TYPE_ANTIMUON
@ GEANT4_TYPE_MUON
@ GEANT4_TYPE_ANTITAU
@ GEANT4_TYPE_ANTIELECTRON
@ GEANT4_TYPE_ELECTRON
@ GEANT4_TYPE_TAU
double getNPE(const Hit &hit)
Get true charge of hit.
Vec getP1(const Evt &evt)
Get momentum vector of the final state of a neutrino interaction.
bool is_meson(const Trk &track)
Test whether given track is a meson (is hadron and third digit of PDG code is zero)
bool has_particleID(const Trk &track, const int type)
Test whether given track corresponds to given particle type.
bool has_leptonic_taudecay(const Evt &event)
Test whether given event contains a leptonic tau-decay.
bool has_hadronic_taudecay(const Evt &event)
Test whether given event contains a hadronic tau-decay.
bool has_muonic_taudecay(const Evt &event)
Test whether given event contains a tau-lepton decay into a muon.
bool is_hadron(const Trk &track)
Test whether given track is a hadron.
@ TRACK_TYPE_PROTON
@ TRACK_TYPE_NUTAU
@ TRACK_TYPE_ANTINEUTRON
@ TRACK_TYPE_SIGMA_PLUS
@ TRACK_TYPE_CHARGED_PION_PLUS
@ TRACK_TYPE_LAMBDA
@ TRACK_TYPE_ANTIPROTON
@ TRACK_TYPE_ANTILAMBDA
@ TRACK_TYPE_NEUTRAL_PION
@ TRACK_TYPE_NUE
@ TRACK_TYPE_ELECTRON
@ TRACK_TYPE_NEUTRAL_SIGMA
@ TRACK_TYPE_PHOTON
@ TRACK_TYPE_MUON
@ TRACK_TYPE_NEUTRON
@ TRACK_TYPE_NUMU
@ TRACK_TYPE_SIGMA_MINUS
bool is_tau(const Trk &track)
Test whether given track is a (anti-)tau.
int count_hadrons(const Evt &evt)
Count the number of hadrons in a given event (not including neutral pions)
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [s]).
std::string to_string(const T &value)
Convert value to string.
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
static const double MASS_MUON
muon mass [GeV]
static const JGeaneWater gWater
Function object for energy loss of muon in sea water.
Definition: JGeane.hh:381
static const double MASS_NEUTRON
neutron mass [GeV]
double getKineticEnergy(const double E, const double m)
Get kinetic energy of particle with given energy and mass.
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [ns]).
static const double MASS_PROTON
proton mass [GeV]
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
const int n
Definition: JPolint.hh:786
Definition: JSTDTypes.hh:14
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:21
std::vector< Hit > mc_hits
MC: list of MC truth hits.
Definition: Evt.hh:48
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:49
int id
offline event identifier
Definition: Evt.hh:22
static int ROOT_IO_VERSION
Streamer version as obtained from ROOT file.
Definition: Evt.hh:274
std::vector< Trk > trks
list of reconstructed tracks (can be several because of prefits,showers, etc).
Definition: Evt.hh:39
Definition: Hit.hh:10
double a
hit amplitude (in p.e.)
Definition: Hit.hh:24
int type
particle type or parametrisation used for hit (mc only)
Definition: Hit.hh:28
double t
hit time (from tdc+calibration or MC truth)
Definition: Hit.hh:23
const JParticle & getPDG(const int pdg) const
Get particle corresponding to given PDG code.
Definition: JPDB.hh:240
static const JPDB & getInstance()
Get particle data book.
Definition: JPDB.hh:131
double mass
mass of particle [GeV]
Definition: JPDB.hh:97
int charge
charge of particle [|e|/3]
Definition: JPDB.hh:98
Primary particle.
Definition: JHead.hh:1174
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition: Trk.hh:15
int status
MC status code, see km3net-dataformat/definitions/trkmembers.csv for values.
Definition: Trk.hh:28
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
int id
track identifier
Definition: Trk.hh:16
std::vector< double > fitinf
place to store additional fit info, see km3net-dataformat/definitions/fitparameters....
Definition: Trk.hh:32
Vec dir
track direction
Definition: Trk.hh:18
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:20
double t
track time [ns] (when the particle is at pos )
Definition: Trk.hh:19
int mother_id
MC id of the parent particle.
Definition: Trk.hh:29
Vec pos
postion [m] of the track at time t
Definition: Trk.hh:17
The Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition: Vec.hh:13
double z
Definition: Vec.hh:14
double x
Definition: Vec.hh:14
double y
Definition: Vec.hh:14
Auxiliary methods for selection of reconstructed tracks.
static const int PDG_MUONBUNDLE
muon bundle reached the can level (mupage)
Definition: trkmembers.hh:38
static const int TRK_ST_FINALSTATE
for MC: the particle must be processed by detector simulation ('track_in' tag in evt files)....
Definition: trkmembers.hh:15
static const int TRK_MOTHER_UNDEFINED
KM3NeT Data Definitions v3.4.0-8-ge14cb17 https://git.km3net.de/common/km3net-dataformat.
Definition: trkmembers.hh:12
static const int TRK_ST_MUONBUNDLE
initial state muon bundle (mupage)
Definition: trkmembers.hh:18
static const int TRK_ST_PRIMARYCOSMIC
initial state cosmic ray ('track_primary' tag in evt files from corant).
Definition: trkmembers.hh:17
static const int TRK_ST_PRIMARYNEUTRINO
initial state neutrino ('neutrino' tag in evt files from gseagen and genhen).
Definition: trkmembers.hh:16
static const int TRK_ST_ININUCLEI
Initial state nuclei (gseagen)
Definition: trkmembers.hh:19
static const int TRK_MOTHER_NONE
mother id of a particle if it has no parent
Definition: trkmembers.hh:13