Jpp 19.3.0-rc.3
the software that should make you happy
Loading...
Searching...
No Matches
JRemovePMT.cc
Go to the documentation of this file.
1
2#include <string>
3#include <iostream>
4#include <iomanip>
5#include <set>
6
11#include "JDAQ/JDAQEventIO.hh"
13
16#include "JSupport/JSupport.hh"
17
21
22#include "Jeep/JParser.hh"
23#include "Jeep/JMessage.hh"
24
25namespace {
26
28
29 /**
30 * Remove PMT(s) from data.
31 *
32 * \param data data (I/O)
33 * \param PMT PMT(s) to be removed
34 */
35 template<class T>
36 inline void remove(std::vector<T>& data, const std::set<JDAQPMTIdentifier>& PMT)
37 {
38 using namespace std;
39
40 vector<T> buffer;
41
42 for (typename vector<T>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
43 if (PMT.find(JDAQPMTIdentifier(hit->getModuleID(), hit->getPMT())) == PMT.end() &&
44 PMT.find(JDAQPMTIdentifier(-1, hit->getPMT())) == PMT.end()) {
45 buffer.push_back(*hit);
46 }
47 }
48
49 data.swap(buffer);
50 }
51}
52
53
54/**
55 * \file
56 *
57 * Example program to remove PMT(s) from data (and set corresponding rate to 0).
58 * \author mdejong
59 */
60int main(int argc, char **argv)
61{
62 using namespace std;
63 using namespace JPP;
64 using namespace KM3NETDAQ;
65
67 JLimit_t& numberOfEvents = inputFile.getLimit();
70 int debug;
71
72 try {
73
74 JParser<> zap("Example program to remove PMT(s) from data (and set corresponding rate to 0).");
75
76 zap['f'] = make_field(inputFile);
77 zap['o'] = make_field(outputFile) = "abc.root";
78 zap['n'] = make_field(numberOfEvents) = JLimit::max();
79 zap['P'] = make_field(PMT) = JPARSER::initialised();
80 zap['d'] = make_field(debug) = 1;
81
82 zap(argc, argv);
83 }
84 catch(const exception& error) {
85 FATAL(error.what() << endl);
86 }
87
88
89 outputFile.open();
90
91 {
93
94 // put base class JDAQTimeslice at end of type list
95
97
99
100 for (counter_type counter = 0; in.hasNext() && counter != inputFile.getLimit(); ++counter) {
101
102 STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
103
104 JDAQTimeslice* timeslice = in.next();
105
106 for (JDAQTimeslice::iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
107
108 bool rm = false;
109
110 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
111
112 if (PMT.find(JDAQPMTIdentifier(frame->getModuleID(), pmt)) != PMT.end() ||
113 PMT.find(JDAQPMTIdentifier(-1, pmt)) != PMT.end()) {
114 rm = true;
115 }
116 }
117
118 if (rm) {
119
121
122 for (JDAQSuperFrame::const_iterator hit = frame->begin(); hit != frame->end(); ++hit) {
123 buffer[hit->getPMT()].push_back(*hit);
124 }
125
126 vector<JDAQHit> data;
127
128 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
129
130 if (PMT.find(JDAQPMTIdentifier(frame->getModuleID(), pmt)) == PMT.end() &&
131 PMT.find(JDAQPMTIdentifier(-1, pmt)) == PMT.end()) {
132 copy(buffer[pmt].begin(), buffer[pmt].end(), back_inserter(data));
133 }
134 }
135
136 sort(data.begin(), data.end());
137
138 frame->clear();
139 frame->add(data.size(), data.data());
140 }
141 }
142
143 out.put(*timeslice);
144 }
145 STATUS(endl);
146 }
147 {
149
150 while (in.hasNext()) {
151
152 STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
153
154 JDAQEvent* event = in.next();
155
156 remove(event->getHits<JDAQTriggeredHit>(), PMT);
157 remove(event->getHits<JDAQSnapshotHit> (), PMT);
158
159 outputFile.put(*event);
160 }
161 STATUS(endl);
162 }
163 {
165
166 while (in.hasNext()) {
167
168 STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
169
170 JDAQSummaryslice* summaryslice = in.next();
171
172 for (JDAQSummaryslice::iterator frame = summaryslice->begin(); frame != summaryslice->end(); ++frame) {
173
174 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
175
176 if (PMT.find(JDAQPMTIdentifier(frame->getModuleID(), pmt)) != PMT.end() ||
177 PMT.find(JDAQPMTIdentifier(-1, pmt)) != PMT.end()) {
178 (*frame)[pmt].setValue(0.0);
179 }
180 }
181 }
182
183 outputFile.put(*summaryslice);
184 }
185 STATUS(endl);
186 }
187
189
190 io >> outputFile;
191
192 outputFile.close();
193}
JDAQPMTIdentifier PMT
Command line options.
string outputFile
KM3NeT DAQ constants, bit handling, etc.
Recording of objects on file according a format that follows from the file name extension.
General purpose messaging.
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define STATUS(A)
Definition JMessage.hh:63
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
Scanning of objects from multiple files according a format that follows from the extension of each fi...
Utility class to parse command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
int main(int argc, char **argv)
Definition JRemovePMT.cc:60
ROOT TTree parameter settings of various packages.
Auxiliary class for multiplexing object iterators.
virtual bool hasNext() override
Check availability of next element.
virtual const pointer_type & next() override
Get next element.
Utility class to parse command line options.
Definition JParser.hh:1698
Object writing to file.
General purpose class for object reading from a list of file names.
virtual bool hasNext() override
Check availability of next element.
counter_type getCounter() const
Get counter.
virtual const pointer_type & next() override
Get next element.
Hit data structure.
Definition JDAQHit.hh:35
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition JHead.cc:163
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Long64_t counter_type
Type definition for counter.
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition JDAQ.hh:26
Append to type list.
Definition JTypeList.hh:62
Auxiliary class for demultiplexing object outputs.
virtual bool put(const JBase_t &object) override
Object output.
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition JParser.hh:68
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