Jpp  18.3.0-209-g56ce19a
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  /**
72  * Enumeration of interaction types based on GENIE codes.
73  */
75  UNKNOWN = 0,
79  };
80 
81 
82  /**
83  * Enumeration of hit types based on km3 codes.
84  */
85  enum JHitType_t {
86  HIT_TYPE_MUON_DIRECT = +5, //!< Direct light from muon
87  HIT_TYPE_MUON_SCATTERED = -5, //!< Scattered light from muon
88  HIT_TYPE_DELTARAYS_DIRECT = +4, //!< Direct light from delta-rays
89  HIT_TYPE_DELTARAYS_SCATTERED = -4, //!< Scattered light from delta-rays
90  HIT_TYPE_BREMSSTRAHLUNG_DIRECT = +3, //!< Direct light from Bremsstrahlung
91  HIT_TYPE_BREMSSTRAHLUNG_SCATTERED = -3, //!< Scattered light from Bremsstrahlung
92  HIT_TYPE_SHOWER_DIRECT = +99, //!< Direct light from primary shower
93  HIT_TYPE_SHOWER_SCATTERED = -99, //!< Scattered light from primary shower
94  HIT_TYPE_NOISE = -1, //!< Random noise
95  HIT_TYPE_UNKNOWN = 0 //!< Unknown source
96  };
97 
98 
99  /**
100  * Get true time of hit.
101  *
102  * \param hit hit
103  * \return true time of hit
104  */
105  inline double getTime(const Hit& hit)
106  {
107  return hit.t;
108  }
109 
110 
111  /**
112  * Get true charge of hit.
113  *
114  * \param hit hit
115  * \return true charge of hit
116  */
117  inline double getNPE(const Hit& hit)
118  {
119  return hit.a;
120  }
121 
122 
123  /**
124  * Verify hit origin.
125  *
126  * \param hit hit
127  * \return true if noise; else false
128  */
129  inline bool is_noise(const Hit& hit)
130  {
131  return hit.type == HIT_TYPE_NOISE;
132  }
133 
134 
135  /**
136  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
137  * Note that the global event time in not taken into account.
138  *
139  * \param event event
140  * \return time range
141  */
142  inline JTimeRange getTimeRange(const Evt& event)
143  {
145 
146  for (std::vector<Hit>::const_iterator hit = event.mc_hits.begin(); hit != event.mc_hits.end(); ++hit) {
147  if (!is_noise(*hit)) {
148  time_range.include(getTime(*hit));
149  }
150  }
151 
152  return time_range;
153  }
154 
155 
156  /**
157  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
158  * The resulting time range is larger than or equal to the given time window.\n
159  * Note that the global event time in not taken into account.
160  *
161  * \param event event
162  * \param T_ns time window [ns]
163  * \return time range
164  */
165  inline JTimeRange getTimeRange(const Evt& event, const JTimeRange& T_ns)
166  {
167  JTimeRange time_range = getTimeRange(event);
168 
169  if (!time_range.is_valid()) {
170  time_range = T_ns;
171  }
172 
173  if (time_range.getLength() < T_ns.getLength()) {
174  time_range.setLength(T_ns.getLength());
175  }
176 
177  return time_range;
178  }
179 
180 
181  /**
182  * Get position.
183  *
184  * \param pos position
185  * \return position
186  */
187  inline JPosition3D getPosition(const Vec& pos)
188  {
189  return JPosition3D(pos.x, pos.y, pos.z);
190  }
191 
192 
193  /**
194  * Get position.
195  *
196  * \param pos position
197  * \return position
198  */
199  inline Vec getPosition(const JPosition3D& pos)
200  {
201  return Vec(pos.getX(), pos.getY(), pos.getZ());
202  }
203 
204 
205  /**
206  * Get position.
207  *
208  * \param track track
209  * \return position
210  */
211  inline JPosition3D getPosition(const Trk& track)
212  {
213  return getPosition(track.pos);
214  }
215 
216 
217  /**
218  * Get direction.
219  *
220  * \param dir direction
221  * \return direction
222  */
223  inline JDirection3D getDirection(const Vec& dir)
224  {
225  return JDirection3D(dir.x, dir.y, dir.z);
226  }
227 
228 
229  /**
230  * Get direction.
231  *
232  * \param dir direction
233  * \return direction
234  */
235  inline Vec getDirection(const JDirection3D& dir)
236  {
237  return Vec(dir.getDX(), dir.getDY(), dir.getDZ());
238  }
239 
240 
241  /**
242  * Get direction.
243  *
244  * \param track track
245  * \return direction
246  */
247  inline JDirection3D getDirection(const Trk& track)
248  {
249  return getDirection(track.dir);
250  }
251 
252 
253  /**
254  * Get axis.
255  *
256  * \param track track
257  * \return axis
258  */
259  inline JAxis3D getAxis(const Trk& track)
260  {
261  return JAxis3D(getPosition(track), getDirection(track));
262  }
263 
264 
265  /**
266  * Get track.
267  *
268  * \param track track
269  * \return track
270  */
271  inline JTrack3E getTrack(const Trk& track)
272  {
273  return JTrack3E(JTrack3D(getAxis(track), track.t), track.E);
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 >= 15) {
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  * Get primary.
516  *
517  * \param evt event
518  * \return primary track
519  */
520  inline const Trk& get_primary(const Evt& evt)
521  {
522  using namespace std;
523  using namespace JPP;
524 
525  for (vector<Trk>::const_iterator i = evt.mc_trks.cbegin(); i != evt.mc_trks.cend(); ++i) {
526  if (is_initialstate(*i) && i->status != TRK_ST_ININUCLEI) {
527  return *i;
528  }
529  }
530 
531  THROW(JIndexOutOfRange, "JAANET::get_primary(): Cannot retrieve primary track for event " << evt.id << ".");
532  }
533 
534 
535  /**
536  * Get vertex.
537  *
538  * \param track track
539  * \return vertex
540  */
541  inline JVertex3D getVertex(const Trk& track)
542  {
543  return JVertex3D(getPosition(track), track.t);
544  }
545 
546  /**
547  * Get event vertex.
548  *
549  * \param event event
550  * \return event vertex
551  */
552  inline JVertex3D getVertex(const Evt& event)
553  {
554  using namespace std;
555  using namespace JPP;
556 
557  if (has_neutrino(event)) {
558 
559  const Trk& neutrino = get_neutrino(event);
560 
561  return getVertex(neutrino);
562 
563  } else if (!event.mc_trks.empty()) {
564 
565  vector<Trk>::const_iterator i = find_if(event.mc_trks.begin(), event.mc_trks.end(), &is_initialstate);
566 
567  if (i != event.mc_trks.end()) {
568  return getVertex(*i);
569  } else {
570  THROW(JValueOutOfRange, "getVertex(): No initial state track found.");
571  }
572 
573  } else if (!event.trks.empty()) { // For reconstructed DAQ events
574 
575  const Trk& track = get_best_reconstructed_track<JPP_RECONSTRUCTION_TYPE>(event);
576 
577  return getVertex(track);
578 
579  } else {
580 
581  THROW(JValueOutOfRange, "getVertex(): Could not find a valid event vertex.");
582  }
583  }
584 
585  /**
586  * Count the number of electrons in a given event
587  *
588  * \param evt event
589  * \return number of electrons
590  */
591  inline int count_electrons(const Evt& evt)
592  {
593  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
594  }
595 
596  /**
597  * Test whether given event has an electron.
598  *
599  * \param evt event
600  * \return true if event has electron; else false
601  */
602  inline bool has_electron(const Evt& evt)
603  {
604  return count_electrons(evt) != 0;
605  }
606 
607  /**
608  * Get first electron from the event tracklist
609  *
610  * \param evt event
611  * \return electron
612  */
613  inline const Trk& get_electron(const Evt& evt)
614  {
615  if (count_electrons(evt) > 0)
616  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
617  else
618  THROW(JIndexOutOfRange, "This event has no electron.");
619  }
620 
621  /**
622  * Test whether given hit was produced by an electron
623  *
624  * Warning: This method will only work with the output of KM3Sim.
625  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
626  *
627  * \param hit hit
628  * \return true if hit was produced by an electron; else false
629  */
630  inline bool from_electron(const Hit& hit)
631  {
632  return hit.type == GEANT4_TYPE_ELECTRON || hit.type == GEANT4_TYPE_ANTIELECTRON;
633  }
634 
635  /**
636  * Count the number of muons in a given event
637  *
638  * \param evt event
639  * \return number of muons
640  */
641  inline int count_muons(const Evt& evt)
642  {
643  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
644  }
645 
646  /**
647  * Test whether given event has a muon.
648  *
649  * \param evt event
650  * \return true if event has muon; else false
651  */
652  inline bool has_muon(const Evt& evt)
653  {
654  return count_muons(evt) != 0;
655  }
656 
657  /**
658  * Get first muon from the event tracklist
659  *
660  * \param evt event
661  * \return muon
662  */
663  inline const Trk& get_muon(const Evt& evt)
664  {
665  if (count_muons(evt) > 0)
666  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
667  else
668  THROW(JIndexOutOfRange, "This event has no muon.");
669  }
670 
671  /**
672  * Test whether given hit was produced by a muon
673  *
674  * Warning: This method will only work with the output of KM3Sim.
675  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
676  *
677  * \param hit hit
678  * \return true if hit was produced by a muon; else false
679  */
680  inline bool from_muon(const Hit& hit)
681  {
682  return hit.type == GEANT4_TYPE_MUON || hit.type == GEANT4_TYPE_ANTIMUON;
683  }
684 
685  /**
686  * Count the number of taus in a given event
687  *
688  * \param evt event
689  * \return number of taus
690  */
691  inline int count_taus(const Evt& evt)
692  {
693  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
694  }
695 
696  /**
697  * Test whether given event has a tau.
698  *
699  * \param evt event
700  * \return true if event has tau; else false
701  */
702  inline bool has_tau(const Evt& evt)
703  {
704  return count_taus(evt) != 0;
705  }
706 
707  /**
708  * Get first tau from the event tracklist
709  *
710  * \param evt event
711  * \return tau
712  */
713  inline const Trk& get_tau(const Evt& evt)
714  {
715  if (count_taus(evt) > 0)
716  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
717  else
718  THROW(JIndexOutOfRange, "This event has no tau.");
719  }
720 
721  /**
722  * Test whether given hit was produced by a tau
723  *
724  * Warning: This method will only work with the output of KM3Sim.
725  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
726  *
727  * \param hit hit
728  * \return true if hit was produced by a tau; else false
729  */
730  inline bool from_tau(const Hit& hit)
731  {
732  return hit.type == GEANT4_TYPE_TAU || hit.type == GEANT4_TYPE_ANTITAU;
733  }
734 
735  /**
736  * Count the number of hadrons in a given event (not including neutral pions)
737  *
738  * \param evt event
739  * \return number of hadrons
740  */
741  inline int count_hadrons(const Evt& evt)
742  {
743  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_hadron);
744  }
745 
746  /**
747  * Test whether given hit was produced by a hadronic shower
748  *
749  * Warning: This method will only work with the output of KM3Sim.
750  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
751  *
752  * \param hit hit
753  * \return true if hit was produced by a hadron; else false
754  */
755  inline bool from_hadron(const Hit& hit)
756  {
757  return (!from_electron(hit) && !from_muon(hit) && !from_tau(hit));
758  }
759 
760 
761  /**
762  * Add time offset to mc event;
763  * according to current definition, the absolute time of the event is defined by the track "t" attribute;
764  * this could change in the future if the global attribute mc_t is assigned to this purpose.
765  *
766  * \param evt event
767  * \param tOff time offset [ns]
768  */
769  inline void add_time_offset(Evt& evt, const double tOff)
770  {
771  for (std::vector<Trk>::iterator p = evt.mc_trks.begin(); p != evt.mc_trks.end(); p++) {
772  p->t += tOff;
773  }
774  }
775 
776 
777  /**
778  * Get kinetic energy of particle with given mass.
779  *
780  * \param E energy [GeV]
781  * \param m mass [GeV]
782  * \return energy [GeV]
783  */
784  inline double getKineticEnergy(const double E, const double m)
785  {
786  if (E > m) {
787  return sqrt((E - m) * (E + m));
788  } else {
789  return 0.0;
790  }
791  }
792 
793 
794  /**
795  * Get track kinetic energy.
796  *
797  * \param trk track
798  * \return kinetic energy [GeV]
799  */
800  inline double getKineticEnergy(const Trk& trk)
801  {
802  const double energy = trk.E;
803  const double mass = JPDB::getInstance().getPDG(trk.type).mass;
804 
805  return getKineticEnergy(energy, mass);
806  }
807 
808 
809  /**
810  * Get initial state energy of a neutrino interaction.\n
811  * This method includes baryon number conservation.
812  *
813  * \param evt event
814  * \return energy [GeV]
815  */
816  inline double getE0(const Evt& evt)
817  {
818  using namespace std;
819 
820  const Trk& neutrino = get_neutrino(evt);
821 
822  double E0 = neutrino.E;
823 
824  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.end(); ++track) {
825 
826  if (!is_finalstate(*track)) { continue; }
827 
828  if (track->type == TRACK_TYPE_NEUTRON ||
829  track->type == TRACK_TYPE_SIGMA_MINUS ||
830  track->type == TRACK_TYPE_NEUTRAL_SIGMA) {
831 
832  E0 += MASS_NEUTRON;
833 
834  } else if (track->type == TRACK_TYPE_PROTON ||
835  track->type == TRACK_TYPE_LAMBDA ||
836  track->type == TRACK_TYPE_SIGMA_PLUS) {
837 
838  E0 += MASS_PROTON;
839 
840  } else if (track->type == TRACK_TYPE_ANTINEUTRON) {
841 
842  E0 -= MASS_NEUTRON;
843 
844  } else if (track->type == TRACK_TYPE_ANTIPROTON ||
845  track->type == TRACK_TYPE_ANTILAMBDA) {
846 
847  E0 -= MASS_PROTON;
848  }
849  }
850 
851  return E0;
852  }
853 
854 
855  /**
856  * Get final state energy of a neutrino interaction.\n
857  * This method includes muon energy loss.
858  *
859  * \param evt event
860  * \return energy [GeV]
861  */
862  inline double getE1(const Evt& evt)
863  {
864  using namespace std;
865  using namespace JPP;
866 
867  double E1 = 0.0;
868 
869  const Trk& neutrino = get_neutrino(evt);
870 
871  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
872 
873  if (!is_finalstate(*track)) { continue; }
874 
875  if (is_muon(*track)) {
876 
877  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
878  evt.mc_trks[track->mother_id] : neutrino );
879 
880  const double distance = (track->pos - mother.pos).len();
881 
882  E1 += gWater(track->E, -distance);
883 
884  } else {
885 
886  E1 += track->E;
887  }
888  }
889 
890  return E1;
891  }
892 
893 
894  /**
895  * Get momentum vector of the initial state of a neutrino interaction.\n
896  * This method assumes neutrino DIS on a stationary nucleus
897  *
898  * \param evt event
899  * \return energy [GeV]
900  */
901  inline Vec getP0(const Evt& evt)
902  {
903  const Trk& neutrino = get_neutrino(evt);
904 
905  return neutrino.E * neutrino.dir;
906  }
907 
908 
909  /**
910  * Get momentum vector of the final state of a neutrino interaction.\n
911  * This method includes muon energy losses.
912  *
913  * \param evt event
914  * \return final state momentum vector [GeV]
915  */
916  inline Vec getP1(const Evt& evt)
917  {
918  using namespace std;
919  using namespace JPP;
920 
921  Vec P1(0,0,0);
922 
923  const Trk& neutrino = get_neutrino(evt);
924 
925  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
926 
927  if (!is_finalstate(*track)) { continue; }
928 
929  double kineticEnergy = 0.0;
930 
931  if (is_muon(*track)) {
932 
933  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
934  evt.mc_trks[track->mother_id] : neutrino );
935 
936  const double distance = (track->pos - mother.pos).len();
937  const double energy = gWater(track->E, -distance);
938 
939  kineticEnergy = getKineticEnergy(energy, MASS_MUON);
940 
941  } else {
942 
943  kineticEnergy = getKineticEnergy(*track);
944  }
945 
946  P1 += kineticEnergy * track->dir;
947  }
948 
949  return P1;
950  }
951 
952 
953  /**
954  * Get final state invariant mass.
955  *
956  * \param event event
957  * \return invariant mass [GeV]
958  */
959  inline double getInvariantMass(const Evt& event)
960  {
961  using namespace std;
962  using namespace JPP;
963 
964  double Minv = 0.0;
965 
966  for (vector<Trk>::const_iterator track = event.mc_trks.cbegin(); track != event.mc_trks.cend(); ++track) {
967 
968  if (!is_finalstate(*track)) { continue; }
969 
970  const double mass = JPDB::getInstance().getPDG(track->type).mass;
971 
972  Minv += mass;
973  }
974 
975  return Minv;
976  }
977 }
978 
979 #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.
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
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.
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: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.
then set_variable PMT_FILE set_variable DAQ_FILE set_variable OUTPUT_FILE set_variable DETECTOR else fatal Wrong number of arguments fi JPrintTree f $DAQ_FILE type
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.3.0-2-g5cc95cf 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
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
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.