Jpp  19.0.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JAAnetToolkit.hh
Go to the documentation of this file.
1 #ifndef __JAANET__JAANETTOOLKIT__
2 #define __JAANET__JAANETTOOLKIT__
3 
4 #include <cstdlib>
5 #include <algorithm>
6 
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 interaction types based on GENIE codes.
75  */
77  UNKNOWN = 0,
81  };
82 
83 
84  /**
85  * Enumeration of hit types based on km3 codes.
86  */
87  enum JHitType_t {
88  HIT_TYPE_MUON_DIRECT = +5, //!< Direct light from muon
89  HIT_TYPE_MUON_SCATTERED = -5, //!< Scattered light from muon
90  HIT_TYPE_DELTARAYS_DIRECT = +4, //!< Direct light from delta-rays
91  HIT_TYPE_DELTARAYS_SCATTERED = -4, //!< Scattered light from delta-rays
92  HIT_TYPE_BREMSSTRAHLUNG_DIRECT = +3, //!< Direct light from Bremsstrahlung
93  HIT_TYPE_BREMSSTRAHLUNG_SCATTERED = -3, //!< Scattered light from Bremsstrahlung
94  HIT_TYPE_SHOWER_DIRECT = +99, //!< Direct light from primary shower
95  HIT_TYPE_SHOWER_SCATTERED = -99, //!< Scattered light from primary shower
96  HIT_TYPE_NOISE = -1, //!< Random noise
97  HIT_TYPE_UNKNOWN = 0 //!< Unknown source
98  };
99 
100 
101  /**
102  * Get true time of hit.
103  *
104  * \param hit hit
105  * \return true time of hit
106  */
107  inline double getTime(const Hit& hit)
108  {
109  return hit.t;
110  }
111 
112 
113  /**
114  * Get true charge of hit.
115  *
116  * \param hit hit
117  * \return true charge of hit
118  */
119  inline double getNPE(const Hit& hit)
120  {
121  return hit.a;
122  }
123 
124 
125  /**
126  * Verify hit origin.
127  *
128  * \param hit hit
129  * \return true if noise; else false
130  */
131  inline bool is_noise(const Hit& hit)
132  {
133  return hit.type == HIT_TYPE_NOISE;
134  }
135 
136 
137  /**
138  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
139  * Note that the global event time in not taken into account.
140  *
141  * \param event event
142  * \return time range
143  */
144  inline JTimeRange getTimeRange(const Evt& event)
145  {
147 
148  for (std::vector<Hit>::const_iterator hit = event.mc_hits.begin(); hit != event.mc_hits.end(); ++hit) {
149  if (!is_noise(*hit)) {
150  time_range.include(getTime(*hit));
151  }
152  }
153 
154  return time_range;
155  }
156 
157 
158  /**
159  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
160  * The resulting time range is larger than or equal to the given time window.\n
161  * Note that the global event time in not taken into account.
162  *
163  * \param event event
164  * \param T_ns time window [ns]
165  * \return time range
166  */
167  inline JTimeRange getTimeRange(const Evt& event, const JTimeRange& T_ns)
168  {
169  JTimeRange time_range = getTimeRange(event);
170 
171  if (!time_range.is_valid()) {
172  time_range = T_ns;
173  }
174 
175  if (time_range.getLength() < T_ns.getLength()) {
176  time_range.setLength(T_ns.getLength());
177  }
178 
179  return time_range;
180  }
181 
182 
183  /**
184  * Get position.
185  *
186  * \param pos position
187  * \return position
188  */
189  inline JPosition3D getPosition(const Vec& pos)
190  {
191  return JPosition3D(pos.x, pos.y, pos.z);
192  }
193 
194 
195  /**
196  * Get position.
197  *
198  * \param pos position
199  * \return position
200  */
201  inline Vec getPosition(const JPosition3D& pos)
202  {
203  return Vec(pos.getX(), pos.getY(), pos.getZ());
204  }
205 
206 
207  /**
208  * Get position.
209  *
210  * \param track track
211  * \return position
212  */
213  inline JPosition3D getPosition(const Trk& track)
214  {
215  return getPosition(track.pos);
216  }
217 
218 
219  /**
220  * Get direction.
221  *
222  * \param dir direction
223  * \return direction
224  */
225  inline JDirection3D getDirection(const Vec& dir)
226  {
227  return JDirection3D(dir.x, dir.y, dir.z);
228  }
229 
230 
231  /**
232  * Get direction.
233  *
234  * \param dir direction
235  * \return direction
236  */
237  inline Vec getDirection(const JDirection3D& dir)
238  {
239  return Vec(dir.getDX(), dir.getDY(), dir.getDZ());
240  }
241 
242 
243  /**
244  * Get direction.
245  *
246  * \param track track
247  * \return direction
248  */
249  inline JDirection3D getDirection(const Trk& track)
250  {
251  return getDirection(track.dir);
252  }
253 
254 
255  /**
256  * Get axis.
257  *
258  * \param track track
259  * \return axis
260  */
261  inline JAxis3D getAxis(const Trk& track)
262  {
263  return JAxis3D(getPosition(track), getDirection(track));
264  }
265 
266 
267  /**
268  * Get track.
269  *
270  * \param track track
271  * \return track
272  */
273  inline JTrack3E getTrack(const Trk& track)
274  {
275  return JTrack3E(JTrack3D(getAxis(track), track.t), track.E);
276  }
277 
278 
279  /**
280  * Get transformation.
281  *
282  * \param track track
283  * \return transformation
284  */
286  {
287  return JTransformation3D(getPosition(track), getDirection(track));
288  }
289 
290 
291  /**
292  * Get track information.
293  *
294  * \param track track
295  * \param index index
296  * \param value default value
297  * \return actual value
298  */
299  inline double getW(const Trk& track, const size_t index, const double value)
300  {
301  if (index < track.fitinf.size())
302  return track.fitinf[index];
303  else
304  return value;
305  }
306 
307 
308  /**
309  * Test whether given track is a photon or neutral pion.
310  *
311  * \param track track
312  * \return true if photon or neutral pion; else false
313  */
314  inline bool is_photon (const Trk& track) { return ( track.type == TRACK_TYPE_PHOTON ||
315  abs(track.type) == TRACK_TYPE_NEUTRAL_PION); }
316 
317  /**
318  * Test whether given track is a neutrino.
319  *
320  * \param track track
321  * \return true if neutrino; else false
322  */
323  inline bool is_neutrino(const Trk& track) { return (abs(track.type) == TRACK_TYPE_NUE ||
324  abs(track.type) == TRACK_TYPE_NUMU ||
325  abs(track.type) == TRACK_TYPE_NUTAU); }
326 
327  /**
328  * Test whether given track is a (anti-)electron.
329  *
330  * \param track track
331  * \return true if (anti-)electron; else false
332  */
333  inline bool is_electron(const Trk& track) { return (abs(track.type) == TRACK_TYPE_ELECTRON); }
334 
335  /**
336  * Test whether given track is a (anti-)muon.
337  *
338  * \param track track
339  * \return true if (anti-)muon; else false
340  */
341  inline bool is_muon (const Trk& track) { return (abs(track.type) == TRACK_TYPE_MUON); }
342 
343  /**
344  * Test whether given track is a (anti-)tau.
345  *
346  * \param track track
347  * \return true if (anti-)tau; else false
348  */
349  inline bool is_tau (const Trk& track) { return (abs(track.type) == TRACK_TYPE_TAU); }
350 
351  /**
352  * Test whether given track is a charged pion.
353  *
354  * \param track track
355  * \return true if charged pion; else false
356  */
357  inline bool is_pion (const Trk& track) { return (abs(track.type) == TRACK_TYPE_CHARGED_PION_PLUS); }
358 
359  /**
360  * Test whether given track is a (anti-)proton.
361  *
362  * \param track track
363  * \return true if (anti-)proton; else false
364  */
365  inline bool is_proton (const Trk& track) { return (abs(track.type) == TRACK_TYPE_PROTON); }
366 
367  /**
368  * Test whether given track is a (anti-)neutron.
369  *
370  * \param track track
371  * \return true if (anti-)neutron; else false
372  */
373  inline bool is_neutron (const Trk& track) { return (abs(track.type) == TRACK_TYPE_NEUTRON); }
374 
375  /**
376  * Test whether given track is a lepton
377  *
378  * \param track track
379  * \return true if lepton; else fails
380  */
381  inline bool is_lepton (const Trk& track) { return (is_neutrino(track) ||
382  is_electron(track) ||
383  is_muon (track) ||
384  is_tau (track)); }
385 
386  /**
387  * Test whether given track is a charged lepton
388  *
389  * \param track track
390  * \return true if lepton; else fails
391  */
392  inline bool is_charged_lepton (const Trk& track) { return (is_electron(track) ||
393  is_muon (track) ||
394  is_tau (track)); }
395 
396  /**
397  * Test whether given track is a hadron
398  *
399  * \param track track
400  * \return true if hadron; else fails
401  */
402  inline bool is_hadron (const Trk& track) { return (abs(track.type) != TRACK_TYPE_PHOTON &&
403  !is_lepton(track)); }
404 
405  /**
406  * Test whether given track is a meson (is hadron and third digit of PDG code is zero)
407  *
408  * \param track track
409  * \return true if hadron; else fails
410  */
411  inline bool is_meson (const Trk& track) { return (is_hadron(track) &&
412  ((int)(track.type / 1000)) % 10 == 0); }
413 
414  /**
415  * Test whether given track is a baryon (is hadron and third digit of PDG code is not zero)
416  *
417  * \param track track
418  * \return true if hadron; else fails
419  */
420  inline bool is_baryon (const Trk& track) { return (is_hadron(track) &&
421  ((int)(track.type / 1000)) % 10 != 0); }
422 
423  /**
424  * Test whether given track corresponds to given particle type.
425  *
426  * \param track track
427  * \param type particle type
428  * \return true if matches the corresponding particle; else false
429  */
430  inline bool has_particleID (const Trk& track, const int type) { return (track.type == type); }
431 
432  /**
433  * Test whether given track corresponds to an initial state particle.
434  *
435  * \param track track
436  * \return true if track corresponds to an initial state particle; else false
437  */
438  inline bool is_initialstate(const Trk& track)
439  {
440  if (Evt::ROOT_IO_VERSION >= 14) {
441 
442  return (track.status == TRK_ST_PRIMARYNEUTRINO ||
443  track.status == TRK_ST_PRIMARYCOSMIC ||
444  track.status == TRK_ST_MUONBUNDLE ||
445  track.status == TRK_ST_ININUCLEI);
446 
447  } else if (Evt::ROOT_IO_VERSION >= 11) {
448 
449  return (track.mother_id == TRK_MOTHER_UNDEFINED ||
450  track.mother_id == TRK_MOTHER_NONE);
451 
452  } else {
453 
454  return is_neutrino(track) && track.id == 0;
455  }
456  }
457 
458  /**
459  * Test whether given track corresponds to a final state particle.
460  *
461  * \param track track
462  * \return true if track corresponds to a final state particle; else false
463  */
464  inline bool is_finalstate(const Trk& track)
465  {
466  if (Evt::ROOT_IO_VERSION >= 15) {
467  return track.status == TRK_ST_FINALSTATE;
468  } else {
469  return !is_initialstate(track);
470  }
471  }
472 
473  /**
474  * Test whether given event has an incoming neutrino.
475  *
476  * \param evt event
477  * \return true if neutrino; else false
478  */
479  inline bool has_neutrino(const Evt& evt)
480  {
481  if (Evt::ROOT_IO_VERSION >= 14) {
482 
483  std::vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
485  return i != evt.mc_trks.cend();
486 
487  } else {
488 
489  return !evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0]);
490  }
491  }
492 
493  /**
494  * Get incoming neutrino.
495  *
496  * \param evt event
497  * \return neutrino
498  */
499  inline const Trk& get_neutrino(const Evt& evt)
500  {
501  if (Evt::ROOT_IO_VERSION >= 14) {
502 
503  std::vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
505 
506  if (i != evt.mc_trks.cend()) { return *i; }
507 
508  } else if (!evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0])) {
509 
510  return evt.mc_trks[0];
511  }
512 
513  THROW(JIndexOutOfRange, "JAANET::get_neutrino(): Event " << evt.id << " has no neutrino.");
514  }
515 
516  /**
517  * Get primary.
518  *
519  * \param evt event
520  * \return primary track
521  */
522  inline const Trk& get_primary(const Evt& evt)
523  {
524  using namespace std;
525  using namespace JPP;
526 
527  for (vector<Trk>::const_iterator i = evt.mc_trks.cbegin(); i != evt.mc_trks.cend(); ++i) {
528  if (is_initialstate(*i) && i->status != TRK_ST_ININUCLEI) {
529  return *i;
530  }
531  }
532 
533  THROW(JIndexOutOfRange, "JAANET::get_primary(): Cannot retrieve primary track for event " << evt.id << ".");
534  }
535 
536 
537  /**
538  * Get vertex.
539  *
540  * \param track track
541  * \return vertex
542  */
543  inline JVertex3D getVertex(const Trk& track)
544  {
545  return JVertex3D(getPosition(track), track.t);
546  }
547 
548  /**
549  * Get event vertex.
550  *
551  * \param event event
552  * \return event vertex
553  */
554  inline JVertex3D getVertex(const Evt& event)
555  {
556  using namespace std;
557  using namespace JPP;
558 
559  if (has_neutrino(event)) {
560 
561  const Trk& neutrino = get_neutrino(event);
562 
563  return getVertex(neutrino);
564 
565  } else if (!event.mc_trks.empty()) {
566 
567  vector<Trk>::const_iterator i = find_if(event.mc_trks.begin(), event.mc_trks.end(), &is_initialstate);
568 
569  if (i != event.mc_trks.end()) {
570  return getVertex(*i);
571  } else {
572  THROW(JValueOutOfRange, "getVertex(): No initial state track found.");
573  }
574 
575  } else if (!event.trks.empty()) { // For reconstructed DAQ events
576 
577  const Trk& track = get_best_reconstructed_track<JPP_RECONSTRUCTION_TYPE>(event);
578 
579  return getVertex(track);
580 
581  } else {
582 
583  THROW(JValueOutOfRange, "getVertex(): Could not find a valid event vertex.");
584  }
585  }
586 
587  /**
588  * Count the number of electrons in a given event
589  *
590  * \param evt event
591  * \return number of electrons
592  */
593  inline int count_electrons(const Evt& evt)
594  {
595  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
596  }
597 
598  /**
599  * Test whether given event has an electron.
600  *
601  * \param evt event
602  * \return true if event has electron; else false
603  */
604  inline bool has_electron(const Evt& evt)
605  {
606  return count_electrons(evt) != 0;
607  }
608 
609  /**
610  * Get first electron from the event tracklist
611  *
612  * \param evt event
613  * \return electron
614  */
615  inline const Trk& get_electron(const Evt& evt)
616  {
617  if (count_electrons(evt) > 0)
618  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
619  else
620  THROW(JIndexOutOfRange, "This event has no electron.");
621  }
622 
623  /**
624  * Test whether given hit was produced by an electron
625  *
626  * Warning: This method will only work with the output of KM3Sim.
627  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
628  *
629  * \param hit hit
630  * \return true if hit was produced by an electron; else false
631  */
632  inline bool from_electron(const Hit& hit)
633  {
634  return hit.type == GEANT4_TYPE_ELECTRON || hit.type == GEANT4_TYPE_ANTIELECTRON;
635  }
636 
637  /**
638  * Count the number of muons in a given event
639  *
640  * \param evt event
641  * \return number of muons
642  */
643  inline int count_muons(const Evt& evt)
644  {
645  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
646  }
647 
648  /**
649  * Test whether given event has a muon.
650  *
651  * \param evt event
652  * \return true if event has muon; else false
653  */
654  inline bool has_muon(const Evt& evt)
655  {
656  return count_muons(evt) != 0;
657  }
658 
659  /**
660  * Get first muon from the event tracklist
661  *
662  * \param evt event
663  * \return muon
664  */
665  inline const Trk& get_muon(const Evt& evt)
666  {
667  if (count_muons(evt) > 0)
668  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
669  else
670  THROW(JIndexOutOfRange, "This event has no muon.");
671  }
672 
673  /**
674  * Test whether given hit was produced by a muon
675  *
676  * Warning: This method will only work with the output of KM3Sim.
677  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
678  *
679  * \param hit hit
680  * \return true if hit was produced by a muon; else false
681  */
682  inline bool from_muon(const Hit& hit)
683  {
684  return hit.type == GEANT4_TYPE_MUON || hit.type == GEANT4_TYPE_ANTIMUON;
685  }
686 
687  /**
688  * Count the number of taus in a given event
689  *
690  * \param evt event
691  * \return number of taus
692  */
693  inline int count_taus(const Evt& evt)
694  {
695  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
696  }
697 
698  /**
699  * Test whether given event has a tau.
700  *
701  * \param evt event
702  * \return true if event has tau; else false
703  */
704  inline bool has_tau(const Evt& evt)
705  {
706  return count_taus(evt) != 0;
707  }
708 
709  /**
710  * Get first tau from the event tracklist
711  *
712  * \param evt event
713  * \return tau
714  */
715  inline const Trk& get_tau(const Evt& evt)
716  {
717  if (count_taus(evt) > 0)
718  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
719  else
720  THROW(JIndexOutOfRange, "This event has no tau.");
721  }
722 
723 
724  /**
725  * Test whether given event contains a tau-lepton decay into a muon.
726  *
727  * \param event event
728  * \return true if event contains a muonic tau-lepton decay; else false
729  */
730  inline bool has_muonic_taudecay(const Evt& event)
731  {
732  using namespace std;
733 
734  if (has_tau(event)) {
735 
736  const Trk& tau = get_tau(event);
737 
738  for (vector<Trk>::const_iterator i = event.mc_trks.cbegin(); i != event.mc_trks.cend(); ++i) {
739  if (is_muon(*i) && i->mother_id == tau.id) {
740  return true;
741  }
742  }
743  }
744 
745  return false;
746  }
747 
748  /**
749  * Test whether given event contains a leptonic tau-decay.
750  *
751  * \param event event
752  * \return true if event contains a leptonic tau-decay; else false
753  */
754  inline bool has_leptonic_taudecay(const Evt& event)
755  {
756  using namespace std;
757 
758  if (has_tau(event)) {
759 
760  const Trk& tau = get_tau(event);
761 
762  for (vector<Trk>::const_iterator i = event.mc_trks.cbegin(); i != event.mc_trks.cend(); ++i) {
763  if ((is_electron(*i) || is_muon(*i)) && i->mother_id == tau.id) {
764  return true;
765  }
766  }
767  }
768 
769  return false;
770  }
771 
772  /**
773  * Test whether given event contains a hadronic tau-decay.
774  *
775  * \param event event
776  * \return true if event contains a hadronic tau-decay; else false
777  */
778  inline bool has_hadronic_taudecay(const Evt& event)
779  {
780  return has_tau(event) && !has_leptonic_taudecay(event);
781  }
782 
783  /**
784  * Count the number of tau-lepton decay prongs in a given event.\n
785  * The number of prongs is defined as the number of charged tau-lepton decay products.
786  *
787  * \param event event
788  * \return number of tau-lepton decay prongs
789  */
790  inline int count_taudecay_prongs(const Evt& event)
791  {
792  using namespace std;
793 
794  int n = 0;
795 
796  if (has_tau(event)) {
797 
798  const Trk& tau = get_tau(event);
799 
800  for (vector<Trk>::const_iterator i = event.mc_trks.cbegin(); i != event.mc_trks.cend(); ++i) {
801  if (i->mother_id == tau.id && JPDB::getInstance().getPDG(i->type).charge != 0) {
802  ++n;
803  }
804  }
805  }
806 
807  return n;
808  }
809 
810  /**
811  * Test whether given hit was produced by a tau
812  *
813  * Warning: This method will only work with the output of KM3Sim.
814  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
815  *
816  * \param hit hit
817  * \return true if hit was produced by a tau; else false
818  */
819  inline bool from_tau(const Hit& hit)
820  {
821  return hit.type == GEANT4_TYPE_TAU || hit.type == GEANT4_TYPE_ANTITAU;
822  }
823 
824  /**
825  * Count the number of hadrons in a given event (not including neutral pions)
826  *
827  * \param evt event
828  * \return number of hadrons
829  */
830  inline int count_hadrons(const Evt& evt)
831  {
832  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_hadron);
833  }
834 
835  /**
836  * Test whether given hit was produced by a hadronic shower
837  *
838  * Warning: This method will only work with the output of KM3Sim.
839  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
840  *
841  * \param hit hit
842  * \return true if hit was produced by a hadron; else false
843  */
844  inline bool from_hadron(const Hit& hit)
845  {
846  return (!from_electron(hit) && !from_muon(hit) && !from_tau(hit));
847  }
848 
849 
850  /**
851  * Add time offset to mc event;
852  * according to current definition, the absolute time of the event is defined by the track "t" attribute;
853  * this could change in the future if the global attribute mc_t is assigned to this purpose.
854  *
855  * \param evt event
856  * \param tOff time offset [ns]
857  */
858  inline void add_time_offset(Evt& evt, const double tOff)
859  {
860  for (std::vector<Trk>::iterator p = evt.mc_trks.begin(); p != evt.mc_trks.end(); p++) {
861  p->t += tOff;
862  }
863  }
864 
865 
866  /**
867  * Get track kinetic energy.
868  *
869  * \param trk track
870  * \return kinetic energy [GeV]
871  */
872  inline double getKineticEnergy(const Trk& trk)
873  {
874  const double energy = trk.E;
875  const double mass = JPDB::getInstance().getPDG(trk.type).mass;
876 
877  return getKineticEnergy(energy, mass);
878  }
879 
880 
881  /**
882  * Get initial state energy of a neutrino interaction.\n
883  * This method includes baryon number conservation.
884  *
885  * \param evt event
886  * \return energy [GeV]
887  */
888  inline double getE0(const Evt& evt)
889  {
890  using namespace std;
891 
892  const Trk& neutrino = get_neutrino(evt);
893 
894  double E0 = neutrino.E;
895 
896  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.end(); ++track) {
897 
898  if (!is_finalstate(*track)) { continue; }
899 
900  if (track->type == TRACK_TYPE_NEUTRON ||
901  track->type == TRACK_TYPE_SIGMA_MINUS ||
902  track->type == TRACK_TYPE_NEUTRAL_SIGMA) {
903 
904  E0 += MASS_NEUTRON;
905 
906  } else if (track->type == TRACK_TYPE_PROTON ||
907  track->type == TRACK_TYPE_LAMBDA ||
908  track->type == TRACK_TYPE_SIGMA_PLUS) {
909 
910  E0 += MASS_PROTON;
911 
912  } else if (track->type == TRACK_TYPE_ANTINEUTRON) {
913 
914  E0 -= MASS_NEUTRON;
915 
916  } else if (track->type == TRACK_TYPE_ANTIPROTON ||
917  track->type == TRACK_TYPE_ANTILAMBDA) {
918 
919  E0 -= MASS_PROTON;
920  }
921  }
922 
923  return E0;
924  }
925 
926 
927  /**
928  * Get final state energy of a neutrino interaction.\n
929  * This method includes muon energy loss.
930  *
931  * \param evt event
932  * \return energy [GeV]
933  */
934  inline double getE1(const Evt& evt)
935  {
936  using namespace std;
937  using namespace JPP;
938 
939  double E1 = 0.0;
940 
941  const Trk& neutrino = get_neutrino(evt);
942 
943  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
944 
945  if (!is_finalstate(*track)) { continue; }
946 
947  if (is_muon(*track)) {
948 
949  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
950  evt.mc_trks[track->mother_id] : neutrino );
951 
952  const double distance = (track->pos - mother.pos).len();
953 
954  E1 += gWater(track->E, -distance);
955 
956  } else {
957 
958  E1 += track->E;
959  }
960  }
961 
962  return E1;
963  }
964 
965 
966  /**
967  * Get momentum vector of the initial state of a neutrino interaction.\n
968  * This method assumes neutrino DIS on a stationary nucleus
969  *
970  * \param evt event
971  * \return energy [GeV]
972  */
973  inline Vec getP0(const Evt& evt)
974  {
975  const Trk& neutrino = get_neutrino(evt);
976 
977  return neutrino.E * neutrino.dir;
978  }
979 
980 
981  /**
982  * Get momentum vector of the final state of a neutrino interaction.\n
983  * This method includes muon energy losses.
984  *
985  * \param evt event
986  * \return final state momentum vector [GeV]
987  */
988  inline Vec getP1(const Evt& evt)
989  {
990  using namespace std;
991  using namespace JPP;
992 
993  Vec P1(0,0,0);
994 
995  const Trk& neutrino = get_neutrino(evt);
996 
997  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
998 
999  if (!is_finalstate(*track)) { continue; }
1000 
1001  double kineticEnergy = 0.0;
1002 
1003  if (is_muon(*track)) {
1004 
1005  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
1006  evt.mc_trks[track->mother_id] : neutrino );
1007 
1008  const double distance = (track->pos - mother.pos).len();
1009  const double energy = gWater(track->E, -distance);
1010 
1011  kineticEnergy = getKineticEnergy(energy, MASS_MUON);
1012 
1013  } else {
1014 
1015  kineticEnergy = getKineticEnergy(*track);
1016  }
1017 
1018  P1 += kineticEnergy * track->dir;
1019  }
1020 
1021  return P1;
1022  }
1023 
1024 
1025  /**
1026  * Get final state invariant mass.
1027  *
1028  * \param event event
1029  * \return invariant mass [GeV]
1030  */
1031  inline double getInvariantMass(const Evt& event)
1032  {
1033  using namespace std;
1034  using namespace JPP;
1035 
1036  double Minv = 0.0;
1037 
1038  for (vector<Trk>::const_iterator track = event.mc_trks.cbegin(); track != event.mc_trks.cend(); ++track) {
1039 
1040  if (!is_finalstate(*track)) { continue; }
1041 
1042  const double mass = JPDB::getInstance().getPDG(track->type).mass;
1043 
1044  Minv += mass;
1045  }
1046 
1047  return Minv;
1048  }
1049 }
1050 
1051 #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.
JInteractionType_t
Enumeration of interaction types based on GENIE codes.
Auxiliary methods for physics calculations.
double getInvariantMass(const Evt &event)
Get final state invariant mass.
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
int charge
charge of particle [|e|/3]
Definition: JPDB.hh:98
bool has_leptonic_taudecay(const Evt &event)
Test whether given event contains a leptonic tau-decay.
int count_taudecay_prongs(const Evt &event)
Count the number of tau-lepton decay prongs in a given event.
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...
double getKineticEnergy(const double E, const double m)
Get kinetic energy of particle with given energy and mass.
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.
const Trk & get_primary(const Evt &evt)
Get primary.
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:240
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.
bool has_muonic_taudecay(const Evt &event)
Test whether given event contains a tau-lepton decay into a muon.
static const int AASHOWER_RECONSTRUCTION_TYPE
AAnet shower fit reconstruction type.
static const JPDB & getInstance()
Get particle data book.
Definition: JPDB.hh:131
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.
const int n
Definition: JPolint.hh:786
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
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
bool has_hadronic_taudecay(const Evt &event)
Test whether given event contains a hadronic tau-decay.
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
double getKineticEnergy(const Trk &trk)
Get track kinetic energy.
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.4.0-8-ge14cb17 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:97
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
std::vector< Trk > trks
list of reconstructed tracks (can be several because of prefits,showers, etc).
Definition: Evt.hh:39
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
then set_variable DETECTOR set_variable OUTPUT_FILE set_variable DAQ_FILE set_variable PMT_FILE else fatal Wrong number of arguments fi JPrintTree f $DAQ_FILE type
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.