Jpp in_tag_pdf_generation
the software that should make you happy
Loading...
Searching...
No Matches
JMonitorMultiplicity.cc
Go to the documentation of this file.
1#include <iostream>
2#include <string>
3
4#include "TFile.h"
5#include "TH1D.h"
6#include "TH2D.h"
7#include "TMath.h"
8#include "TROOT.h"
9#include "TVectorD.h"
10
11#include "JAAnet/JHead.hh"
14#include "JDAQ/JDAQEvaluator.hh"
23#include "Jeep/JParser.hh"
24#include "Jeep/JMessage.hh"
26#include "Jeep/JTimer.hh"
28#include "JLang/JString.hh"
29#include "JROOT/JManager.hh"
32#include "JSupport/JMeta.hh"
35#include "JSupport/JSupport.hh"
38#include "JTrigger/JHit.hh"
39#include "JTrigger/JHitR0.hh"
40#include "JTrigger/JMatchL0.hh"
44#include "JTools/JRange.hh"
45
46/**
47 * \file
48 *
49 * Auxiliary program to analyse coincidence multiplicity;
50 * \author mlincetto, mkarel
51 */
52int main(int argc, char **argv)
53{
54 using namespace std;
55 using namespace KM3NETDAQ;
56 using namespace JPP;
57
58 //-----------------------------------------------------------
59 // parameter interface
60 //-----------------------------------------------------------
61
63 JLimit_t& numberOfEvents = inputFile.getLimit();
64 string outputFile;
65 string detectorFile;
66 Double_t TMax_ns;
67 JRange<double> rateVeto_Hz;
68 JRange<double> totVeto_ns;
69 JROOTClassSelector selector;
70 int debug;
71 int filterLevel;
72 int muteChannel;
73 bool twoDim;
74 bool monitorOccupancy;
75 bool notJoin;
76
77 try {
78
79 JParser<> zap("Auxiliary program to estimate PMT and hit multiplicities.");
80
81 zap['f'] = make_field(inputFile);
82 zap['o'] = make_field(outputFile) = "monitormultiplicity.root";
83 zap['a'] = make_field(detectorFile);
84 zap['n'] = make_field(numberOfEvents) = JLimit::max();
85 zap['T'] = make_field(TMax_ns) = 10.0; // [ns]
86 zap['V'] = make_field(rateVeto_Hz) = JRange<double>(0, 1e5); // [Hz]
87 zap['t'] = make_field(totVeto_ns) = JRange<double>(0, 1e4); // min and max value of ToT of hits to be monitored [ns]
89 zap['d'] = make_field(debug) = 1;
90 zap['M'] = make_field(muteChannel) = -1;
91 zap['F'] = make_field(filterLevel, "0 = raw data; 1 = filtered data; 2+ = only clean frames") = 1;
92 zap['D'] = make_field(twoDim, "2D mode for background subtraction");
93 zap['O'] = make_field(monitorOccupancy, "produces 2D PMT vs multiplicity plots");
94 zap['j'] = make_field(notJoin, "do not join consecutive hits on same PMT (diagnostic purpose)");
95
96 zap(argc, argv);
97 }
98 catch(const exception &error) {
99 FATAL(error.what() << endl);
100 }
101
102
103 //-----------------------------------------------------------
104 // check the input parameters
105 //-----------------------------------------------------------
106
107 JTriggerParameters parameters;
108 bool hasTriggerParameters = true;
109
110 try {
111 parameters = getTriggerParameters(inputFile);
112 NOTICE("Get trigger parameters from input." << endl);
113 }
114 catch(const JException& error) {
115 ERROR("No trigger parameters from input." << endl);
116 hasTriggerParameters = false;
117
118 }
119
120 // keep track of prescale factor for proper MC normalization.
121 int prescale = 1;
122
123 if (hasTriggerParameters) {
124
125 // for L1 (most common for sea analysis) check that time window is not larger than in trigger
126
127 if (parameters.writeL1.prescale > 0 && (selector == "JDAQTimeslice" || selector == "JDAQTimesliceL1") &&
128 parameters.TMaxLocal_ns < TMax_ns) {
129 FATAL("TMax_ns (-T) " << TMax_ns << " ns is larger than in the trigger " << parameters.TMaxLocal_ns << " ns." << endl);
130 }
131
132 if (selector == "JDAQTimeslice") {
133 // Jpp v8 legacy compatibility
134 prescale = (parameters.writeL1.prescale != 0) ? parameters.writeL1.prescale : parameters.writeL0.prescale;
135 }
136
137 if (selector == "JDAQTimesliceL0") {
138 prescale = parameters.writeL0.prescale;
139 }
140
141 if (selector == "JDAQTimesliceL1") {
142 prescale = parameters.writeL1.prescale;
143 }
144
145 if (selector == "JDAQTimesliceSN") {
146 prescale = parameters.writeSN.prescale;
147 }
148
149 }
150
151 if (prescale == 0) {
152 FATAL("[R] According to trigger parameters, the " << selector << " stream should be empty.");
153 } else {
154 NOTICE("[R] Prescale factor is " << prescale << endl);
155 }
156
157 if (!totVeto_ns.is_valid()) {
158 FATAL("Invalid time over threshold " << totVeto_ns << endl);
159 }
160
161 //-----------------------------------------------------------
162 // load the detector file
163 //-----------------------------------------------------------
164
166
167 try {
168 load(detectorFile, detector);
169 }
170 catch(const JException& error) {
171 FATAL(error);
172 }
173
174 if (detector.empty())
175 FATAL("Empty detector." << endl);
176
177 const JModuleRouter router(detector);
178 JSummaryRouter summaryRouter;
179
180 int detSize = getNumberOfModules(detector);
181
182 int detID = detector.getID();
183
184 JDetectorAddressMap& detectorAddressMap = getDetectorAddressMap(detID);
185
186
187 //-----------------------------------------------------------
188 // summary of input parameters
189 //-----------------------------------------------------------
190
191 NOTICE("[R] Frame time [ms] " << getFrameTime() * 1e-6 << endl);
192 NOTICE("[R] Tmax [ns] " << TMax_ns << endl);
193 NOTICE("[R] Rate range [Hz] " << rateVeto_Hz << endl);
194 NOTICE("[R] ToT range [ns] " << totVeto_ns << endl);
195 NOTICE("[R] Timeslice stream " << selector << endl);
196 NOTICE("[R] Detector file has " << detSize << " DOMs." << endl);
197
198 //-----------------------------------------------------------
199 // initialise histograms
200 //-----------------------------------------------------------
201
202 // nominal livetime + effective livetime for each DOM
203 TH1D* h_livetime = new TH1D("LiveTime", "Livetime", 1 + detSize, 0.0, 1.0 + detSize);
204
205 int ibin = 1;
206 h_livetime->GetXaxis()->SetBinLabel(ibin++, "Nominal");
207 for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
208 if (!module->empty()) {
209 h_livetime->GetXaxis()->SetBinLabel(ibin++, TString(getLabel(*module)));
210 }
211 }
212
213 h_livetime->GetYaxis()->SetTitle("Livetime (s)");
214
215
216 // multiplicity histograms
217 const int nx = 1 + NUMBER_OF_PMTS;
218 const double xmin = -0.5;
219 const double xmax = nx - 0.5;
220
221 typedef JManager<TString, TH1D> JManager_t;
222 typedef JManager<TString, TH2D> JManager2D_t;
223 JManager_t H(new TH1D("%_H", NULL, nx, xmin, xmax)); H->Sumw2();
224 JManager_t P(new TH1D("%_P", NULL, nx, xmin, xmax)); P->Sumw2();
225
226 JManager2D_t HT(new TH2D("%_HT", NULL, nx, xmin, xmax, 100, -TMax_ns, +TMax_ns));
227
228 TH2D* CO_proto = new TH2D("%_CO", NULL, NUMBER_OF_PMTS, 0.5, 0.5 + NUMBER_OF_PMTS, NUMBER_OF_PMTS, -0.5, -0.5 + NUMBER_OF_PMTS);
229
230 if (monitorOccupancy) {
232 setAxisLabels(CO_proto->GetYaxis(), memo);
233 }
234
235 JManager2D_t CO(CO_proto);
236
237 // benchmark histogram for time distribution of high-multiplicities
238 const int sn_mul_min = 6;
239 TH1D* time_distr = new TH1D("T_distr", "Coincidence (M > 6) distribution in timeslice", 100, 0, 1e8);
240
241 //-----------------------------------------------------------
242 // output file and auxiliary data structures
243 //-----------------------------------------------------------
244
245 // output file
246 TFile out(outputFile.c_str(), "recreate");
247
248 // auxiliary counters
249
250 int count_HRV = 0; // high-rate-veto counter
251 int count_FAF = 0; // fifo-almost-full counter
252 int count_URV = 0; // user-rate-veto counter
253
254 int treeMismatchCount = 0;
255
256 //-----------------------------------------------------------
257 // file pre-processing
258 //-----------------------------------------------------------
259
261
263
264 vector<JHitR0> filteredData;
265
266 typedef JSuperFrame2D<JHitR0> JSuperFrame2D_t;
267 typedef JSuperFrame1D<JHitR0> JSuperFrame1D_t;
268 const JMatchL0<JHitR0> match(TMax_ns);
269
270 counter_type counter = 0;
271
272 // process Monte Carlo aanet header
273 bool hasMCHeader = true ;
274 bool isMupage = false;
275 double liveTimeMC = 0 ;
276 double liveTimeMCErr = 0 ;
277
278 JHead mc_header;
279
280 try {
281 mc_header = (JHead)getHeader(inputFile);
282 }
283 catch (const JException& error) {
284 hasMCHeader = false;
285 }
286
287 NOTICE("[R] File source is " << (hasMCHeader ? "MC" : "DAQ") << "." << endl);
288 if (hasMCHeader) {
289 liveTimeMC = mc_header.livetime.numberOfSeconds;
290 liveTimeMCErr = mc_header.livetime.errorOfSeconds;
291 if (liveTimeMC > 0) {
292 isMupage = true;
293 NOTICE("[R] mupage live time is (" << liveTimeMC << " +/- " << liveTimeMCErr << ") s." << endl);
294 } else {
295 NOTICE("[R] aanet header found but mupage live time is zero. Likely not a mupage file.");
296 }
297 }
298
299 // retrieving and exporting run number to tree metadata;
300
301 int runNumber = 0;
302
303 if (summaryFile.hasNext()) {
304 runNumber = summaryFile.getEntry(0)->getRunNumber();
305 NOTICE("[R] Processing run #" << runNumber << endl);
306 summaryFile.rewind();
307 }
308
309 //-----------------------------------------------------------
310 // file processing: loop over timeslices
311 //-----------------------------------------------------------
312
313 for ( ; in.hasNext() && counter != inputFile.getLimit(); ++counter) {
314
315 STATUS("Entry: " << setw(10) << counter << '\r');
316
317 const JDAQTimeslice* timeslice = in.next();
318 const Long64_t index = summaryFile.find(*timeslice);
319 JDAQSummaryslice* summary = (index != -1 ? summaryFile.getEntry(index) : NULL);
320
321 summaryRouter.update(summary);
322
323 // check correspondence between summary data and timeslice data
324 if (summary != NULL) {
325 if (timeslice->getFrameIndex() != summary->getFrameIndex()) {
326 ++treeMismatchCount;
327 ERROR("[R>T] Frame indices do not match [counter=" << counter << ", timeslice=" << timeslice->getFrameIndex() << ", summaryslice=" << summary->getFrameIndex() << "]" << endl);
328 if (treeMismatchCount > 1) {
329 FATAL("ROOT TTrees not aligned. Run #" << runNumber << endl);
330 }
331 continue;
332 }
333 }
334
335 //-----------------------------------------------------------
336 // timeslice processing: loop over superframes
337 //-----------------------------------------------------------
338 for (JDAQTimeslice::const_iterator super_frame = timeslice->begin(); super_frame != timeslice->end(); ++super_frame) {
339
340 // empty superframes are not skipped because SN timeslices could have empty superframes even if the module was active;
341 // if DAQ error or FIFO full -> frame counts as dead time;
342 // if frame is missing but no errors -> frame counts as live time (SN stream has many empty frames);
343
344 int moduleID = super_frame->getModuleID();
345
346 if (!router.hasModule(moduleID)) {
347 // module is not mapped in detector file, nothing can be done
348 continue;
349 }
350
351 JModuleAddressMap moduleAddressMap = detectorAddressMap.get(moduleID);
352
353 //----------------------------------------------------------
354 // retrieval of summary data from superframe or summaryframe
355 //----------------------------------------------------------
356
357 vector<double> rate_Hz(NUMBER_OF_PMTS, 0.0);
358 JDAQFrameStatus status;
359
360 if (summary == NULL) {
361 WARNING("[R>T>F] Missing summary for timeslice." << endl);
362 // PMT rates estimated from timeslice data
363 for (JDAQSuperFrame::const_iterator i = super_frame->begin(); i != super_frame->end(); ++i) {
364 rate_Hz[i->getPMT()] += 1e9 / getFrameTime();
365 }
366
367 // frame status retrieved from superframe
368 status = super_frame->getDAQFrameStatus();
369
370 } else {
371
372 // summary frame lookup
373
374 JDAQSummaryFrame summary_frame;
375
376 if (!summaryRouter.hasSummaryFrame(super_frame->getModuleIdentifier())) {
377 summaryRouter.print(cout);
378 FATAL("[R>D>F] Summary frame not found, module = " << super_frame->getModuleID() << endl);
379 } else {
380 summary_frame = summaryRouter.getSummaryFrame(super_frame->getModuleID());
381 }
382
383 // PMT rates recovered from summaryslice data
384 for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
385 rate_Hz[i] = summary_frame.getRate(i);
386 }
387
388 // frame status retrieved from summaryframe
389 status = summary_frame.getDAQFrameStatus();
390
391 }
392
393 //---------------------------
394 // processing of summary data
395 //---------------------------
396
397 // DAQ error (e.g. missing UDP trailer) -> discard frame
398 if (!status.testDAQStatus()) {
399 WARNING("[R>D>F] DAQ status not okay in timeslice #" << timeslice->getFrameIndex() << ", DOM=" << super_frame->getModuleID() << endl);
400 continue;
401 }
402
403 vector<bool> veto(NUMBER_OF_PMTS, false);
404 bool moduleLiveTime = 0;
405
406 // default status of module is alive if at least one PMT is active
407 for (vector<double>::const_iterator r = rate_Hz.begin(); !moduleLiveTime && r != rate_Hz.end(); ++r) {
408 if ((*r) > 0) {
409 moduleLiveTime = 1;
410 }
411 }
412
413 int frame_HRV_count = status.countHighRateVeto();
414 int frame_FAF_count = status.countFIFOStatus();
415 int frame_URV_count = 0; // user rate veto
416
417 for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
418 if ((veto[i] = (veto[i] || !rateVeto_Hz(rate_Hz[i])))) {
419 frame_URV_count++;
420 }
421 }
422
423 count_HRV += frame_HRV_count;
424 count_FAF += frame_FAF_count;
425 count_URV += frame_URV_count;
426
427 int frame_VETO_count = frame_HRV_count + frame_FAF_count + frame_URV_count;
428
429
430 if (filterLevel >= 2 && frame_VETO_count > 0) {
431 moduleLiveTime = 0;
432 fill(veto.begin(), veto.end(), true);
433 } else if (filterLevel >= 1 && frame_HRV_count + frame_FAF_count > 0) {
434 for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
435 veto[i] = (veto[i] || (status.testFIFOStatus(i) || status.testHighRateVeto(i)));
436 }
437 }
438
439 // recover DOM in case the dead channel corresponds to the muted channel
440
441 if (muteChannel != -1) {
442
443 // veto conditions are not exclusive
444 int mute_VETO_count =
445 status.testHighRateVeto(muteChannel) +
446 status.testFIFOStatus(muteChannel ) +
447 !rateVeto_Hz(rate_Hz[muteChannel] );
448
449 if (mute_VETO_count == frame_VETO_count) {
450 moduleLiveTime = 1;
451 fill(veto.begin(), veto.end(), false);
452 }
453
454 veto[muteChannel] = true;
455
456 }
457
458 // get address and module for current superframe DOM;
459 const JModuleAddress& address = router.getAddress(super_frame->getModuleID());
460 const JModule& module = detector.getModule(address);
461 const TString moduleLabel = getLabel(module);
462
463 // acknowledge DOM is live in this timeslice for analysis purpose;
464 if (moduleLiveTime) {
465 h_livetime->Fill(moduleLabel, getFrameTime() * 1e-9 * moduleLiveTime);
466 }
467
468 //--------------------------------------------------------------------
469 // superframe processing: hit joining, calibration, filter, processing
470 //--------------------------------------------------------------------
471
472 JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*super_frame, router.getModule(super_frame->getModuleID()));
473
474 if (!notJoin) {
475 for (JSuperFrame2D_t::iterator i = buffer.begin(); i != buffer.end(); ++i) {
476 i->join(match);
477 }
478 }
479
480 JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer);
481
482 filteredData.clear();
483
484 for (vector<JHitR0>::const_iterator h = data.begin(); h != data.end(); ++h) {
485
486 const int pmt = h->getPMT();
487
488 // veto and calibration;
489 if (!veto[pmt] &&
490 rateVeto_Hz(rate_Hz[pmt]) &&
491 totVeto_ns(h->getToT()) ) {
492 filteredData.push_back(*h);
493 }
494 }
495
496 // processing signal;
497
498 if (filteredData.size() > 1) {
499
500 TH1D* h_hit = H[moduleLabel];
501 TH1D* h_pmt = P[moduleLabel];
502
503 TH2D* h2_hit = HT[moduleLabel];
504 TH2D* h2_co = CO[moduleLabel];
505
506 sort(filteredData.begin(), filteredData.end());
507
508 //-----------------------------------------------------------
509 // superframe processing: loop on calibrated hits
510 //-----------------------------------------------------------
511 for (vector<JHitR0>::const_iterator p = filteredData.begin(); p != filteredData.end(); ) {
512
513 vector<JHitR0>::const_iterator q = p;
514
515 // coincidence building using a fixed time window
516 while (++q != filteredData.end() && q->getT() - p->getT() <= TMax_ns ) {}
517
518 // multiplicity estimation
519 int hit_multiplicity = distance(p,q);
520 set<int> pmthit;
521 for (vector<JHitR0>::const_iterator __p = p ; __p != q; ++__p ) { pmthit.insert(__p->getPMT()) ; }
522 int pmt_multiplicity = pmthit.size() ;
523
524 // histogram filling
525 h_pmt->Fill(pmt_multiplicity);
526 h_hit->Fill(hit_multiplicity);
527
528 if (twoDim && hit_multiplicity > 1) {
529 const double W = factorial(hit_multiplicity, 2);
530 for (vector<JHitR0>::const_iterator __p = p; __p != q; ++__p) {
531 for (vector<JHitR0>::const_iterator __q = __p; ++__q != q; ) {
532 double dt = JCombinatorics::getSign(__p->getPMT(), __q->getPMT()) * (__q->getT() - __p->getT());
533 h2_hit->Fill(hit_multiplicity, dt, 1.0/W);
534 }
535 }
536 }
537
538 if (monitorOccupancy) {
539 for (set<int>::const_iterator r = pmthit.begin(); r != pmthit.end(); ++r) {
540 JPMTPhysicalAddress pmtAddress = moduleAddressMap.getAddressTranslator(*r);
541 TString pmtLabel(pmtAddress.toString());
542 h2_co->Fill(pmt_multiplicity, pmtLabel, 1.0/pmt_multiplicity);
543 }
544 }
545
546 // benchmarking
547 if (hit_multiplicity >= sn_mul_min) {
548 time_distr->Fill(p->getT());
549 }
550
551 p = q;
552
553 }
554
555 }
556
557 }
558
559 }
560
561 NOTICE("[R] " << counter << " timeslices processed." << endl);
562
563 if (isMupage) {
564 h_livetime->Scale(liveTimeMC / (getFrameTime() * 1e-9 * counter * prescale));
565 }
566
567 // nominal live time;
568 double nominalLiveTime = (isMupage ? liveTimeMC : counter * getFrameTime() * 1e-9) / prescale;
569 h_livetime->Fill("Nominal", nominalLiveTime);
570
571 out.cd();
572
573 int nLiveDOMs = H.size();
574
575 h_livetime->Write();
576
577 // writing hit multiplicities is redundant given new options
578 // H.Write(out);
579
580 P.Write(out);
581
582 if (twoDim) {
583 HT.Write(out);
584 }
585
586 if (monitorOccupancy) {
587 CO.Write(out);
588 }
589
590 time_distr->Write();
591
592 double rate_HRV = count_HRV / (1.0 * counter);
593 double rate_FAF = count_FAF / (1.0 * counter);
594 double rate_URV = count_URV / (1.0 * counter);
595 double biolumIndex = (rate_HRV + rate_FAF) / (nLiveDOMs * NUMBER_OF_PMTS);
596 NOTICE("[R] Average HRV count per timeslice: " << rate_HRV << endl);
597 NOTICE("[R] Average FIFO-almost-full count per timeslice: " << rate_FAF << endl);
598 NOTICE("[R] Average user rate veto count per timeslice: " << rate_URV << endl);
599 NOTICE("[R] " << H.size() << " DOMs were active in the run." << endl);
600 NOTICE("[R] Bioluminescence proxy index: " << biolumIndex << endl);
601
602 putObject(out, JMeta(argc,argv));
603
604 out.Close();
605
606 NOTICE(endl);
607
608 return (counter ? 0 : 1);
609
610}
string outputFile
KM3NeT DAQ constants, bit handling, etc.
Detector support kit.
Data structure for detector geometry and calibration.
Basic data structure for L0 hit.
Dynamic ROOT object management.
Match operator for consecutive hits.
Auxiliary methods for mathematics.
General purpose messaging.
#define STATUS(A)
Definition JMessage.hh:63
#define NOTICE(A)
Definition JMessage.hh:64
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define WARNING(A)
Definition JMessage.hh:65
ROOT I/O of application specific meta data.
Direct access to module in detector data structure.
int main(int argc, char **argv)
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
Auxiliary class to define a range between two values.
ROOT TTree parameter settings of various packages.
Basic data structure for time and time over threshold information of hit.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Monte Carlo run header.
Definition JHead.hh:1236
JAANET::livetime livetime
Definition JHead.hh:1604
Lookup table for PMT addresses in detector.
const JModuleAddressMap & get(const int id) const
Get module address map.
virtual const JModuleAddressMap & getDefaultModuleAddressMap() const =0
Get default module address map.
Detector data structure.
Definition JDetector.hh:96
Lookup table for PMT addresses in optical module.
const JPMTAddressTranslator & getAddressTranslator(const int tdc) const
Get PMT address translator.
Address of module in detector data structure.
Router for direct addressing of module data in detector data structure.
bool hasModule(const JObjectID &id) const
Has module.
const JModuleAddress & getAddress(const JObjectID &id) const
Get address of module.
const JModule & getModule(const JObjectID &id) const
Get module parameters.
Data structure for a composite optical module.
Definition JModule.hh:75
Data structure for PMT physical address.
std::string toString() const
Convert PMT physical address to string.
General exception.
Definition JException.hh:24
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
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys.
Definition JManager.hh:47
General purpose class for object reading from a list of file names.
Router for fast addressing of summary data in KM3NETDAQ::JDAQSummaryslice data structure as a functio...
bool hasSummaryFrame(const JDAQModuleIdentifier &module) const
Has summary frame.
const JDAQSummaryFrame & getSummaryFrame(const JDAQModuleIdentifier &module) const
Get summary frame.
void update(const JDAQSummaryslice *ps)
Update router.
Template definition for direct access of elements in ROOT TChain.
static int getSign(const int first, const int second)
Sign of pair of indices.
Range of values.
Definition JRange.hh:42
bool is_valid() const
Check validity of range.
Definition JRange.hh:311
L0 match criterion.
Definition JMatchL0.hh:29
1-dimensional frame with time calibrated data from one optical module.
2-dimensional frame with time calibrated data from one optical module.
int getFrameIndex() const
Get frame index.
bool testFIFOStatus() const
Test FIFO status.
int countHighRateVeto() const
Count high-rate veto status.
const JDAQFrameStatus & getDAQFrameStatus() const
Get DAQ frame status.
int countFIFOStatus() const
Count FIFO status.
bool testHighRateVeto() const
Test high-rate veto status.
bool testDAQStatus() const
Test DAQ status of packets.
Hit data structure.
Definition JDAQHit.hh:35
Data storage class for rate measurements of all PMTs in one module.
double getRate(const int tdc, const double factor=1.0) const
Get count rate.
std::string getLabel(const JLocation &location)
Get module label for monitoring and other applications.
Definition JLocation.hh:247
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
JDetectorAddressMap & getDetectorAddressMap()
Get detector address map.
void setAxisLabels(TAxis *axis, const JModuleAddressMap &memo)
Set axis with PMT address labels.
long long int factorial(const long long int n)
Determine factorial.
static const double H
Planck constant [eV s].
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
bool putObject(TDirectory &dir, const TObject &object)
Write object to ROOT directory.
std::set< JROOTClassSelector > getROOTClassSelection(const bool option=false)
Get ROOT class selection.
Long64_t counter_type
Type definition for counter.
JTriggerParameters getTriggerParameters(const JMultipleFileScanner_t &file_list)
Get trigger parameters.
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
static const JModuleCounter getNumberOfModules
Function object to count unique modules.
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
double getFrameTime()
Get frame time duration.
Definition JDAQClock.hh:162
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition JDAQ.hh:26
Detector file.
Definition JHead.hh:227
double numberOfSeconds
Live time [s].
Definition JHead.hh:896
double errorOfSeconds
Uncertainty on live time [s].
Definition JHead.hh:897
Auxiliary class to select ROOT class based on class name.
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 class for ROOT I/O of application specific meta data.
Definition JMeta.hh:72