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