Jpp  17.1.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JAAnetToolkit.hh
Go to the documentation of this file.
1 #ifndef __JAANET__JAANETTOOLKIT__
2 #define __JAANET__JAANETTOOLKIT__
3 
4 #include <cstdlib>
5 #include <algorithm>
6 
14 
15 #include "JLang/JException.hh"
16 #include "JLang/JPredicate.hh"
17 
20 #include "JGeometry3D/JAxis3D.hh"
21 #include "JGeometry3D/JTrack3D.hh"
22 #include "JGeometry3D/JTrack3E.hh"
23 #include "JGeometry3D/JVertex3D.hh"
25 
26 #include "JPhysics/JConstants.hh"
27 #include "JPhysics/JTimeRange.hh"
28 #include "JPhysics/JGeane.hh"
29 
30 #include "JAAnet/JParticleTypes.hh"
31 #include "JAAnet/JPDB.hh"
32 
33 
34 /**
35  * \file
36  *
37  * Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
38  * \author mdejong
39  */
40 namespace JAANET {}
41 namespace JPP { using namespace JAANET; }
42 
43 /**
44  * Extensions to Evt data format.
45  */
46 namespace JAANET {
47 
51 
59 
63 
64  /**
65  * AAnet shower fit reconstruction type.
66  */
67  static const int AASHOWER_RECONSTRUCTION_TYPE = 101;
68 
69  /**
70  * Enumeration of hit types based on km3 codes.
71  */
72  enum JHitType_t {
73  HIT_TYPE_MUON_DIRECT = +5, //!< Direct light from muon
74  HIT_TYPE_MUON_SCATTERED = -5, //!< Scattered light from muon
75  HIT_TYPE_DELTARAYS_DIRECT = +4, //!< Direct light from delta-rays
76  HIT_TYPE_DELTARAYS_SCATTERED = -4, //!< Scattered light from delta-rays
77  HIT_TYPE_BREMSSTRAHLUNG_DIRECT = +3, //!< Direct light from Bremsstrahlung
78  HIT_TYPE_BREMSSTRAHLUNG_SCATTERED = -3, //!< Scattered light from Bremsstrahlung
79  HIT_TYPE_SHOWER_DIRECT = +99, //!< Direct light from primary shower
80  HIT_TYPE_SHOWER_SCATTERED = -99, //!< Scattered light from primary shower
81  HIT_TYPE_NOISE = -1, //!< Random noise
82  HIT_TYPE_UNKNOWN = 0 //!< Unknown source
83  };
84 
85 
86  /**
87  * Get true time of hit.
88  *
89  * \param hit hit
90  * \return true time of hit
91  */
92  inline double getTime(const Hit& hit)
93  {
94  return hit.t;
95  }
96 
97 
98  /**
99  * Get true charge of hit.
100  *
101  * \param hit hit
102  * \return true charge of hit
103  */
104  inline double getNPE(const Hit& hit)
105  {
106  return hit.a;
107  }
108 
109 
110  /**
111  * Verify hit origin.
112  *
113  * \param hit hit
114  * \return true if noise; else false
115  */
116  inline bool is_noise(const Hit& hit)
117  {
118  return hit.type == HIT_TYPE_NOISE;
119  }
120 
121 
122  /**
123  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
124  * Note that the global event time in not taken into account.
125  *
126  * \param event event
127  * \return time range
128  */
129  inline JTimeRange getTimeRange(const Evt& event)
130  {
132 
133  for (std::vector<Hit>::const_iterator hit = event.mc_hits.begin(); hit != event.mc_hits.end(); ++hit) {
134  if (!is_noise(*hit)) {
135  time_range.include(getTime(*hit));
136  }
137  }
138 
139  return time_range;
140  }
141 
142 
143  /**
144  * Get time range (i.e.\ time between earliest and latest hit) of Monte Carlo event.\n
145  * The resulting time range is larger than or equal to the given time window.\n
146  * Note that the global event time in not taken into account.
147  *
148  * \param event event
149  * \param T_ns time window [ns]
150  * \return time range
151  */
152  inline JTimeRange getTimeRange(const Evt& event, const JTimeRange& T_ns)
153  {
154  JTimeRange time_range = getTimeRange(event);
155 
156  if (!time_range.is_valid()) {
157  time_range = T_ns;
158  }
159 
160  if (time_range.getLength() < T_ns.getLength()) {
161  time_range.setLength(T_ns.getLength());
162  }
163 
164  return time_range;
165  }
166 
167 
168  /**
169  * Get position.
170  *
171  * \param pos position
172  * \return position
173  */
174  inline JPosition3D getPosition(const Vec& pos)
175  {
176  return JPosition3D(pos.x, pos.y, pos.z);
177  }
178 
179 
180  /**
181  * Get position.
182  *
183  * \param pos position
184  * \return position
185  */
186  inline Vec getPosition(const JPosition3D& pos)
187  {
188  return Vec(pos.getX(), pos.getY(), pos.getZ());
189  }
190 
191 
192  /**
193  * Get position.
194  *
195  * \param track track
196  * \return position
197  */
198  inline JPosition3D getPosition(const Trk& track)
199  {
200  return getPosition(track.pos);
201  }
202 
203 
204  /**
205  * Get direction.
206  *
207  * \param dir direction
208  * \return direction
209  */
210  inline JDirection3D getDirection(const Vec& dir)
211  {
212  return JDirection3D(dir.x, dir.y, dir.z);
213  }
214 
215 
216  /**
217  * Get direction.
218  *
219  * \param dir direction
220  * \return direction
221  */
222  inline Vec getDirection(const JDirection3D& dir)
223  {
224  return Vec(dir.getDX(), dir.getDY(), dir.getDZ());
225  }
226 
227 
228  /**
229  * Get direction.
230  *
231  * \param track track
232  * \return direction
233  */
234  inline JDirection3D getDirection(const Trk& track)
235  {
236  return getDirection(track.dir);
237  }
238 
239 
240  /**
241  * Get axis.
242  *
243  * \param track track
244  * \return axis
245  */
246  inline JAxis3D getAxis(const Trk& track)
247  {
248  return JAxis3D(getPosition(track), getDirection(track));
249  }
250 
251 
252  /**
253  * Get track.
254  *
255  * \param track track
256  * \return track
257  */
258  inline JTrack3E getTrack(const Trk& track)
259  {
260  return JTrack3E(JTrack3D(getAxis(track), track.t), track.E);
261  }
262 
263 
264  /**
265  * Get vertex.
266  *
267  * \param track track
268  * \return vertex
269  */
270  inline JVertex3D getVertex(const Trk& track)
271  {
272  return JVertex3D(getPosition(track), track.t);
273  }
274 
275 
276  /**
277  * Get transformation.
278  *
279  * \param track track
280  * \return transformation
281  */
283  {
284  return JTransformation3D(getPosition(track), getDirection(track));
285  }
286 
287 
288  /**
289  * Get track information.
290  *
291  * \param track track
292  * \param index index
293  * \param value default value
294  * \return actual value
295  */
296  inline double getW(const Trk& track, const size_t index, const double value)
297  {
298  if (index < track.fitinf.size())
299  return track.fitinf[index];
300  else
301  return value;
302  }
303 
304 
305  /**
306  * Test whether given track is a photon or neutral pion.
307  *
308  * \param track track
309  * \return true if photon or neutral pion; else false
310  */
311  inline bool is_photon (const Trk& track) { return ( track.type == TRACK_TYPE_PHOTON ||
312  abs(track.type) == TRACK_TYPE_NEUTRAL_PION); }
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 corresponds to given particle type.
422  *
423  * \param track track
424  * \param type particle type
425  * \return true if matches the corresponding particle; else false
426  */
427  inline bool has_particleID (const Trk& track, const int type) { return (track.type == type); }
428 
429  /**
430  * Test whether given event has an incoming neutrino.
431  *
432  * \param evt event
433  * \return true if neutrino; else false
434  */
435  inline bool has_neutrino(const Evt& evt)
436  {
437  return !evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0]);
438  }
439 
440  /**
441  * Get incoming neutrino.
442  *
443  * \param evt event
444  * \return neutrino
445  */
446  inline const Trk& get_neutrino(const Evt& evt)
447  {
448  if (has_neutrino(evt)) {
449 
450  std::vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
452 
453  return (i != evt.mc_trks.cend() ? *i : *find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
455  } else {
456 
457  THROW(JIndexOutOfRange, "JAANET::get_neutrino(): Event " << evt.id << " has no neutrino.");
458  }
459  }
460 
461  /**
462  * Count the number of electrons in a given event
463  *
464  * \param evt event
465  * \return number of electrons
466  */
467  inline int count_electrons(const Evt& evt)
468  {
469  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
470  }
471 
472  /**
473  * Test whether given event has an electron.
474  *
475  * \param evt event
476  * \return true if event has electron; else false
477  */
478  inline bool has_electron(const Evt& evt)
479  {
480  return count_electrons(evt) != 0;
481  }
482 
483  /**
484  * Get first electron from the event tracklist
485  *
486  * \param evt event
487  * \return electron
488  */
489  inline const Trk& get_electron(const Evt& evt)
490  {
491  if (count_electrons(evt) > 0)
492  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
493  else
494  THROW(JIndexOutOfRange, "This event has no electron.");
495  }
496 
497  /**
498  * Test whether given hit was produced by an electron
499  *
500  * Warning: This method will only work with the output of KM3Sim.
501  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
502  *
503  * \param hit hit
504  * \return true if hit was produced by an electron; else false
505  */
506  inline bool from_electron(const Hit& hit)
507  {
508  return hit.type == GEANT4_TYPE_ELECTRON || hit.type == GEANT4_TYPE_ANTIELECTRON;
509  }
510 
511  /**
512  * Count the number of muons in a given event
513  *
514  * \param evt event
515  * \return number of muons
516  */
517  inline int count_muons(const Evt& evt)
518  {
519  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
520  }
521 
522  /**
523  * Test whether given event has a muon.
524  *
525  * \param evt event
526  * \return true if event has muon; else false
527  */
528  inline bool has_muon(const Evt& evt)
529  {
530  return count_muons(evt) != 0;
531  }
532 
533  /**
534  * Get first muon from the event tracklist
535  *
536  * \param evt event
537  * \return muon
538  */
539  inline const Trk& get_muon(const Evt& evt)
540  {
541  if (count_muons(evt) > 0)
542  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
543  else
544  THROW(JIndexOutOfRange, "This event has no muon.");
545  }
546 
547  /**
548  * Test whether given hit was produced by a muon
549  *
550  * Warning: This method will only work with the output of KM3Sim.
551  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
552  *
553  * \param hit hit
554  * \return true if hit was produced by a muon; else false
555  */
556  inline bool from_muon(const Hit& hit)
557  {
558  return hit.type == GEANT4_TYPE_MUON || hit.type == GEANT4_TYPE_ANTIMUON;
559  }
560 
561  /**
562  * Count the number of taus in a given event
563  *
564  * \param evt event
565  * \return number of taus
566  */
567  inline int count_taus(const Evt& evt)
568  {
569  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
570  }
571 
572  /**
573  * Test whether given event has a tau.
574  *
575  * \param evt event
576  * \return true if event has tau; else false
577  */
578  inline bool has_tau(const Evt& evt)
579  {
580  return count_taus(evt) != 0;
581  }
582 
583  /**
584  * Get first tau from the event tracklist
585  *
586  * \param evt event
587  * \return tau
588  */
589  inline const Trk& get_tau(const Evt& evt)
590  {
591  if (count_taus(evt) > 0)
592  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
593  else
594  THROW(JIndexOutOfRange, "This event has no tau.");
595  }
596 
597  /**
598  * Test whether given hit was produced by a tau
599  *
600  * Warning: This method will only work with the output of KM3Sim.
601  * For JSirene or KM3, you should check if track.id == hit.origin instead.
602  *
603  * \param hit hit
604  * \return true if hit was produced by a tau; else false
605  */
606  inline bool from_tau(const Hit& hit)
607  {
608  return hit.type == GEANT4_TYPE_TAU || hit.type == GEANT4_TYPE_ANTITAU;
609  }
610 
611  /**
612  * Count the number of hadrons in a given event (not including neutral pions)
613  *
614  * \param evt event
615  * \return number of hadrons
616  */
617  inline int count_hadrons(const Evt& evt)
618  {
619  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_hadron);
620  }
621 
622  /**
623  * Test whether given hit was produced by a hadronic shower
624  *
625  * Warning: This method will only work with the output of KM3Sim.
626  * For JSirene or KM3, you should check if track.id == hit.origin instead.
627  *
628  * \param hit hit
629  * \return true if hit was produced by a hadron; else false
630  */
631  inline bool from_hadron(const Hit& hit)
632  {
633  return (!from_electron(hit) && !from_muon(hit) && !from_tau(hit));
634  }
635 
636 
637  /**
638  * Add time offset to mc event;
639  * according to current definition, the absolute time of the event is defined by the track "t" attribute;
640  * this could change in the future if the global attribute mc_t is assigned to this purpose.
641  *
642  * \param evt event
643  * \param tOff time offset [ns]
644  */
645  inline void add_time_offset(Evt& evt, const double tOff)
646  {
647  for (std::vector<Trk>::iterator p = evt.mc_trks.begin(); p != evt.mc_trks.end(); p++) {
648  p->t += tOff;
649  }
650  }
651 
652 
653  /**
654  * Get kinetic energy of particle with given mass.
655  *
656  * \param E energy [GeV]
657  * \param m mass [GeV]
658  * \return energy [GeV]
659  */
660  inline double getKineticEnergy(const double E, const double m)
661  {
662  if (E > m) {
663  return sqrt((E - m) * (E + m));
664  } else {
665  return 0.0;
666  }
667  }
668 
669 
670  /**
671  * Get track kinetic energy.
672  *
673  * \param trk track
674  * \return kinetic energy [GeV]
675  */
676  inline double getKineticEnergy(const Trk& trk)
677  {
678  const double energy = trk.E;
679  const double mass = JPDB::getInstance().getPDG(trk.type).mass;
680 
681  return getKineticEnergy(energy, mass);
682  }
683 
684 
685  /**
686  * Get initial state energy of a neutrino interaction.\n
687  * This method includes baryon number conservation.
688  *
689  * Note: If the given event does not correspond to a neutrino interaction\n
690  * this method returns 0.0 GeV.
691  *
692  * \param evt event
693  * \return energy [GeV]
694  */
695  inline double getE0(const Evt& evt)
696  {
697  double E0 = 0.0;
698 
699  if (!evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0])) {
700 
701  const Trk& neutrino = evt.mc_trks[0];
702 
703  E0 += neutrino.E;
704 
705  for (size_t i = 1; i != evt.mc_trks.size(); ++i) {
706 
707  const Trk& track = evt.mc_trks[i];
708 
709  if (track.type == TRACK_TYPE_NEUTRON ||
710  track.type == TRACK_TYPE_SIGMA_MINUS) {
711  E0 += MASS_NEUTRON;
712  }
713 
714  if (track.type == TRACK_TYPE_PROTON ||
715  track.type == TRACK_TYPE_LAMBDA ||
716  track.type == TRACK_TYPE_SIGMA_PLUS) {
717  E0 += MASS_PROTON;
718  }
719 
720  if (track.type == TRACK_TYPE_ANTINEUTRON) {
721  E0 -= MASS_NEUTRON;
722  }
723 
724  if (track.type == TRACK_TYPE_ANTIPROTON ||
725  track.type == TRACK_TYPE_ANTILAMBDA) {
726  E0 -= MASS_PROTON;
727  }
728  }
729  }
730 
731  return E0;
732  }
733 
734 
735  /**
736  * Get final state energy of a neutrino interaction.\n
737  * This method includes muon energy loss.
738  *
739  * \param evt event
740  * \return energy [GeV]
741  */
742  inline double getE1(const Evt& evt)
743  {
744  double E1 = 0.0;
745 
746  if (!evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0])) {
747 
748  const Trk& neutrino = evt.mc_trks[0];
749 
750  for (size_t i = 1; i != evt.mc_trks.size(); ++i) {
751 
752  const Trk& track = evt.mc_trks[i];
753 
754  if (is_muon(track)) {
755 
756  const double distance = (track.pos - neutrino.pos).len();
757 
758  E1 += JPHYSICS::gWater(track.E, -distance);
759 
760  } else {
761 
762  E1 += track.E;
763  }
764  }
765  }
766 
767  return E1;
768  }
769 
770 
771  /**
772  * Get momentum vector of the initial state of a neutrino interaction.\n
773  * This method assumes neutrino DIS on a stationary nucleus
774  *
775  * \param evt event
776  * \return energy [GeV]
777  */
778  inline Vec getP0(const Evt& evt)
779  {
780  Vec P0(0,0,0);
781 
782  if (!evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0])) {
783 
784  const Trk& neutrino = evt.mc_trks[0];
785 
786  P0 = neutrino.E * neutrino.dir;
787  }
788 
789  return P0;
790  }
791 
792 
793  /**
794  * Get momentum vector of the final state of a neutrino interaction.\n
795  * This method includes muon energy losses.
796  *
797  * Note: If the given event does not correspond to a neutrino interaction\n
798  * a null vector will be returned.
799  *
800  * \param evt event
801  * \return final state momentum vector [GeV]
802  */
803  inline Vec getP1(const Evt& evt)
804  {
805  Vec P1(0,0,0);
806 
807  if (!evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0])) {
808 
809  const Trk& neutrino = evt.mc_trks[0];
810 
811  for (size_t i = 1; i != evt.mc_trks.size(); ++i) {
812 
813  const Trk& track = evt.mc_trks[i];
814 
815  double kineticEnergy = 0.0;
816 
817  if (is_muon(track)) {
818 
819  const double distance = (track.pos - neutrino.pos).len();
820  const double energy = JPHYSICS::gWater(track.E, -distance);
821 
822  kineticEnergy = getKineticEnergy(energy, JPHYSICS::MASS_MUON);
823 
824  } else {
825 
826  kineticEnergy = getKineticEnergy(track);
827  }
828 
829  P1 += kineticEnergy * track.dir;
830  }
831  }
832 
833  return P1;
834  }
835 }
836 
837 #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.
static const JRange< double, std::less< double > > DEFAULT_RANGE
Default range.
Definition: JRange.hh:555
bool is_electron(const Trk &track)
Test whether given track is a (anti-)electron.
int count_electrons(const Evt &evt)
Count the number of electrons in a given event.
JHitType_t
Enumeration of hit types based on km3 codes.
Exceptions.
Scattered light from primary shower.
JPredicate< JResult_t T::*, JComparison::eq > make_predicate(JResult_t T::*member, const JResult_t value)
Helper method to create predicate for data member.
Definition: JPredicate.hh:128
bool is_proton(const Trk &track)
Test whether given track is a (anti-)proton.
const Trk & get_muon(const Evt &evt)
Get first muon from the event tracklist.
Data structure for direction in three dimensions.
Definition: JDirection3D.hh:33
then usage E
Definition: JMuonPostfit.sh:35
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.
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_neutrino() const
Check if this is a netrino.
Definition: Trk.hh:100
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:696
double getE1(const Evt &evt)
Get final state energy of a neutrino interaction.
static const double MASS_MUON
muon mass [GeV]
bool has_particleID(const Trk &track, const int type)
Test whether given track corresponds to given particle type.
bool is_charged_lepton(const Trk &track)
Test whether given track is a charged lepton.
Scattered light from Bremsstrahlung.
static const JGeaneWater gWater
Function object for energy loss of muon in sea water.
Definition: JGeane.hh:381
const JParticle & getPDG(const int pdg) const
Get particle corresponding to given PDG code.
Definition: JPDB.hh:234
Direct light from Bremsstrahlung.
double getTime(const Hit &hit)
Get true time of hit.
bool is_baryon(const Trk &track)
Test whether given track is a baryon (is hadron and third digit of PDG code is not zero) ...
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:20
double y
Definition: Vec.hh:14
3D track with energy.
Definition: JTrack3E.hh:30
Direct light from primary shower.
static const int AASHOWER_RECONSTRUCTION_TYPE
AAnet shower fit reconstruction type.
static const JPDB & getInstance()
Get particle data book.
Definition: JPDB.hh:125
bool is_noise(const Hit &hit)
Verify hit origin.
JTimeRange getTimeRange(const Evt &event)
Get time range (i.e. time between earliest and latest hit) of Monte Carlo event.
static const double MASS_NEUTRON
neutron mass [GeV]
bool is_neutrino(const Trk &track)
Test whether given track is a neutrino.
Axis object.
Definition: JAxis3D.hh:38
bool is_lepton(const Trk &track)
Test whether given track is a lepton.
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, for jgandalf, see JFitParameters.hh
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
double getY() const
Get y position.
Definition: JVector3D.hh:104
bool is_photon(const Trk &track)
Test whether given track is a photon or neutral pion.
bool from_hadron(const Hit &hit)
Test whether given hit was produced by a hadronic shower.
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
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
bool is_pion(const Trk &track)
Test whether given track is a charged pion.
Vec pos
postion of the track at time t [m]
Definition: Trk.hh:17
double mass
mass of particle [GeV]
Definition: JPDB.hh:92
int count_taus(const Evt &evt)
Count the number of taus in a given event.
Definition: Hit.hh:8
const Trk & get_tau(const Evt &evt)
Get first tau from the event tracklist.
Direct light from muon.
then P1
bool has_muon(const Evt &evt)
Test whether given event has a muon.
std::vector< Hit > mc_hits
MC: list of MC truth hits.
Definition: Evt.hh:48
double t
hit time (from tdc+calibration or MC truth)
Definition: Hit.hh:23
Scattered light from delta-rays.
static const double MASS_PROTON
proton mass [GeV]
Definition of particle types.
double getX() const
Get x position.
Definition: JVector3D.hh:94
bool has_tau(const Evt &evt)
Test whether given event has a tau.
int id
offline event identifier
Definition: Evt.hh:22
double getNPE(const Hit &hit)
Get true charge of hit.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
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:162
Exception for accessing an index in a collection that is outside of its range.
Definition: JException.hh:90
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
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.