Jpp
JAAnetToolkit.hh
Go to the documentation of this file.
1 #ifndef __JAANET__JAANETTOOLKIT__
2 #define __JAANET__JAANETTOOLKIT__
3 
4 #include <stdlib.h>
5 #include <limits>
6 #include <algorithm>
7 
8 #include "evt/Hit.hh"
9 #include "evt/Vec.hh"
10 #include "evt/Trk.hh"
11 #include "evt/Evt.hh"
12 
13 #include "JLang/JException.hh"
14 #include "JTools/JRange.hh"
17 #include "JGeometry3D/JAxis3D.hh"
18 #include "JGeometry3D/JTrack3D.hh"
19 #include "JGeometry3D/JTrack3E.hh"
20 #include "JGeometry3D/JVertex3D.hh"
22 #include "JTrigger/JHitL0.hh"
23 #include "JDetector/JTimeRange.hh"
25 #include "JDetector/JPMTRouter.hh"
26 #include "JFit/JFitApplications.hh"
27 #include "JAAnet/JParticleTypes.hh"
28 
29 
30 /**
31  * \file
32  *
33  * Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
34  * \author mdejong
35  */
36 namespace JAANET {}
37 namespace JPP { using namespace JAANET; }
38 
39 /**
40  * Extensions to AAnet data format.
41  */
42 namespace JAANET {
43 
52  using JTRIGGER::JHitL0;
57  using JTOOLS::JRange;
58 
59  /**
60  * AAshower reconstruction type for AAnet.
61  */
62  static const int AASHOWER_RECONSTRUCTION_TYPE = 101;
63 
64  /**
65  * Enumeration of hit types based on km3 codes.
66  */
67  enum JHitType_t {
68  HIT_TYPE_MUON_DIRECT = +5, //!< Direct light from muon
69  HIT_TYPE_MUON_SCATTERED = -5, //!< Scattered light from muon
70  HIT_TYPE_DELTARAYS_DIRECT = +4, //!< Direct light from delta-rays
71  HIT_TYPE_DELTARAYS_SCATTERED = -4, //!< Scattered light from delta-rays
72  HIT_TYPE_BREMSSTRAHLUNG_DIRECT = +3, //!< Direct light from Bremsstrahlung
73  HIT_TYPE_BREMSSTRAHLUNG_SCATTERED = -3, //!< Scattered light from Bremsstrahlung
74  HIT_TYPE_SHOWER_DIRECT = +99, //!< Direct light from primary shower
75  HIT_TYPE_SHOWER_SCATTERED = -99, //!< Scattered light from primary shower
76  HIT_TYPE_NOISE = -1, //!< Random noise
77  HIT_TYPE_UNKNOWN = 0 //!< Unknown source
78  };
79 
80 
81  /**
82  * Get true time of hit.
83  *
84  * \param hit hit
85  * \return true time of hit
86  */
87  inline double getTime(const Hit& hit)
88  {
89  // if IGAIN = 0 in km3 program, then pure_xx values are not set.
90 
91  if (hit.pure_a == 0.0)
92  return hit.t;
93  else
94  return hit.pure_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  // if IGAIN = 0 in km3 program, then pure_xx values are not set.
107 
108  if (hit.pure_a == 0.0)
109  return hit.a;
110  else
111  return hit.pure_a;
112  }
113 
114 
115  /**
116  * Verify hit origin.
117  *
118  * \param hit hit
119  * \return true if noise; else false
120  */
121  inline bool is_noise(const Hit& hit)
122  {
123  return hit.type == HIT_TYPE_NOISE;
124  }
125 
126 
127  /**
128  * Get time range (i.e. earlist and latest hit time) of Monte Carlo event.
129  * Note that the global event time in not taken into account.
130  *
131  * \param event event
132  * \return time range
133  */
134  inline JTimeRange getTimeRange(const Evt& event)
135  {
137 
138  for (std::vector<Hit>::const_iterator hit = event.mc_hits.begin(); hit != event.mc_hits.end(); ++hit) {
139  if (!is_noise(*hit)) {
140  timeRange.include(getTime(*hit));
141  }
142  }
143 
144  return timeRange;
145  }
146 
147 
148  /**
149  * Get time range (i.e. earlist and latest hit time) of Monte Carlo event.
150  * The given time window is offset so that the lower limit is equal to the earliest hit time.
151  * The time window is then applied to all hits.
152  * Note that the global event time in not taken into account.
153  *
154  * \param event event
155  * \param T_ns time window [ns]
156  * \return time range
157  */
158  inline JTimeRange getTimeRange(const Evt& event, const JTimeRange& T_ns)
159  {
160  using namespace std;
161 
163 
164  double t0 = numeric_limits<double>::max();
165 
166  for (std::vector<Hit>::const_iterator hit = event.mc_hits.begin(); hit != event.mc_hits.end(); ++hit) {
167  if (!is_noise(*hit)) {
168  if (getTime(*hit) < t0) {
169  t0 = getTime(*hit);
170  }
171  }
172  }
173 
174  if (t0 != numeric_limits<double>::max()) {
175 
176  JTimeRange window(T_ns);
177 
178  window.fixLowerLimit(t0);
179 
180  for (std::vector<Hit>::const_iterator hit = event.mc_hits.begin(); hit != event.mc_hits.end(); ++hit) {
181  if (!is_noise(*hit) && window(getTime(*hit))) {
182  timeRange.include(getTime(*hit));
183  }
184  }
185  }
186 
187  return timeRange;
188  }
189 
190 
191  /**
192  * Get position.
193  *
194  * \param v vector
195  * \return position
196  */
197  inline JPosition3D getPosition(const Vec& v)
198  {
199  return JPosition3D(v.x, v.y, v.z);
200  }
201 
202 
203  /**
204  * Get position.
205  *
206  * \param track track
207  * \return position
208  */
209  inline JPosition3D getPosition(const Trk& track)
210  {
211  return getPosition(track.pos);
212  }
213 
214 
215  /**
216  * Get direction.
217  *
218  * \param v vector
219  * \return direction
220  */
221  inline JDirection3D getDirection(const Vec& v)
222  {
223  return JDirection3D(v.x, v.y, v.z);
224  }
225 
226  /**
227  * Get direction.
228  *
229  * \param track track
230  * \return direction
231  */
232  inline JDirection3D getDirection(const Trk& track)
233  {
234  return getDirection(track.dir);
235  }
236 
237 
238  /**
239  * Get axis.
240  *
241  * \param track track
242  * \return axis
243  */
244  inline JAxis3D getAxis(const Trk& track)
245  {
246  return JAxis3D(getPosition(track), getDirection(track));
247  }
248 
249 
250  /**
251  * Get track.
252  *
253  * \param track track
254  * \return track
255  */
256  inline JTrack3E getTrack(const Trk& track)
257  {
258  return JTrack3E(JTrack3D(getAxis(track), track.t), track.E);
259  }
260 
261 
262  /**
263  * Get vertex.
264  *
265  * \param track track
266  * \return vertex
267  */
268  inline JVertex3D getVertex(const Trk& track)
269  {
270  return JVertex3D(getPosition(track), track.t);
271  }
272 
273 
274  /**
275  * Get transformation.
276  *
277  * \param track track
278  * \return transformation
279  */
280  inline JTransformation3D getTransformation(const Trk& track)
281  {
282  return JTransformation3D(getPosition(track), getDirection(track));
283  }
284 
285 
286  /**
287  * Get transformation.
288  *
289  * \param hit hit
290  * \return hit
291  */
292  inline JHitL0 getHit(const Hit& hit)
293  {
294  using namespace JPP;
295 
296  return JHitL0(JDAQPMTIdentifier(hit.dom_id, hit.channel_id), JAxis3D(getPosition(hit.pos), getDirection(hit.dir)), JHit(hit.t, hit.tot) );
297  }
298 
299 
300  /**
301  * Get transformation.
302  *
303  * \param hit hit
304  * \param router module router
305  * \return hit
306  */
307  inline JHitL0 getHit(const Hit& hit, const JModuleRouter& router)
308  {
309  using namespace JPP;
310 
311  const JPMT& pmt = router.getModule(hit.dom_id).getPMT(hit.channel_id);
312 
313  return JHitL0(JDAQPMTIdentifier(hit.dom_id, hit.channel_id), pmt.getAxis(), JHit(getTime(hit.tdc, pmt.getCalibration()), getToT(hit.tot, pmt.getCalibration())));
314  }
315 
316 
317  /**
318  * Get transformation.
319  *
320  * \param hit hit
321  * \param router pmt router
322  * \return hit
323  */
324  inline JHitL0 getHit(const Hit& hit, const JPMTRouter& router)
325  {
326  using namespace JPP;
327 
328  const JPMT& pmt = router.getPMT(hit.pmt_id);
329 
330  return JHitL0(JDAQPMTIdentifier(hit.dom_id, hit.channel_id), pmt.getAxis(), JHit(getTime(hit.tdc, pmt.getCalibration()), getToT(hit.tot, pmt.getCalibration())));
331  }
332 
333 
334  /**
335  * Test whether given track is a photon or neutral pion.
336  *
337  * \param track track
338  * \return true if photon or neutral pion; else false
339  */
340  inline bool is_photon (const Trk& track) { return ( track.type == TRACK_TYPE_PHOTON ||
341  abs(track.type) == TRACK_TYPE_NEUTRAL_PION); }
342 
343  /**
344  * Test whether given track is a neutrino.
345  *
346  * \param track track
347  * \return true if neutrino; else false
348  */
349  inline bool is_neutrino(const Trk& track) { return (abs(track.type) == TRACK_TYPE_NUE ||
350  abs(track.type) == TRACK_TYPE_NUMU ||
351  abs(track.type) == TRACK_TYPE_NUTAU); }
352 
353  /**
354  * Test whether given track is a (anti-)electron.
355  *
356  * \param track track
357  * \return true if (anti-)electron; else false
358  */
359  inline bool is_electron(const Trk& track) { return (abs(track.type) == TRACK_TYPE_ELECTRON); }
360 
361  /**
362  * Test whether given track is a (anti-)muon.
363  *
364  * \param track track
365  * \return true if (anti-)muon; else false
366  */
367  inline bool is_muon (const Trk& track) { return (abs(track.type) == TRACK_TYPE_MUON); }
368 
369  /**
370  * Test whether given track is a (anti-)tau.
371  *
372  * \param track track
373  * \return true if (anti-)tau; else false
374  */
375  inline bool is_tau (const Trk& track) { return (abs(track.type) == TRACK_TYPE_TAU); }
376 
377  /**
378  * Test whether given track is a charged pion.
379  *
380  * \param track track
381  * \return true if charged pion; else false
382  */
383  inline bool is_pion (const Trk& track) { return (abs(track.type) == TRACK_TYPE_CHARGED_PION_PLUS); }
384 
385  /**
386  * Test whether given track is a (anti-)proton.
387  *
388  * \param track track
389  * \return true if (anti-)proton; else false
390  */
391  inline bool is_proton (const Trk& track) { return (abs(track.type) == TRACK_TYPE_PROTON); }
392 
393  /**
394  * Test whether given track is a (anti-)neutron.
395  *
396  * \param track track
397  * \return true if (anti-)neutron; else false
398  */
399  inline bool is_neutron (const Trk& track) { return (abs(track.type) == TRACK_TYPE_NEUTRON); }
400 
401  /**
402  * Test whether given track is a hadron
403  *
404  * \param track track
405  * \return true if hadron; else fails
406  */
407  inline bool is_hadron (const Trk& track) { return (!is_neutrino(track) &&
408  !is_electron(track) &&
409  !is_muon (track) &&
410  !is_tau (track)); }
411 
412  /**
413  * Test whether given track corresponds to given particle type.
414  *
415  * \param track track
416  * \param type particle type
417  * \return true if matches the corresponding particle; else false
418  */
419  inline bool has_particleID (const Trk& track, const int type) { return (track.type == type); }
420 
421  /**
422  * Test whether given event has an incoming neutrino.
423  *
424  * \param evt event
425  * \return true if neutrino; else false
426  */
427  inline bool has_neutrino(const Evt& evt)
428  {
429  return !evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0]);
430  }
431 
432  /**
433  * Get incoming neutrino.
434  *
435  * \param evt event
436  * \return neutrino
437  */
438  inline const Trk& get_neutrino(const Evt& evt)
439  {
440  if (has_neutrino(evt))
441  return evt.mc_trks[0];
442  else
443  THROW(JIndexOutOfRange, "This event has no neutrino.");
444  }
445 
446  /**
447  * Count the number of electrons in a given event
448  *
449  * \param evt event
450  * \return number of electrons
451  */
452  inline int count_electrons(const Evt& evt)
453  {
454  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
455  }
456 
457  /**
458  * Test whether given event has an electron.
459  *
460  * \param evt event
461  * \return true if event has electron; else false
462  */
463  inline bool has_electron(const Evt& evt)
464  {
465  return count_electrons(evt) != 0;
466  }
467 
468  /**
469  * Get first electron from the event tracklist
470  *
471  * \param evt event
472  * \return electron
473  */
474  inline const Trk& get_electron(const Evt& evt)
475  {
476  if (count_electrons(evt) > 0)
477  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
478  else
479  THROW(JIndexOutOfRange, "This event has no electron.");
480  }
481 
482  /**
483  * Test whether given hit was produced by an electron
484  *
485  * Warning: This method will only work with the output of KM3Sim.
486  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
487  *
488  * \param hit hit
489  * \return true if hit was produced by an electron; else false
490  */
491  inline bool from_electron(const Hit& hit)
492  {
493  return hit.type == GEANT4_TYPE_ELECTRON || hit.type == GEANT4_TYPE_ANTIELECTRON;
494  }
495 
496  /**
497  * Count the number of muons in a given event
498  *
499  * \param evt event
500  * \return number of muons
501  */
502  inline int count_muons(const Evt& evt)
503  {
504  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
505  }
506 
507  /**
508  * Test whether given event has a muon.
509  *
510  * \param evt event
511  * \return true if event has muon; else false
512  */
513  inline bool has_muon(const Evt& evt)
514  {
515  return count_muons(evt) != 0;
516  }
517 
518  /**
519  * Get first muon from the event tracklist
520  *
521  * \param evt event
522  * \return muon
523  */
524  inline const Trk& get_muon(const Evt& evt)
525  {
526  if (count_muons(evt) > 0)
527  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
528  else
529  THROW(JIndexOutOfRange, "This event has no muon.");
530  }
531 
532  /**
533  * Test whether given hit was produced by a muon
534  *
535  * Warning: This method will only work with the output of KM3Sim.
536  * For JSirene or KM3, you should check if <tt>track.id == hit.origin</tt> instead.
537  *
538  * \param hit hit
539  * \return true if hit was produced by a muon; else false
540  */
541  inline bool from_muon(const Hit& hit)
542  {
543  return hit.type == GEANT4_TYPE_MUON || hit.type == GEANT4_TYPE_ANTIMUON;
544  }
545 
546  /**
547  * Count the number of taus in a given event
548  *
549  * \param evt event
550  * \return number of taus
551  */
552  inline int count_taus(const Evt& evt)
553  {
554  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
555  }
556 
557  /**
558  * Test whether given event has a tau.
559  *
560  * \param evt event
561  * \return true if event has tau; else false
562  */
563  inline bool has_tau(const Evt& evt)
564  {
565  return count_taus(evt) != 0;
566  }
567 
568  /**
569  * Get first tau from the event tracklist
570  *
571  * \param evt event
572  * \return tau
573  */
574  inline const Trk& get_tau(const Evt& evt)
575  {
576  if (count_taus(evt) > 0)
577  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
578  else
579  THROW(JIndexOutOfRange, "This event has no tau.");
580  }
581 
582  /**
583  * Test whether given hit was produced by a tau
584  *
585  * Warning: This method will only work with the output of KM3Sim.
586  * For JSirene or KM3, you should check if track.id == hit.origin instead.
587  *
588  * \param hit hit
589  * \return true if hit was produced by a tau; else false
590  */
591  inline bool from_tau(const Hit& hit)
592  {
593  return hit.type == GEANT4_TYPE_TAU || hit.type == GEANT4_TYPE_ANTITAU;
594  }
595 
596  /**
597  * Count the number of hadrons in a given event (not including neutral pions)
598  *
599  * \param evt event
600  * \return number of hadrons
601  */
602  inline int count_hadrons(const Evt& evt)
603  {
604  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_hadron);
605  }
606 
607  /**
608  * Test whether given hit was produced by a hadronic shower
609  *
610  * Warning: This method will only work with the output of KM3Sim.
611  * For JSirene or KM3, you should check if track.id == hit.origin instead.
612  *
613  * \param hit hit
614  * \return true if hit was produced by a hadron; else false
615  */
616  inline bool from_hadron(const Hit& hit)
617  {
618  return (!from_electron(hit) && !from_muon(hit) && !from_tau(hit));
619  }
620 
621  /**
622  * Auxiliary function to get average direction and total energy of a hadronic cascade
623  *
624  * \param evt event
625  * \return cascade (default if there are no hadrons)
626  */
627  inline Trk get_hadronic_cascade(const Evt& evt)
628  {
629  Trk cascade;
630 
631  for (vector<Trk>::const_iterator track = evt.mc_trks.begin(); track != evt.mc_trks.end(); ++track) {
632  if (is_hadron(*track)) {
633  cascade.E += track->E;
634  cascade.dir += track->dir * track->E;
635  cascade.t = track->t;
636  cascade.pos = track->pos;
637  }
638  }
639 
640  cascade.dir.normalize();
641 
642  return cascade;
643  }
644 
645  /**
646  * Add time offset to mc event;
647  * according to current definition, the absolute time of the event is defined by the track "t" attribute;
648  * this could change in the future if the global attribute mc_t is assigned to this purpose.
649  *
650  * \param evt event
651  * \param tOff time offset [ns]
652  */
653  inline void add_time_offset(Evt& evt, const double tOff)
654  {
655  for (vector<Trk>::iterator p = evt.mc_trks.begin(); p != evt.mc_trks.end(); p++) {
656  p->t += tOff;
657  }
658  }
659 
660  /**
661  * Reconstruction type dependent comparison of track quality.
662  *
663  * The specialisation of this class should implement the function object operator
664  * <pre>
665  * inline bool operator()(const Trk& first, const Trk& second) const
666  * </pre>
667  * and return true if the first track is better then the second track.
668  */
669  template<int reconstruction_type>
670  struct quality_sorter {
671  /**
672  * The default comparison is based on:
673  * -# number of reconstruction stages, i.e. <tt>Trk::rec_stages.size()</tt> (the larger the better);
674  * -# likelihood of the final track, i.e. <tt>Trk::lik</tt> (the larger the better).
675  *
676  * \param first first track
677  * \param second second track
678  * \return true if first track has better quality than second; else false
679  */
680  inline bool operator()(const Trk& first, const Trk& second) const
681  {
682  if (first.rec_stages.size() == second.rec_stages.size())
683  return first.lik > second.lik;
684  else
685  return first.rec_stages.size() > second.rec_stages.size();
686  }
687  };
688 
689  /**
690  * Auxiliary class to test whether given track has specified history.\n
691  * The history of a track consists of a reconstruction type and a range of application types.
692  */
693  struct has_history {
694  /**
695  * Constructor.
696  *
697  * \param type reconstruction type
698  * \param range range of application types
699  */
700  has_history(const int type, const JRange<int> range) :
701  type (type),
702  range(range)
703  {}
704 
705  /**
706  * \param track track
707  * \return true if given track has specified history; else false
708  */
709  inline bool operator()(const Trk& track) const
710  {
711  if (track.rec_type == type)
712  return std::find_if(track.rec_stages.begin(), track.rec_stages.end(), range) != track.rec_stages.end();
713  else
714  return false;
715  }
716 
717  const int type; //!< reconstruction type
718  const JRange<int> range; //!< range of application types
719  };
720 
721  /**
722  * Test whether given track has muon prefit in history.
723  *
724  * \param track track
725  * \return true if muon prefit in history; else false
726  */
727  inline bool has_muon_prefit(const Trk& track)
728  {
730  }
731 
732  /**
733  * Test whether given track has muon simplex fit in history.
734  *
735  * \param track track
736  * \return true if muon simplex fit in history; else false
737  */
738  inline bool has_muon_simplex(const Trk& track)
739  {
741  }
742 
743  /**
744  * Test whether given track has muon gandalf fit in history.
745  *
746  * \param track track
747  * \return true if muon gandalf fit in history; else false
748  */
749  inline bool has_muon_gandalf(const Trk& track)
750  {
752  }
753 
754  /**
755  * Test whether given track has muon energy fit in history.
756  *
757  * \param track track
758  * \return true if muon energy fit in history; else false
759  */
760  inline bool has_muon_energy(const Trk& track)
761  {
763  }
764 
765  /**
766  * Test whether given track has muon start fit in history.
767  *
768  * \param track track
769  * \return true if muon start fit in history; else false
770  */
771  inline bool has_muon_start(const Trk& track)
772  {
774  }
775 
776  /**
777  * Test whether given track has default muon fit in history.
778  *
779  * \param track track
780  * \return true if muon fit in history; else false
781  */
782  inline bool has_muon_fit(const Trk& track)
783  {
785  }
786 
787  /**
788  * Test whether given track has shower prefit in history.
789  *
790  * \param track track
791  * \return true if shower prefit in history; else false
792  */
793  inline bool has_shower_prefit(const Trk& track)
794  {
796  }
797 
798  /**
799  * Test whether given track has shower position fit in history.
800  *
801  * \param track track
802  * \return true if shower position fit in history; else false
803  */
804  inline bool has_shower_positionfit(const Trk& track)
805  {
807  }
808 
809  /**
810  * Test whether given track has shower complete fit in history.
811  *
812  * \param track track
813  * \return true if shower complete fit in history; else false
814  */
815  inline bool has_shower_completefit(const Trk& track)
816  {
818  }
819 
820  /**
821  * Test whether given track has default shower fit in history.
822  *
823  * \param track track
824  * \return true if shower fit in history; else false
825  */
826  inline bool has_shower_fit(const Trk& track)
827  {
829  }
830 
831  /**
832  * Test whether given track has default shower fit in history.
833  *
834  * \param track track
835  * \return true if shower fit in history; else false
836  */
837  inline bool has_aashower_fit(const Trk& track)
838  {
840  }
841 
842  /**
843  * Test whether given event has a track according selection.\n
844  * The track selector corresponds to the function operator <tt>bool selector(const Trk&);</tt>.
845  *
846  * \param evt event
847  * \param selector track selector
848  * \return true if at least one corresponding track; else false
849  */
850  template<class JTrackSelector_t>
851  inline bool has_reconstructed_track(const Evt& evt, JTrackSelector_t selector)
852  {
853  return std::find_if(evt.trks.begin(), evt.trks.end(), selector) != evt.trks.end();
854  }
855 
856  /**
857  * Test whether given event has a track according selection.
858  *
859  * \param evt event
860  * \param range range of application types
861  * \return true if at least one corresponding track; else false
862  */
863  template<int reconstruction_type>
864  inline bool has_reconstructed_track(const Evt& evt, const JRange<int> range = JRange<int>())
865  {
866  return has_reconstructed_track(evt, has_history(reconstruction_type, range));
867  }
868 
869  /**
870  * Test whether given event has a track with muon reconstruction.
871  *
872  * \param evt event
873  * \return true if at least one reconstructed muon; else false
874  */
875  inline bool has_reconstructed_muon(const Evt& evt)
876  {
878  }
879 
880  /**
881  * Test whether given event has a track with shower reconstruction.
882  *
883  * \param evt event
884  * \return true if at least one reconstructed shower; else false
885  */
886  inline bool has_reconstructed_shower(const Evt& evt)
887  {
889  }
890 
891  /**
892  * Test whether given event has a track with aashower reconstruction.
893  *
894  * \param evt event
895  * \return true if at least one reconstructed shower; else false
896  */
897  inline bool has_reconstructed_aashower(const Evt& evt)
898  {
900  }
901 
902  /**
903  * Get best reconstructed track.\n
904  * The track selector corresponds to the function operator <tt>bool selector(const Trk&);</tt> and
905  * the track comparator to <tt>bool comparator(const Trk&, const Trk&);</tt>.
906  *
907  * \param evt event
908  * \param selector track selector
909  * \param comparator track comparator
910  * \return track
911  */
912  template<class JTrackSelector_t, class JQualitySorter_t>
913  inline const Trk& get_best_reconstructed_track(const Evt& evt,
914  JTrackSelector_t selector,
915  JQualitySorter_t comparator)
916  {
917  std::vector<Trk>::const_iterator p = std::find_if(evt.trks.begin(), evt.trks.end(), selector);
918 
919  for (std::vector<Trk>::const_iterator i = p; i != evt.trks.end(); ++i) {
920  if (selector(*i) && comparator(*i, *p)) {
921  p = i;
922  }
923  }
924 
925  if (p != evt.trks.end())
926  return *p;
927  else
928  THROW(JIndexOutOfRange, "This event has no reconstructed track with given selector.");
929  }
930 
931  /**
932  * Get best reconstructed track.
933  *
934  * \param evt event
935  * \param range range of application types
936  * \return track
937  */
938  template<int reconstruction_type>
939  inline const Trk& get_best_reconstructed_track(const Evt& evt, const JRange<int> range = JRange<int>())
940  {
941  return get_best_reconstructed_track(evt, has_history(reconstruction_type, range), quality_sorter<reconstruction_type>());
942  }
943 
944 
945  /**
946  * Get best reconstructed muon.
947  *
948  * \param evt event
949  * \return track
950  */
951  inline const Trk& get_best_reconstructed_muon(const Evt& evt)
952  {
954  }
955 
956  /**
957  * Get best reconstructed shower.
958  *
959  * \param evt event
960  * \return track
961  */
962  inline const Trk& get_best_reconstructed_shower(const Evt& evt)
963  {
965  }
966 
967  /**
968  * Get best reconstructed aashower.
969  *
970  * \param evt event
971  * \return track
972  */
973  inline const Trk& get_best_reconstructed_aashower(const Evt& evt)
974  {
976  }
977 }
978 
979 #endif
JException.hh
JAANET::has_shower_positionfit
bool has_shower_positionfit(const Trk &track)
Test whether given track has shower position fit in history.
Definition: JAAnetToolkit.hh:804
JAANET::has_muon_energy
bool has_muon_energy(const Trk &track)
Test whether given track has muon energy fit in history.
Definition: JAAnetToolkit.hh:760
JFIT::JMUONSIMPLEX
JSimplex.cc.
Definition: JFitApplications.hh:24
JAANET::get_electron
const Trk & get_electron(const Evt &evt)
Get first electron from the event tracklist.
Definition: JAAnetToolkit.hh:474
JDETECTOR::getToT
double getToT(const T &tot, const JCalibration &cal)
Get calibrated time-over-threshold of hit.
Definition: JDetector/JCalibration.hh:230
JAANET::HIT_TYPE_UNKNOWN
Unknown source.
Definition: JAAnetToolkit.hh:77
JFIT::JMUONENERGY
JEnergy.cc.
Definition: JFitApplications.hh:26
JAANET::GEANT4_TYPE_ANTIELECTRON
Definition: JParticleTypes.hh:19
JLANG::JIndexOutOfRange
Exception for accessing an index in a collection that is outside of its range.
Definition: JException.hh:90
JTransformation3D.hh
JFitApplications.hh
JVertex3D.hh
JGEOMETRY3D::JTransformation3D
Transformation.
Definition: JTransformation3D.hh:24
JFIT::JMUONGANDALF
JGandalf.cc.
Definition: JFitApplications.hh:25
JTrack3D.hh
JAANET::is_hadron
bool is_hadron(const Trk &track)
Test whether given track is a hadron.
Definition: JAAnetToolkit.hh:407
JGEOMETRY3D::JVertex3D
3D vertex.
Definition: JVertex3D.hh:31
JAANET::from_electron
bool from_electron(const Hit &hit)
Test whether given hit was produced by an electron.
Definition: JAAnetToolkit.hh:491
JAANET::getAxis
JAxis3D getAxis(const Trk &track)
Get axis.
Definition: JAAnetToolkit.hh:244
JPosition3D.hh
JDirection3D.hh
JAANET::has_tau
bool has_tau(const Evt &evt)
Test whether given event has a tau.
Definition: JAAnetToolkit.hh:563
JAANET::has_muon_fit
bool has_muon_fit(const Trk &track)
Test whether given track has default muon fit in history.
Definition: JAAnetToolkit.hh:782
JAANET::TRACK_TYPE_PROTON
Definition: JParticleTypes.hh:77
JAANET::get_tau
const Trk & get_tau(const Evt &evt)
Get first tau from the event tracklist.
Definition: JAAnetToolkit.hh:574
JAANET::GEANT4_TYPE_ANTITAU
Definition: JParticleTypes.hh:58
JAANET::is_photon
bool is_photon(const Trk &track)
Test whether given track is a photon or neutral pion.
Definition: JAAnetToolkit.hh:340
JAANET::is_neutron
bool is_neutron(const Trk &track)
Test whether given track is a (anti-)neutron.
Definition: JAAnetToolkit.hh:399
JAANET::AASHOWER_RECONSTRUCTION_TYPE
static const int AASHOWER_RECONSTRUCTION_TYPE
AAshower reconstruction type for AAnet.
Definition: JAAnetToolkit.hh:62
JDETECTOR::JModuleRouter::getModule
const JModule & getModule(const JObjectID &id) const
Get module parameters.
Definition: JModuleRouter.hh:89
JGEOMETRY3D::JAxis3D
Axis object.
Definition: JAxis3D.hh:38
JAANET::getHit
JHitL0 getHit(const Hit &hit, const JPMTRouter &router)
Get transformation.
Definition: JAAnetToolkit.hh:324
JFIT::JSHOWERPOSITIONFIT
JShowerPositionFit.cc.
Definition: JFitApplications.hh:33
JAANET::has_muon_gandalf
bool has_muon_gandalf(const Trk &track)
Test whether given track has muon gandalf fit in history.
Definition: JAAnetToolkit.hh:749
JAANET::is_tau
bool is_tau(const Trk &track)
Test whether given track is a (anti-)tau.
Definition: JAAnetToolkit.hh:375
JAANET::TRACK_TYPE_NUMU
Definition: JParticleTypes.hh:67
JFIT::JMUONSTART
JStart.cc.
Definition: JFitApplications.hh:27
JAANET::is_proton
bool is_proton(const Trk &track)
Test whether given track is a (anti-)proton.
Definition: JAAnetToolkit.hh:391
JAANET::TRACK_TYPE_NEUTRAL_PION
Definition: JParticleTypes.hh:71
JAANET::getDirection
JDirection3D getDirection(const Trk &track)
Get direction.
Definition: JAAnetToolkit.hh:232
std::vector
Definition: JSTDTypes.hh:12
JAANET::HIT_TYPE_DELTARAYS_DIRECT
Direct light from delta-rays.
Definition: JAAnetToolkit.hh:70
JTOOLS::JRange
Range of values.
Definition: JRange.hh:34
JFIT::JSHOWERPREFIT
JShowerPrefit.cc.
Definition: JFitApplications.hh:32
JGEOMETRY3D::JDirection3D
Data structure for direction in three dimensions.
Definition: JDirection3D.hh:32
JAANET::getNPE
double getNPE(const Hit &hit)
Get true charge of hit.
Definition: JAAnetToolkit.hh:104
JTOOLS::JRange< double >::DEFAULT_RANGE
static const JRange< double, std::less< double > > DEFAULT_RANGE
Default range.
Definition: JRange.hh:558
JAANET::getTrack
JTrack3E getTrack(const Trk &track)
Get track.
Definition: JAAnetToolkit.hh:256
JAANET::has_reconstructed_aashower
bool has_reconstructed_aashower(const Evt &evt)
Test whether given event has a track with aashower reconstruction.
Definition: JAAnetToolkit.hh:897
JAANET::getTransformation
JTransformation3D getTransformation(const Trk &track)
Get transformation.
Definition: JAAnetToolkit.hh:280
JDETECTOR::JPMTRouter
Router for direct addressing of PMT data in detector data structure.
Definition: JPMTRouter.hh:33
JAANET::HIT_TYPE_NOISE
Random noise.
Definition: JAAnetToolkit.hh:76
JAANET::from_muon
bool from_muon(const Hit &hit)
Test whether given hit was produced by a muon.
Definition: JAAnetToolkit.hh:541
JAANET::has_history::operator()
bool operator()(const Trk &track) const
Definition: JAAnetToolkit.hh:709
JAANET::is_pion
bool is_pion(const Trk &track)
Test whether given track is a charged pion.
Definition: JAAnetToolkit.hh:383
JTimeRange.hh
JAANET::is_muon
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
Definition: JAAnetToolkit.hh:367
JAANET::TRACK_TYPE_NUTAU
Definition: JParticleTypes.hh:69
JAANET::HIT_TYPE_MUON_SCATTERED
Scattered light from muon.
Definition: JAAnetToolkit.hh:69
JAANET::HIT_TYPE_SHOWER_DIRECT
Direct light from primary shower.
Definition: JAAnetToolkit.hh:74
JFIT::JSHOWERBEGIN
Start of shower fit applications.
Definition: JFitApplications.hh:31
JTOOLS::JTimeRange
JRange< double > JTimeRange
Type definition for time range.
Definition: JTools/JTimeRange.hh:19
JAANET::TRACK_TYPE_ELECTRON
Definition: JParticleTypes.hh:64
JAxis3D.hh
JPP
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Definition: JAAnetToolkit.hh:37
JAANET::GEANT4_TYPE_MUON
Definition: JParticleTypes.hh:23
JAANET::has_shower_prefit
bool has_shower_prefit(const Trk &track)
Test whether given track has shower prefit in history.
Definition: JAAnetToolkit.hh:793
JRange.hh
JAANET::has_history::has_history
has_history(const int type, const JRange< int > range)
Constructor.
Definition: JAAnetToolkit.hh:700
JAANET::HIT_TYPE_MUON_DIRECT
Direct light from muon.
Definition: JAAnetToolkit.hh:68
JAANET::has_history::type
const int type
reconstruction type
Definition: JAAnetToolkit.hh:717
JAANET::has_aashower_fit
bool has_aashower_fit(const Trk &track)
Test whether given track has default shower fit in history.
Definition: JAAnetToolkit.hh:837
JDETECTOR::JModule::getPMT
const JPMT & getPMT(const int index) const
Get PMT.
Definition: JModule.hh:173
JAANET::count_taus
int count_taus(const Evt &evt)
Count the number of taus in a given event.
Definition: JAAnetToolkit.hh:552
JAANET::has_history::range
const JRange< int > range
range of application types
Definition: JAAnetToolkit.hh:718
JFIT::JMUONPREFIT
JPrefit.cc.
Definition: JFitApplications.hh:23
JAANET::get_best_reconstructed_track
const Trk & get_best_reconstructed_track(const Evt &evt, const JRange< int > range=JRange< int >())
Get best reconstructed track.
Definition: JAAnetToolkit.hh:939
JAANET::has_electron
bool has_electron(const Evt &evt)
Test whether given event has an electron.
Definition: JAAnetToolkit.hh:463
JParticleTypes.hh
JAANET::from_hadron
bool from_hadron(const Hit &hit)
Test whether given hit was produced by a hadronic shower.
Definition: JAAnetToolkit.hh:616
JGEOMETRY3D::JPosition3D
Data structure for position in three dimensions.
Definition: JPosition3D.hh:35
JAANET::get_neutrino
const Trk & get_neutrino(const Evt &evt)
Get incoming neutrino.
Definition: JAAnetToolkit.hh:438
JAANET::get_best_reconstructed_shower
const Trk & get_best_reconstructed_shower(const Evt &evt)
Get best reconstructed shower.
Definition: JAAnetToolkit.hh:962
JAANET::TRACK_TYPE_TAU
Definition: JParticleTypes.hh:68
JAANET::get_best_reconstructed_muon
const Trk & get_best_reconstructed_muon(const Evt &evt)
Get best reconstructed muon.
Definition: JAAnetToolkit.hh:951
JAANET::has_neutrino
bool has_neutrino(const Evt &evt)
Test whether given event has an incoming neutrino.
Definition: JAAnetToolkit.hh:427
JAANET::HIT_TYPE_DELTARAYS_SCATTERED
Scattered light from delta-rays.
Definition: JAAnetToolkit.hh:71
JTrack3E.hh
THROW
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:669
JAANET::has_history
Auxiliary class to test whether given track has specified history.
Definition: JAAnetToolkit.hh:693
JAANET::getTime
double getTime(const Hit &hit)
Get true time of hit.
Definition: JAAnetToolkit.hh:87
JModuleRouter.hh
JAANET::from_tau
bool from_tau(const Hit &hit)
Test whether given hit was produced by a tau.
Definition: JAAnetToolkit.hh:591
JGEOMETRY3D::JTrack3E
3D track with energy.
Definition: JTrack3E.hh:29
JFIT::JPP_RECONSTRUCTION_TYPE
static const int JPP_RECONSTRUCTION_TYPE
Jpp reconstruction type for AAnet.
Definition: JFitApplications.hh:18
JAANET::TRACK_TYPE_CHARGED_PION_PLUS
Definition: JParticleTypes.hh:72
JAANET::has_reconstructed_track
bool has_reconstructed_track(const Evt &evt, const JRange< int > range=JRange< int >())
Test whether given event has a track according selection.
Definition: JAAnetToolkit.hh:864
JAANET::has_shower_completefit
bool has_shower_completefit(const Trk &track)
Test whether given track has shower complete fit in history.
Definition: JAAnetToolkit.hh:815
JAANET::get_best_reconstructed_aashower
const Trk & get_best_reconstructed_aashower(const Evt &evt)
Get best reconstructed aashower.
Definition: JAAnetToolkit.hh:973
JAANET::JHitType_t
JHitType_t
Enumeration of hit types based on km3 codes.
Definition: JAAnetToolkit.hh:67
JAANET::TRACK_TYPE_NUE
Definition: JParticleTypes.hh:65
JAANET::has_reconstructed_shower
bool has_reconstructed_shower(const Evt &evt)
Test whether given event has a track with shower reconstruction.
Definition: JAAnetToolkit.hh:886
JFIT::JMUONEND
End of muon fit applications.
Definition: JFitApplications.hh:29
JAANET::quality_sorter::operator()
bool operator()(const Trk &first, const Trk &second) const
The default comparison is based on:
Definition: JAAnetToolkit.hh:680
JTOOLS::v
data_type v[N+1][M+1]
Definition: JPolint.hh:707
JAANET
Extensions to AAnet data format.
Definition: JAAnetToolkit.hh:36
JAANET::get_hadronic_cascade
Trk get_hadronic_cascade(const Evt &evt)
Auxiliary function to get average direction and total energy of a hadronic cascade.
Definition: JAAnetToolkit.hh:627
JDETECTOR::JModuleRouter
Router for direct addressing of module data in detector data structure.
Definition: JModuleRouter.hh:34
JAANET::count_electrons
int count_electrons(const Evt &evt)
Count the number of electrons in a given event.
Definition: JAAnetToolkit.hh:452
JAANET::getVertex
JVertex3D getVertex(const Trk &track)
Get vertex.
Definition: JAAnetToolkit.hh:268
JAANET::has_particleID
bool has_particleID(const Trk &track, const int type)
Test whether given track corresponds to given particle type.
Definition: JAAnetToolkit.hh:419
JAANET::has_muon_start
bool has_muon_start(const Trk &track)
Test whether given track has muon start fit in history.
Definition: JAAnetToolkit.hh:771
JFIT::has_history
bool has_history(const JFit &fit, const int type)
Test whether given fit has specified history.
Definition: JEvtToolkit.hh:262
JAANET::has_shower_fit
bool has_shower_fit(const Trk &track)
Test whether given track has default shower fit in history.
Definition: JAAnetToolkit.hh:826
JAANET::has_reconstructed_muon
bool has_reconstructed_muon(const Evt &evt)
Test whether given event has a track with muon reconstruction.
Definition: JAAnetToolkit.hh:875
JAANET::is_neutrino
bool is_neutrino(const Trk &track)
Test whether given track is a neutrino.
Definition: JAAnetToolkit.hh:349
JTRIGGER::JHitL0
Data structure for L0 hit.
Definition: JHitL0.hh:27
JAANET::get_muon
const Trk & get_muon(const Evt &evt)
Get first muon from the event tracklist.
Definition: JAAnetToolkit.hh:524
JPMTRouter.hh
JAANET::count_muons
int count_muons(const Evt &evt)
Count the number of muons in a given event.
Definition: JAAnetToolkit.hh:502
std
Definition: jaanetDictionary.h:36
JAANET::is_electron
bool is_electron(const Trk &track)
Test whether given track is a (anti-)electron.
Definition: JAAnetToolkit.hh:359
JAANET::GEANT4_TYPE_ANTIMUON
Definition: JParticleTypes.hh:22
JGEOMETRY3D::JTrack3D
3D track.
Definition: JTrack3D.hh:30
JFIT::JSHOWERCOMPLETEFIT
JShowerFit.cc.
Definition: JFitApplications.hh:34
JAANET::has_muon_prefit
bool has_muon_prefit(const Trk &track)
Test whether given track has muon prefit in history.
Definition: JAAnetToolkit.hh:727
JAANET::HIT_TYPE_BREMSSTRAHLUNG_DIRECT
Direct light from Bremsstrahlung.
Definition: JAAnetToolkit.hh:72
JAANET::count_hadrons
int count_hadrons(const Evt &evt)
Count the number of hadrons in a given event (not including neutral pions)
Definition: JAAnetToolkit.hh:602
JAANET::GEANT4_TYPE_ELECTRON
Definition: JParticleTypes.hh:20
JAANET::add_time_offset
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...
Definition: JAAnetToolkit.hh:653
JAANET::TRACK_TYPE_MUON
Definition: JParticleTypes.hh:66
JAANET::TRACK_TYPE_PHOTON
Definition: JParticleTypes.hh:70
JFIT::JMUONBEGIN
Start of muon fit applications.
Definition: JFitApplications.hh:22
JAANET::HIT_TYPE_SHOWER_SCATTERED
Scattered light from primary shower.
Definition: JAAnetToolkit.hh:75
JAANET::has_muon_simplex
bool has_muon_simplex(const Trk &track)
Test whether given track has muon simplex fit in history.
Definition: JAAnetToolkit.hh:738
JAANET::HIT_TYPE_BREMSSTRAHLUNG_SCATTERED
Scattered light from Bremsstrahlung.
Definition: JAAnetToolkit.hh:73
JFIT::JSHOWEREND
End of shower fit applications.
Definition: JFitApplications.hh:36
JHitL0.hh
JAANET::getTimeRange
JTimeRange getTimeRange(const Evt &event, const JTimeRange &T_ns)
Get time range (i.e.
Definition: JAAnetToolkit.hh:158
JAANET::GEANT4_TYPE_TAU
Definition: JParticleTypes.hh:59
JAANET::TRACK_TYPE_NEUTRON
Definition: JParticleTypes.hh:78
JDETECTOR::JPMTRouter::getPMT
const JPMT & getPMT(const JPMTAddress &address) const
Get PMT.
Definition: JPMTRouter.hh:90
JAANET::has_muon
bool has_muon(const Evt &evt)
Test whether given event has a muon.
Definition: JAAnetToolkit.hh:513
JAANET::getPosition
JPosition3D getPosition(const Trk &track)
Get position.
Definition: JAAnetToolkit.hh:209
JAANET::quality_sorter
Reconstruction type dependent comparison of track quality.
Definition: JAAnetToolkit.hh:670
JAANET::is_noise
bool is_noise(const Hit &hit)
Verify hit origin.
Definition: JAAnetToolkit.hh:121