Jpp 20.0.0-rc.2
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
30
32
33#include "Jeep/JParser.hh"
34
35
36namespace {
37
38 /**
39 * Auxiliary data structure for sorting of hits.
40 */
41 static const struct {
42 /**
43 * Compare hits by PMT identifier and time.
44 *
45 * \param first first hit
46 * \param second second hit
47 * \return true if first before second; else false
48 */
49 template<class T>
50 bool operator()(const T& first, const T& second) const
51 {
52 if (first.getPMTIdentifier() == second.getPMTIdentifier())
53 return first.getT() < second.getT();
54 else
55 return first.getPMTIdentifier() < second.getPMTIdentifier();
56 }
57 } compare;
58}
59
60
61
62/**
63 * \file
64 *
65 * Program to extract information about hits and all PMTs in the detector for water property studies using I/O of JFIT::JEvt data.
66 */
67int main(int argc, char **argv)
68{
69
70 using namespace std;
71 using namespace JPP;
72 using namespace KM3NETDAQ;
73
75 JACOUSTICS::JEvt>::typelist calibration_types;
76 typedef JMultipleFileScanner<calibration_types> JCalibration_t;
77
79 typedef JSingleFileScanner_t::pointer_type pointer_type;
80
81
83 string outputFile;
84 JLimit_t& numberOfEvents = inputFile.getLimit();
85 string detectorFile;
86 JCalibration_t calibrationFile;
87 int debug;
88 double Tmax_s;
89 double R_Hz = 6000;
90 JPMTParametersMap pmtParameters;
91
92 try {
93
94 JParser<> zap("Program to extract information about hits and all PMTs in the detector for water property studies");
95
96 zap['f'] = make_field(inputFile);
97 zap['o'] = make_field(outputFile) = "hits_pmts.root";
98 zap['a'] = make_field(detectorFile);
99 zap['+'] = make_field(calibrationFile) = JPARSER::initialised();
100 zap['P'] = make_field(pmtParameters, "PMT simulation data (or corresponding file name)") = JPARSER::initialised();
101 zap['n'] = make_field(numberOfEvents) = JLimit::max();
102 zap['t'] = make_field(Tmax_s) = 100;
103 zap['d'] = make_field(debug) = 1;
104
105 zap(argc, argv);
106 }
107 catch(const exception& error) {
108 FATAL(error.what() << endl);
109 }
110
111
113
114 try {
115 load(detectorFile, detector);
116 }
117 catch(const JException& error) {
118 FATAL(error);
119 }
120
121 unique_ptr<JDynamics> dynamics;
122
123 try {
124 dynamics.reset(new JDynamics(detector, Tmax_s));
125 dynamics->load(calibrationFile);
126 }
127 catch(const exception& error) {
128 if (!calibrationFile.empty()) {
129 FATAL(error.what());
130 }
131 }
132
133 TFile *out;
134 try {
135 out = new TFile(outputFile.c_str(), "RECREATE");
136 }
137 catch(const exception& error){
138 FATAL(error.what() << endl);
139 }
140
141
142
143
144 // 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
145 TTree *hits_tree = new TTree("Hits","Tree with hits");
146 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;
147 vector<int> hit_first, hit_HRV, hit_status, hit_dom_id, hit_channel_id;
148 hits_tree->Branch("first",&hit_first);
149 hits_tree->Branch("QE",&hit_QE);
150 hits_tree->Branch("rate",&hit_rate);
151 hits_tree->Branch("status",&hit_status);
152 hits_tree->Branch("HRV",&hit_HRV);
153 hits_tree->Branch("pos_x",&hit_pos_x);
154 hits_tree->Branch("pos_y",&hit_pos_y);
155 hits_tree->Branch("pos_z",&hit_pos_z);
156 hits_tree->Branch("dir_x",&hit_dir_x);
157 hits_tree->Branch("dir_y",&hit_dir_y);
158 hits_tree->Branch("dir_z",&hit_dir_z);
159 hits_tree->Branch("time",&hit_time);
160 hits_tree->Branch("tot",&hit_tot);
161 hits_tree->Branch("channel_id",&hit_channel_id);
162 hits_tree->Branch("dom_id",&hit_dom_id);
163
164 // Create tree to store information about all PMTs of the detector
165 TTree *pmts_tree = new TTree("PMTs","Tree with PMTs");
166 vector<float> pmt_pos_x, pmt_pos_y, pmt_pos_z, pmt_dir_x, pmt_dir_y, pmt_dir_z, pmt_QE, pmt_rate;
167 vector<int> pmt_channel_id, pmt_dom_id, pmt_status, pmt_HRV;
168 pmts_tree->Branch("pos_x",&pmt_pos_x);
169 pmts_tree->Branch("pos_y",&pmt_pos_y);
170 pmts_tree->Branch("pos_z",&pmt_pos_z);
171 pmts_tree->Branch("dir_x",&pmt_dir_x);
172 pmts_tree->Branch("dir_y",&pmt_dir_y);
173 pmts_tree->Branch("dir_z",&pmt_dir_z);
174 pmts_tree->Branch("QE",&pmt_QE);
175 pmts_tree->Branch("channel_id",&pmt_channel_id);
176 pmts_tree->Branch("dom_id",&pmt_dom_id);
177 pmts_tree->Branch("rate",&pmt_rate);
178 pmts_tree->Branch("status",&pmt_status);
179 pmts_tree->Branch("HRV",&pmt_HRV);
180
181
182 const JModuleRouter router(dynamics ? dynamics->getDetector() : detector);
183
184 JSummaryFileRouter summary(inputFile);
185
186 JTreeScanner<JDAQEvent, JDAQEvaluator> in(inputFile, inputFile.getLimit());
187
188 typedef vector<JHitL0> JDataL0_t;
189
190 const JBuildL0<JHitL0> buildL0;
191
192 while (in.hasNext()) {
193
194 STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
195
196 pointer_type ps = in.next();
197
198 JDAQEvent* tev = ps;
199
200 summary.update(*tev);
201
202 if (dynamics) {
203 dynamics->update(*tev);
204 }
205
206 JDataL0_t dataL0;
207
208 buildL0(*tev, router, true, back_inserter(dataL0));
209
210 sort(dataL0.begin(), dataL0.end(), JHitL0::compare);
211
212 int module_id_last = -1;
213 int channel_id_last = -1;
214
215 // Loop over hits
216 for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
217
218 JHitW0 hit(*i, summary.getRate(*i, R_Hz));
219
220 // Positions
221 hit_pos_x.push_back(hit.getX());
222 hit_pos_y.push_back(hit.getY());
223 hit_pos_z.push_back(hit.getZ());
224 // Directions
225 hit_dir_x.push_back(hit.getDX());
226 hit_dir_y.push_back(hit.getDY());
227 hit_dir_z.push_back(hit.getDZ());
228 // Time and ToT
229 hit_time.push_back(hit.getT());
230 hit_tot.push_back(hit.getToT());
231 // DOM and channel IDs
232 hit_dom_id.push_back(hit.getModuleID());
233 hit_channel_id.push_back(hit.getPMTAddress());
234
235
236 const JPMT& pmt = router.getPMT(JPMTIdentifier(hit.getModuleID(), hit.getPMTAddress()));
237
238 hit_status.push_back(pmt.getStatus()); // PMT status
239 hit_QE.push_back(pmtParameters.getQE(JPMTIdentifier(hit.getModuleID(), hit.getPMTAddress()))); // PMT QE
240
241 const JDAQSummaryFrame& SummaryFrame = summary.getSummaryFrame(hit.getModuleID());
242 hit_HRV.push_back(SummaryFrame.testHighRateVeto(hit.getPMTAddress())); // Flag for PMT in HRV
243
244 hit_rate.push_back(hit.getR()); // PMT rate
245
246 // Check if the hit is the first one on the channel (1 means first one)
247 if ((hit.getModuleID()==module_id_last) && (hit.getPMTAddress()==channel_id_last)){
248 hit_first.push_back(0);
249 }
250 else{
251 hit_first.push_back(1);
252 }
253 module_id_last = hit.getModuleID();
254 channel_id_last = hit.getPMTAddress();
255
256
257 }
258
259
260 // Loop over all PMTs
261 for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
262
263 JModule module = *i;
264
265 for (size_t pmt = 0; pmt != module.size(); ++pmt) {
266
267 pmt_pos_x.push_back(module[pmt].getX());
268 pmt_pos_y.push_back(module[pmt].getY());
269 pmt_pos_z.push_back(module[pmt].getZ());
270
271 pmt_dir_x.push_back(module[pmt].getDX());
272 pmt_dir_y.push_back(module[pmt].getDY());
273 pmt_dir_z.push_back(module[pmt].getDZ());
274
275 pmt_status.push_back(router.getPMT(JPMTIdentifier(module.getID(), pmt)).getStatus());
276
277 pmt_QE.push_back(pmtParameters.getQE(JPMTIdentifier(module.getID(), pmt)));
278
279 pmt_channel_id.push_back(pmt);
280 pmt_dom_id.push_back(module.getID());
281
282 const JDAQSummaryFrame& SummaryFrame = summary.getSummaryFrame(module.getID());
283 pmt_HRV.push_back(SummaryFrame.testHighRateVeto(pmt));
284
285 pmt_rate.push_back(summary.getRate(JDAQPMTIdentifier(module.getID(), pmt)));
286
287 }
288 }
289
290
291 hits_tree->Fill();
292 hit_first.clear();
293 hit_QE.clear();
294 hit_rate.clear();
295 hit_status.clear();
296 hit_HRV.clear();
297 hit_pos_x.clear();
298 hit_pos_y.clear();
299 hit_pos_z.clear();
300 hit_dir_x.clear();
301 hit_dir_y.clear();
302 hit_dir_z.clear();
303 hit_time.clear();
304 hit_tot.clear();
305 hit_dom_id.clear();
306 hit_channel_id.clear();
307
308
309 pmts_tree->Fill();
310 pmt_pos_x.clear();
311 pmt_pos_y.clear();
312 pmt_pos_z.clear();
313 pmt_dir_x.clear();
314 pmt_dir_y.clear();
315 pmt_dir_z.clear();
316 pmt_channel_id.clear();
317 pmt_dom_id.clear();
318 pmt_QE.clear();
319 pmt_rate.clear();
320 pmt_status.clear();
321 pmt_HRV.clear();
322
323 }
324
325 STATUS(endl);
326 out->Write();
327 out->Close();
328}
329
330
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:72
Direct access to module in detector data structure.
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:75
Auxiliary class for map of PMT parameters.
double getQE(const JPMTIdentifier &id) const
Get QE of given PMT.
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:23
double getR() const
Get rate.
Definition JHitW0.hh:52
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.
double getT() const
Get calibrated time of hit.
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
Acoustic event fit.
Orientation of module.
Dynamic detector calibration.
Definition JDynamics.hh:86
int getStatus() const
Get status.
Definition JStatus.hh:63
Auxiliary class for recursive type list generation.
Definition JTypeList.hh:351
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