Jpp 20.0.0-195-g190c9e876
the software that should make you happy
Loading...
Searching...
No Matches
JHitsExtractor.cc
Go to the documentation of this file.
1#include <string>
2#include <iostream>
3#include <iomanip>
4#include <vector>
5#include <set>
6#include <algorithm>
7#include <limits>
8#include <map>
9#include <memory>
10
11#include "TFile.h"
12#include "TTree.h"
13
14#include "JDAQ/JDAQEventIO.hh"
17
18#include "JTrigger/JHitL0.hh"
19#include "JTrigger/JBuildL0.hh"
20
25
27
31
33
34#include "Jeep/JParser.hh"
35
36
37namespace {
38
39 /**
40 * Auxiliary data structure for sorting of hits.
41 */
42 static const struct {
43 /**
44 * Compare hits by PMT identifier and time.
45 *
46 * \param first first hit
47 * \param second second hit
48 * \return true if first before second; else false
49 */
50 template<class T>
51 bool operator()(const T& first, const T& second) const
52 {
53 if (first.getPMTIdentifier() == second.getPMTIdentifier())
54 return first.getT() < second.getT();
55 else
56 return first.getPMTIdentifier() < second.getPMTIdentifier();
57 }
58 } compare;
59}
60
61
62
63/**
64 * \file
65 *
66 * Program to extract information about hits and all PMTs in the detector for water property studies using I/O of JFIT::JEvt data.
67 */
68int main(int argc, char **argv)
69{
70
71 using namespace std;
72 using namespace JPP;
73 using namespace KM3NETDAQ;
74
75 typedef JMultipleFileScanner<calibration_types> JCalibration_t;
77 typedef JSingleFileScanner_t::pointer_type pointer_type;
78
79
81 string outputFile;
82 JLimit_t& numberOfEvents = inputFile.getLimit();
83 string detectorFile;
84 JCalibration_t calibrationFile;
85 double Tmax_s;
86 double rate_Hz;
87 JPMTParametersMap pmtParameters;
88 int debug;
89
90 try {
91
92 JParser<> zap("Program to extract information about hits and all PMTs in the detector for water property studies");
93
94 zap['f'] = make_field(inputFile);
95 zap['o'] = make_field(outputFile) = "hits_pmts.root";
96 zap['n'] = make_field(numberOfEvents) = JLimit::max();
97 zap['a'] = make_field(detectorFile);
98 zap['+'] = make_field(calibrationFile) = JPARSER::initialised();
99 zap['T'] = make_field(Tmax_s) = 100;
100 zap['P'] = make_field(pmtParameters) = JPARSER::initialised();
101 zap['B'] = make_field(rate_Hz) = 6.0e3;
102 zap['d'] = make_field(debug) = 1;
103
104 zap(argc, argv);
105 }
106 catch(const exception& error) {
107 FATAL(error.what() << endl);
108 }
109
110
112
113 try {
114 load(detectorFile, detector);
115 }
116 catch(const JException& error) {
117 FATAL(error);
118 }
119
120 unique_ptr<JDynamics> dynamics;
121
122 try {
123 dynamics.reset(new JDynamics(detector, Tmax_s));
124 dynamics->load(calibrationFile);
125 }
126 catch(const exception& error) {
127 if (!calibrationFile.empty()) {
128 FATAL(error.what());
129 }
130 }
131
132 TFile *out;
133 try {
134 out = new TFile(outputFile.c_str(), "RECREATE");
135 }
136 catch(const exception& error){
137 FATAL(error.what() << endl);
138 }
139
140
141
142
143 // Create tree to store information about hits: first hit or not, PMT QE for this hit, PMT rate for this hit, PMT status for this hit, HRV
144 TTree *hits_tree = new TTree("Hits","Tree with hits");
145 vector<float> hit_rate, hit_QE, hit_pos_x, hit_pos_y, hit_pos_z, hit_dir_x, hit_dir_y, hit_dir_z, hit_time, hit_tot;
146 vector<int> hit_first, hit_HRV, hit_status, hit_dom_id, hit_channel_id;
147 hits_tree->Branch("first",&hit_first);
148 hits_tree->Branch("QE",&hit_QE);
149 hits_tree->Branch("rate",&hit_rate);
150 hits_tree->Branch("status",&hit_status);
151 hits_tree->Branch("HRV",&hit_HRV);
152 hits_tree->Branch("pos_x",&hit_pos_x);
153 hits_tree->Branch("pos_y",&hit_pos_y);
154 hits_tree->Branch("pos_z",&hit_pos_z);
155 hits_tree->Branch("dir_x",&hit_dir_x);
156 hits_tree->Branch("dir_y",&hit_dir_y);
157 hits_tree->Branch("dir_z",&hit_dir_z);
158 hits_tree->Branch("time",&hit_time);
159 hits_tree->Branch("tot",&hit_tot);
160 hits_tree->Branch("channel_id",&hit_channel_id);
161 hits_tree->Branch("dom_id",&hit_dom_id);
162
163 // Create tree to store information about all PMTs of the detector
164 TTree *pmts_tree = new TTree("PMTs","Tree with PMTs");
165 vector<float> pmt_pos_x, pmt_pos_y, pmt_pos_z, pmt_dir_x, pmt_dir_y, pmt_dir_z, pmt_QE, pmt_rate;
166 vector<int> pmt_channel_id, pmt_dom_id, pmt_status, pmt_HRV;
167 pmts_tree->Branch("pos_x",&pmt_pos_x);
168 pmts_tree->Branch("pos_y",&pmt_pos_y);
169 pmts_tree->Branch("pos_z",&pmt_pos_z);
170 pmts_tree->Branch("dir_x",&pmt_dir_x);
171 pmts_tree->Branch("dir_y",&pmt_dir_y);
172 pmts_tree->Branch("dir_z",&pmt_dir_z);
173 pmts_tree->Branch("QE",&pmt_QE);
174 pmts_tree->Branch("channel_id",&pmt_channel_id);
175 pmts_tree->Branch("dom_id",&pmt_dom_id);
176 pmts_tree->Branch("rate",&pmt_rate);
177 pmts_tree->Branch("status",&pmt_status);
178 pmts_tree->Branch("HRV",&pmt_HRV);
179
180
181 const JModuleRouter router(dynamics ? dynamics->getDetector() : detector);
182
183 JSummaryFileRouter summary(inputFile);
184
185 JTreeScanner<JDAQEvent, JDAQEvaluator> in(inputFile, inputFile.getLimit());
186
187 typedef vector<JHitL0> JDataL0_t;
188
189 const JBuildL0<JHitL0> buildL0;
190
191 while (in.hasNext()) {
192
193 STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
194
195 pointer_type ps = in.next();
196
197 JDAQEvent* tev = ps;
198
199 summary.update(*tev);
200
201 if (dynamics) {
202 dynamics->update(*tev);
203 }
204
205 JDataL0_t dataL0;
206
207 buildL0(*tev, router, true, back_inserter(dataL0));
208
209 sort(dataL0.begin(), dataL0.end(), JHitL0::compare);
210
211 int module_id_last = -1;
212 int channel_id_last = -1;
213
214 // Loop over hits
215 for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
216
217 const JPMTIdentifier id(i->getModuleID(), i->getPMTAddress());
218
219 const JPMTParameters& wip = pmtParameters.getPMTParameters(id);
220
221 const int type = wip.getType();
222 const double QE = wip.QE;
223 const double R_Hz = summary.getRate(i->getPMTIdentifier(), rate_Hz);
224
225 JHitW0 hit(*i, type, QE, R_Hz);
226
227 // Positions
228 hit_pos_x.push_back(hit.getX());
229 hit_pos_y.push_back(hit.getY());
230 hit_pos_z.push_back(hit.getZ());
231 // Directions
232 hit_dir_x.push_back(hit.getDX());
233 hit_dir_y.push_back(hit.getDY());
234 hit_dir_z.push_back(hit.getDZ());
235 // Time and ToT
236 hit_time.push_back(hit.getT());
237 hit_tot.push_back(hit.getToT());
238 // DOM and channel IDs
239 hit_dom_id.push_back(hit.getModuleID());
240 hit_channel_id.push_back(hit.getPMTAddress());
241
242
243 const JPMT& pmt = router.getPMT(JPMTIdentifier(hit.getModuleID(), hit.getPMTAddress()));
244
245 hit_status.push_back(pmt.getStatus<JStatus::status_type>()); // PMT status
246 hit_QE.push_back(hit.getQE()); // PMT QE
247
248 const JDAQSummaryFrame& SummaryFrame = summary.getSummaryFrame(hit.getModuleID());
249
250 hit_HRV.push_back(SummaryFrame.testHighRateVeto(hit.getPMTAddress()) ||
251 SummaryFrame.testFIFOStatus(hit.getPMTAddress())); // Flag for PMT in HRV
252
253 hit_rate.push_back(hit.getR()); // PMT rate
254
255 // Check if the hit is the first one on the channel (1 means first one)
256 if ((hit.getModuleID()==module_id_last) && (hit.getPMTAddress()==channel_id_last)){
257 hit_first.push_back(0);
258 }
259 else{
260 hit_first.push_back(1);
261 }
262 module_id_last = hit.getModuleID();
263 channel_id_last = hit.getPMTAddress();
264
265
266 }
267
268
269 // Loop over all PMTs
270 for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
271
272 JModule module = *i;
273
274 for (size_t pmt = 0; pmt != module.size(); ++pmt) {
275
276 pmt_pos_x.push_back(module[pmt].getX());
277 pmt_pos_y.push_back(module[pmt].getY());
278 pmt_pos_z.push_back(module[pmt].getZ());
279
280 pmt_dir_x.push_back(module[pmt].getDX());
281 pmt_dir_y.push_back(module[pmt].getDY());
282 pmt_dir_z.push_back(module[pmt].getDZ());
283
284 pmt_status.push_back(router.getPMT(JPMTIdentifier(module.getID(), pmt)).getStatus<JStatus::status_type>());
285
286 pmt_QE.push_back(pmtParameters.getQE(JPMTIdentifier(module.getID(), pmt)));
287
288 pmt_channel_id.push_back(pmt);
289 pmt_dom_id.push_back(module.getID());
290
291 const JDAQSummaryFrame& SummaryFrame = summary.getSummaryFrame(module.getID());
292 pmt_HRV.push_back(SummaryFrame.testHighRateVeto(pmt));
293
294 pmt_rate.push_back(summary.getRate(JDAQPMTIdentifier(module.getID(), pmt)));
295
296 }
297 }
298
299
300 hits_tree->Fill();
301 hit_first.clear();
302 hit_QE.clear();
303 hit_rate.clear();
304 hit_status.clear();
305 hit_HRV.clear();
306 hit_pos_x.clear();
307 hit_pos_y.clear();
308 hit_pos_z.clear();
309 hit_dir_x.clear();
310 hit_dir_y.clear();
311 hit_dir_z.clear();
312 hit_time.clear();
313 hit_tot.clear();
314 hit_dom_id.clear();
315 hit_channel_id.clear();
316
317
318 pmts_tree->Fill();
319 pmt_pos_x.clear();
320 pmt_pos_y.clear();
321 pmt_pos_z.clear();
322 pmt_dir_x.clear();
323 pmt_dir_y.clear();
324 pmt_dir_z.clear();
325 pmt_channel_id.clear();
326 pmt_dom_id.clear();
327 pmt_QE.clear();
328 pmt_rate.clear();
329 pmt_status.clear();
330 pmt_HRV.clear();
331
332 }
333
334 STATUS(endl);
335 out->Write();
336 out->Close();
337}
338
339
string outputFile
Data structure for detector geometry and calibration.
Dynamic detector calibration.
Basic data structure for L0 hit.
int main(int argc, char **argv)
#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:74
Direct access to module in detector data structure.
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
Scanning of objects from a single file according a format that follows from the extension of each fil...
Detector data structure.
Definition JDetector.hh:96
Router for direct addressing of module data in detector data structure.
const JPMT & getPMT(const JPMTIdentifier &id) const
Get PMT parameters.
Data structure for a composite optical module.
Definition JModule.hh:76
Auxiliary class for map of PMT parameters.
const JPMTParameters & getPMTParameters(const JPMTIdentifier &id) const
Get PMT parameters.
double getQE(const JPMTIdentifier &id) const
Get QE of given PMT.
Data structure for PMT parameters.
double QE
relative quantum efficiency
int getType() const
Get type for for time-slewing correction.
Data structure for PMT geometry, calibration and status.
Definition JPMT.hh:49
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
General exception.
Definition JException.hh:24
int getID() const
Get identifier.
Definition JObjectID.hh:50
Utility class to parse command line options.
Definition JParser.hh:1698
Auxiliary class for a hit with background rate value.
Definition JHitW0.hh:25
double getQE() const
Get relative quantum efficiency.
Definition JHitW0.hh:84
double getT() const
Get calibrated time of hit.
Definition JHitW0.hh:62
double getR() const
Get rate.
Definition JHitW0.hh:95
General purpose class for object reading from a list of file names.
Object reading from a list of files.
File router for fast addressing of summary data.
void update(const JDAQHeader &header)
Update router.
double getRate(const JDAQPMTIdentifier &id) const
Get rate.
const JDAQSummaryFrame & getSummaryFrame(const JDAQModuleIdentifier &module) const
Get summary frame.
Template definition for direct access of elements in ROOT TChain.
Template L0 hit builder.
Definition JBuildL0.hh:38
double getToT() const
Get calibrated time over threshold of hit.
bool testFIFOStatus() const
Test FIFO status.
bool testHighRateVeto() const
Test high-rate veto status.
int getModuleID() const
Get module identifier.
int getPMTAddress() const
Get PMT identifier.
Data storage class for rate measurements of all PMTs in one module.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
Detector file.
Definition JHead.hh:227
Dynamic detector calibration.
Definition JDynamics.hh:81
status_type getStatus(const JType< status_type > &type) const
Get status.
Definition JStatus.hh:60
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
Auxiliary base class for file name.
Auxiliary data structure for sorting of hits.
Definition JHitL0.hh:85