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