Jpp in_tag_pdf_generation
the software that should make you happy
Loading...
Searching...
No Matches
JDAQPrint.cc File Reference

Example program to print DAQ events with associated information. More...

#include <string>
#include <iostream>
#include <iomanip>
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JSingleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JSupport/JSummaryFileRouter.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Example program to print DAQ events with associated information.

Author
mdejong

Definition in file JDAQPrint.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 29 of file JDAQPrint.cc.

30{
31 using namespace std;
32 using namespace JPP;
33 using namespace KM3NETDAQ;
34
35 JSingleFileScanner<> inputFile;
36 JLimit_t& numberOfEvents = inputFile.getLimit();
37 string detectorFile;
38 int debug;
39
40 try {
41
42 JParser<> zap("Example program to print DAQ events with associated information.");
43
44 zap['f'] = make_field(inputFile);
45 zap['a'] = make_field(detectorFile) = "";
46 zap['n'] = make_field(numberOfEvents) = JLimit::max();
47 zap['d'] = make_field(debug) = 1;
48
49 zap(argc, argv);
50 }
51 catch(const exception& error) {
52 FATAL(error.what() << endl);
53 }
54
55
57
58 if (detectorFile != "") {
59
60 try {
61 load(detectorFile, detector);
62 }
63 catch(const JException& error) {
64 FATAL(error);
65 }
66 }
67
68 const JModuleRouter router(detector);
69
70
71 JSummaryFileRouter summary(inputFile);
72
73 for (JTreeScanner<JDAQEvent, JDAQEvaluator> in(inputFile, inputFile.getLimit()); in.hasNext(); ) {
74
75 const JDAQEvent* tev = in.next();
76
77 summary.update(*tev);
78
79 cout << "Detector identifier " << tev->getDetectorID() << endl;
80 cout << "Run number " << tev->getRunNumber() << endl;
81 cout << "Frame index " << tev->getFrameIndex() << endl;
82 cout << "UTC " << tev->getTimesliceStart() << endl;
83 cout << "Trigger counter " << tev->getCounter() << endl;
84 cout << "Trigger mask (hex) " << hex << tev->getTriggerMask() << dec << endl;
85 cout << "Overlays " << tev->getOverlays() << endl;
86 cout << "Triggered hits " << tev->size<JDAQTriggeredHit>() << endl;
87 cout << "Snapshot hits " << tev->size<JDAQSnapshotHit> () << endl;
88 cout << endl;
89
90 cout << " module PMT TDC ToT trigger rate [kHz] serial x y z dx dy dz t0 status" << endl;
91 {
93
94 for (JDAQEvent::const_iterator<JHit_t> hit = tev->begin<JHit_t>(); hit != tev->end<JHit_t>(); ++hit) {
95
96 cout << setw(8) << hit->getModuleID() << ' '
97 << setw(2) << (int) hit->getPMT() << ' '
98 << setw(10) << (int) hit->getT() << ' '
99 << setw(3) << (int) hit->getToT() << ' '
100 << setw(10) << hex << hit->getTriggerMask() << dec;
101
102 if (summary.hasSummaryFrame(hit->getModuleID()))
103 cout << ' ' << FIXED(9,2) << summary.getRate(JDAQPMTIdentifier(hit->getModuleID(), hit->getPMT())) * 1.0e-3;
104 else
105 cout << ' ' << setw(9) << " ? ";
106
107 if (router.hasModule(hit->getModuleID())) {
108 cout << ' ' << router.getPMT(JPMTIdentifier(hit->getModuleID(), hit->getPMT()));
109 }
110
111 cout << endl;
112 }
113 }
114 cout << endl;
115
116 if (debug >= debug_t) {
117
118 typedef JDAQSnapshotHit JHit_t;
119
120 for (JDAQEvent::const_iterator<JHit_t> hit = tev->begin<JHit_t>(); hit != tev->end<JHit_t>(); ++hit) {
121
122 cout << setw(8) << hit->getModuleID() << ' '
123 << setw(2) << (int) hit->getPMT() << ' '
124 << setw(10) << (int) hit->getT() << ' '
125 << setw(3) << (int) hit->getToT() << ' '
126 << setw(10) << ' ';
127
128 if (summary.hasSummaryFrame(hit->getModuleID()))
129 cout << ' ' << FIXED(9,2) << summary.getRate(JDAQPMTIdentifier(hit->getModuleID(), hit->getPMT())) * 1.0e-3;
130 else
131 cout << ' ' << setw(9) << " ? ";
132
133 if (router.hasModule(hit->getModuleID())) {
134 cout << ' ' << router.getPMT(JPMTIdentifier(hit->getModuleID(), hit->getPMT()));
135 }
136
137 cout << endl;
138 }
139 }
140 cout << endl;
141 }
142}
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Detector data structure.
Definition JDetector.hh:96
Router for direct addressing of module data in detector data structure.
General exception.
Definition JException.hh:24
Utility class to parse command line options.
Definition JParser.hh:1698
File router for fast addressing of summary data.
Template definition for direct access of elements in ROOT TChain.
int getDetectorID() const
Get detector identifier.
int getRunNumber() const
Get run number.
JDAQUTCExtended getTimesliceStart() const
Get start of timeslice.
int getFrameIndex() const
Get frame index.
unsigned int getOverlays() const
Get number of overlays.
Template const_iterator.
Definition JDAQEvent.hh:68
static JTriggerMask_t getTriggerMask(const JDAQTriggeredHit &hit)
Get trigger mask of given hit.
Definition JDAQEvent.hh:226
const_iterator< T > end() const
Get end of data.
unsigned int size() const
Get number of hits.
const_iterator< T > begin() const
Get begin of data.
JTriggerCounter_t getCounter() const
Get trigger counter.
JTriggerCounter_t next()
Increment trigger counter.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
@ debug_t
debug
Definition JMessage.hh:29
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
Auxiliary data structure for floating point format specification.
Definition JManip.hh:448
Detector file.
Definition JHead.hh:227
Auxiliary class to set-up Hit.
Definition JSirene.hh:58
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45
static counter_type max()
Get maximum counter value.
Definition JLimit.hh:128