Jpp 20.0.0-195-g190c9e876
the software that should make you happy
Loading...
Searching...
No Matches
JCalibrateK40.cc
Go to the documentation of this file.
1#include <string>
2#include <iostream>
3#include <iomanip>
4#include <limits>
5#include <array>
6
7#include "TROOT.h"
8#include "TFile.h"
9#include "TH1D.h"
10#include "TH2D.h"
11
14
17#include "JDAQ/JDAQEvaluator.hh"
20#include "JTrigger/JHitR0.hh"
21#include "JTrigger/JMatchL0.hh"
34#include "JSupport/JSupport.hh"
35#include "JSupport/JMeta.hh"
37#include "JTools/JRange.hh"
39
41#include "JROOT/JRootToolkit.hh"
45
46#include "Jeep/JPrint.hh"
47#include "Jeep/JParser.hh"
48#include "Jeep/JMessage.hh"
49
50namespace {
51
52 using JLANG::JObjectID;
56
57 // Methods for background estimation
58
59 static const char* const rates_t = "rates"; //!< singles rates from summary data
60 static const char* const count_t = "counts"; //!< singles rates from L0 data
61 static const char* const tails_t = "tails"; //!< tails of time residual distribution
62 static const char* const rndms_t = "randoms"; //!< random time offsets (only L0 data)
63
64 /**
65 * Maximal time-over-threshold.
66 */
67 static const double ToTmax_ns = 250;
68
69 /**
70 * Auxiliary data structure for histogram management.
71 */
72 struct JHistogram :
73 public JObjectID,
74 public JCombinatorics
75 {
76 /**
77 * Default constructor.
78 */
79 JHistogram() :
80 h2s(NULL),
81 h1b(NULL),
82 h1L(NULL)
83 {}
84
85
86 /**
87 * Constructor.
88 *
89 * \param module module
90 * \param axis axis
91 */
92 JHistogram(const JModule& module, const JAbstractHistogram<Double_t>& axis) :
93 JObjectID(module.getID())
94 {
95 using namespace JPP;
96
97 // sort the PMT pairs by opening angle in the order largest angle -> smallest angle
98
99 this->configure(module.size());
100
101 this->sort(JPairwiseComparator(module));
102
103 const Int_t nx = this->getNumberOfPairs();
104 const Double_t xmin = -0.5;
105 const Double_t xmax = nx - 0.5;
106
107 h2s = new TH2D(MAKE_CSTRING(module.getID() << _2S), NULL, nx, xmin, xmax, axis.getNumberOfBins(), axis.getLowerLimit(), axis.getUpperLimit());
108 h1b = new TH1D(MAKE_CSTRING(module.getID() << _1B), NULL, nx, xmin, xmax);
109 h1L = new TH1D(MAKE_CSTRING(module.getID() << _1L), NULL, nx, xmin, xmax);
110
111 h2s->Sumw2();
112 h1b->Sumw2();
113 h1L->Sumw2();
114 }
115
116 /**
117 * Check validity.
118 *
119 * \return true if valid; else false
120 */
121 bool is_valid() const
122 {
123 return (h2s != NULL &&
124 h1b != NULL &&
125 h1L != NULL);
126 }
127
128 TH2D* h2s; //!< time distribution per PMT pair
129 TH1D* h1b; //!< background rates per PMT pair
130 TH1D* h1L; //!< livetimes per PMT pair
131 };
132}
133
134
135/**
136 * \file
137 *
138 * Auxiliary program to determine PMT parameters from K40 data.
139 *
140 * By default, the combinatorial background is estimated from the singles rate.\n
141 * In case L0 data are taken (and no summary data are available),
142 * the singles rate can be determined from the amount of L0 data.\n
143 * One can also estimate the background from the tails in the time residual distribution.\n
144 * In that case, option -B can be used to specify the minimal and maximal time offset in ns.\n
145 * For example -B "10 20".\n
146 * Note that the minimal and maximal time should be larger that the width of
147 * the time residual distribution and less than the time window of the L1 coincidence, respectively.\n
148 * The time window is applied symmetrically around a time offset of zero.\n
149 * Finally, if L0 data are available, one can estimate the background using
150 * the same procedure but with a random (read wrong) time calibration.
151 *
152 * The option -t can be used to specify three time-over-threshold values.\n
153 * The time-over-threshold of each hit in the coincidence should be greater or equal the first value and less or equal the last value.\n
154 * Furthermore, the time-over-threshold of one of the two hits should be greater or equals the second value.\n
155 * If you change the (default) values, you may also want to change the corresponding values in JEditPMTParameters.cc.
156 *
157 * The option -M can be used to specify a range of multiplicities.\n
158 * For example -M "2 2" will strictly select two-fold coincidences.
159 * Note that such a selection can bias the result.
160 *
161 * Also consult <a href="https://common.pages.km3net.de/jpp/Detector_calibration.PDF">documentation</a>.
162 * \author mdejong
163 */
164int main(int argc, char **argv)
165{
166 using namespace std;
167 using namespace JPP;
168 using namespace KM3NETDAQ;
169
170 JMultipleFileScanner_t inputFile;
171 counter_type numberOfEvents;
172 string outputFile;
173 string detectorFile;
174 Double_t TMax_ns;
175 JRange<double> rateRange_Hz;
176 array <double, 3> totRange_ns;
177 JRange<double> Tail_ns;
178 JRange<int> multiplicity;
179 double deadTime_us;
180 JROOTClassSelector selector;
181 string background;
182 JPreprocessor option;
183 int debug;
184
185 try {
186
187 JParser<> zap("Auxiliary program to determine PMT parameters from K40 data.");
188
189 zap['f'] = make_field(inputFile, "input file.");
190 zap['o'] = make_field(outputFile, "output file.") = "calibrate_k40.root";
191 zap['n'] = make_field(numberOfEvents) = JLimit::max();
192 zap['a'] = make_field(detectorFile, "detector file.");
193 zap['T'] = make_field(TMax_ns, "time window [ns].") = 20.0;
194 zap['V'] = make_field(rateRange_Hz, "PMT rate range [Hz].") = JRange<double>(0.0, 20.0e3);
195 zap['t'] = make_field(totRange_ns, "PMT time-over-threshold range [ns].") = { 4.0, 7.0, ToTmax_ns };
196 zap['b'] = make_field(background, "background estimation method.") = rates_t, count_t, tails_t, rndms_t;
197 zap['B'] = make_field(Tail_ns, "time window used for background estimation.") = JRange<double>(15.0, 20.0);
198 zap['M'] = make_field(multiplicity, "multiplicity range of hits on DOM.") = JRange<int>(2, 31);
199 zap['D'] = make_field(deadTime_us, "L1 dead time (us)") = 0.0;
200 zap['C'] = make_field(selector, "timeslice selector, e.g. JDAQTimesliceL1.") = getROOTClassSelection<JDAQTimesliceTypes_t>();
201 zap['O'] = make_field(option, "hit pre-processing option.") = JPreprocessor::getOptions(JPreprocessor::filter_t);
202 zap['d'] = make_field(debug, "debug flag.") = 1;
203
204 zap(argc, argv);
205 }
206 catch(const exception &error) {
207 FATAL(error.what() << endl);
208 }
209
210 //-----------------------------------------------------------
211 // check the input parameters
212 //-----------------------------------------------------------
213
214
215 if (selector == JDAQTimesliceL2::Class_Name() ||
216 selector == JDAQTimesliceSN::Class_Name()) {
217 FATAL("Option -C <selector> " << selector << " not compatible with calibration method." << endl);
218 }
219
220 if (selector == JDAQTimeslice ::Class_Name() ||
221 selector == JDAQTimesliceL1::Class_Name()) {
222
223 JTriggerParameters parameters;
224
225 try {
226 parameters = getTriggerParameters(inputFile);
227 }
228 catch(const JException& error) {
229 FATAL("No trigger parameters from input:" << error.what() << endl);
230 }
231
232 if ((selector == JDAQTimeslice ::Class_Name() && parameters.writeL1.prescale > 0) ||
233 (selector == JDAQTimesliceL1::Class_Name())) {
234
235 if (parameters.TMaxLocal_ns < TMax_ns) {
236 FATAL("Option -T <TMax_ns> = " << TMax_ns << " is larger than in the trigger " << parameters.TMaxLocal_ns << endl);
237 }
238
239 if (background == rndms_t ||
240 background == count_t) {
241 FATAL("Option -C <selector> " << selector << " incompatible with option -b <background> " << background << endl);
242 }
243 }
244 }
245
246 if (!multiplicity.is_valid() || multiplicity.getLowerLimit() < 2 || multiplicity.getUpperLimit() > NUMBER_OF_PMTS) {
247 FATAL("Invalid option -M <multiplicity> " << multiplicity << endl);
248 }
249
250 if (totRange_ns[0] > totRange_ns[1] ||
251 totRange_ns[1] > totRange_ns[2]) {
252 FATAL("Invalid option -t <totRange_ns> " << JEEPZ() << totRange_ns << endl);
253 }
254
255 if (background == tails_t) {
256
257 if (!Tail_ns.is_valid()) {
258 FATAL("Invalid option -B <Tail_ns> " << Tail_ns << endl);
259 }
260
261 if (Tail_ns.getUpperLimit() > TMax_ns) {
262 FATAL("Invalid option -B <Tail_ns> " << Tail_ns << "; upper limit larger than option -T <TMax_ns> " << TMax_ns << endl);
263 }
264 }
265
266 //-----------------------------------------------------------
267 // load detector file
268 //-----------------------------------------------------------
269
271
272 try {
273 load(detectorFile, detector);
274 }
275 catch(const JException& error) {
276 FATAL(error);
277 }
278
279 if (detector.empty()) {
280 FATAL("Empty detector." << endl);
281 }
282
283 const JModuleRouter router(detector);
284
285 //-----------------------------------------------------------
286 // determine time offsets for background method rndms_t
287 //-----------------------------------------------------------
288
289 vector<double> t0; // time offsets for random background evaluation [ns]
290
291 for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
292 t0.push_back(i * 2 * TMax_ns);
293 }
294
295 //-----------------------------------------------------------
296 // correction factor for rate due to specified time-over-threshold range
297 //-----------------------------------------------------------
298
300
301 const int NPE = 1;
302 const double RTU = (cpu.getIntegralOfChargeProbability(cpu.getNPE(totRange_ns[0]), cpu.getNPE(totRange_ns[2]), NPE)
303 /
304 cpu.getIntegralOfChargeProbability(cpu.getNPE(0.0), cpu.getNPE(ToTmax_ns), NPE));
305
306
307 //-----------------------------------------------------------
308 // initialise histograms
309 //-----------------------------------------------------------
310
311 vector<JHistogram> zmap(detector.size());
312
313 const Double_t ymin = -floor(TMax_ns) + 0.5;
314 const Double_t ymax = +floor(TMax_ns) - 0.5;
315 const Int_t ny = (Int_t) ((ymax - ymin) / 1.0);
316
317 for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
318
319 const JModuleAddress& address = router.getAddress(module->getID());
320
321 if (!module->empty()) {
322
323 NOTICE("Booking histograms for module " << module->getID() << endl);
324
325 zmap[address.first] = JHistogram(*module, JAbstractHistogram<Double_t>(ny, ymin, ymax));
326 }
327 }
328
329
330 typedef JHitR0 hit_type;
331 typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
332 typedef JSuperFrame1D<hit_type> JSuperFrame1D_t;
333
334 const JMatchL0<hit_type> match(TMax_ns); // time window self-coincidences [ns]
335
336 const double deadTime_ns = deadTime_us * 1e3;
337
338
339 TFile out(outputFile.c_str(), "recreate");
340
341 putObject(out, JMeta(argc, argv));
342
343 counter_type counter = 0;
344
345 for (JMultipleFileScanner_t::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
346
347 STATUS("Processing: " << *i << endl);
348
349 JSummaryFileRouter summary(*i);
350
353
354 for (JDAQHeader header; in.hasNext() && counter != numberOfEvents; ++counter) {
355
356 STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
357
358 const JDAQTimeslice* timeslice = in.next();
359
360 if (header.getDetectorID() != timeslice->getDetectorID() ||
361 header.getRunNumber () != timeslice->getRunNumber ()) {
362
363 header = timeslice->getDAQHeader();
364
365 putObject(out, header);
366 }
367
368 if (background == rates_t) {
369
370 summary.update(timeslice->getDAQHeader());
371
372 if (timeslice->getFrameIndex() != summary.getFrameIndex()) {
373
374 ERROR("Frame indices do not match at "
375 << "[counter = " << counter << "]"
376 << "[timeslice = " << timeslice->getFrameIndex() << "]"
377 << "[summaryslice = " << summary.getFrameIndex() << "]" << endl);
378
379 continue;
380 }
381 }
382
383 for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
384
385 if (router.hasModule(frame->getModuleID())) {
386
387 const JModule& module = router.getModule(frame->getModuleID());
388
389 //-----------------------------------------------------------
390 // determine background rates and veto per PMT
391 //-----------------------------------------------------------
392
393 vector<double> rate_Hz(NUMBER_OF_PMTS, 0.0);
394 vector<bool> veto (NUMBER_OF_PMTS, false);
395
396 if (background == count_t) {
397
398 for (JDAQSuperFrame::const_iterator i = frame->begin(); i != frame->end(); ++i) {
399 rate_Hz[i->getPMT()] += RTU * 1e9 / getFrameTime();
400 }
401
402 } else if (background == tails_t) {
403
404 // see below at: background == tails_t
405
406
407 } else if (background == rndms_t) {
408
409 // see below at: background == rndms_t
410
411 } else if (background == rates_t) {
412
413 if (summary.hasSummaryFrame(frame->getModuleID())) {
414
415 const JDAQSummaryFrame& sum = summary.getSummaryFrame(frame->getModuleID());
416
417 for (int i = 0; i != NUMBER_OF_PMTS; ++i){
418 rate_Hz[i] = sum.getRate(i) * RTU;
419 }
420
421 } else {
422
423 FATAL("Summary frame of module " << frame->getModuleID() << " not found at frame index " << timeslice->getFrameIndex() << endl);
424 }
425 }
426
427 const JDAQFrameStatus status = frame->getDAQFrameStatus();
428
429 // Set veto according DAQ or PMT status and rate;
430 // Hits of PMTs with veto will not be counted in livetime nor coincidences nor background
431
432 for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
433 if (!getDAQStatus(status, module, i) ||
434 !getPMTStatus(status, module, i) ||
435 !rateRange_Hz(rate_Hz[i])) {
436 veto[i] = true;
437 }
438 }
439
440 const JHistogram& histogram = zmap[router.getAddress(frame->getModuleID()).first];
441 const JCombinatorics& combinatorics = histogram.getCombinatorics();
442
443 TH2D* h2s = histogram.h2s;
444 TH1D* h1b = histogram.h1b;
445 TH1D* h1L = histogram.h1L;
446
447 //-----------------------------------------------------------
448 // process data
449 //-----------------------------------------------------------
450
451 JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*frame, module);
452
453 buffer.preprocess(option, match);
454
455 for (JSuperFrame2D_t::iterator i = buffer.begin(); i != buffer.end(); ++i) {
456 if (veto[i->getPMTAddress()]) {
457 i->reset();
458 }
459 }
460
461 JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer);
462
463 if (debug >= debug_t) {
464
465 size_t ns = 0;
466
467 for (JSuperFrame2D_t::iterator i = buffer.begin(); i != buffer.end(); ++i) {
468 ns += i->size();
469 }
470
471 cout << "Module: "
472 << setw(10) << frame->getModuleID() << ' '
473 << setw(6) << timeslice->getFrameIndex() << ' '
474 << setw(4) << frame->size() << ' '
475 << setw(4) << ns << ' '
476 << setw(4) << data.size() << endl;
477 }
478
479 // Signal;
480 // Hits from PMTs that do not comply with rate cuts have been exluded above
481
482 size_t numberOfSignalEvents = 0;
483
484 double t1 = numeric_limits<double>::lowest();
485
486 for (vector<hit_type>::const_iterator p = data.begin(); p != data.end(); ) {
487
488 vector<hit_type>::const_iterator q = p;
489
490 while (++q != data.end() && q->getT() - p->getT() <= TMax_ns) {}
491
492 const int N = distance(p,q);
493
494 if (multiplicity(N)) {
495
496 numberOfSignalEvents += 1;
497
498 if (p->getT() > t1 + deadTime_ns) {
499
500 for (vector<hit_type>::const_iterator __p = p; __p != q; ++__p) {
501 for (vector<hit_type>::const_iterator __q = __p; ++__q != q; ) {
502
503 if ((__p->getToT() >= totRange_ns[0]) &&
504 (__q->getToT() >= totRange_ns[0]) &&
505 (__p->getToT() >= totRange_ns[1] ||
506 __q->getToT() >= totRange_ns[1]) &&
507 (__p->getToT() <= totRange_ns[2]) &&
508 (__q->getToT() <= totRange_ns[2])) {
509
510 h2s->Fill((double) combinatorics.getIndex(__p->getPMT(),__q->getPMT()),
511 JCombinatorics::getSign(__p->getPMT(),__q->getPMT()) * (__q->getT() - __p->getT()));
512 }
513 }
514 }
515
516 t1 = p->getT();
517 }
518 }
519
520 p = q;
521 }
522
523 // Background;
524 // Note that rate cuts and veto have already been accounted for when data buffer was filled
525
526 if (background == rndms_t) {
527
528 for (vector<hit_type>::iterator i = data.begin(); i != data.end(); ++i) {
529 *i = hit_type(i->getPMT(), JHit(i->getT() + t0[i->getPMT()], i->getToT()));
530 }
531
532 sort(data.begin(), data.end());
533
534 double t1 = numeric_limits<double>::lowest();
535
536 for (vector<hit_type>::const_iterator p = data.begin(); p != data.end(); ) {
537
538 vector<hit_type>::const_iterator q = p;
539
540 while (++q != data.end() && q->getT() - p->getT() <= TMax_ns) {}
541
542 const int N = distance(p,q);
543
544 if (multiplicity(N)) {
545
546 if (p->getT() > t1 + deadTime_ns) {
547
548 for (vector<hit_type>::const_iterator __p = p; __p != q; ++__p) {
549 for (vector<hit_type>::const_iterator __q = __p; ++__q != q; ) {
550
551 if ((__p->getToT() >= totRange_ns[0]) &&
552 (__q->getToT() >= totRange_ns[0]) &&
553 (__p->getToT() >= totRange_ns[1] ||
554 __q->getToT() >= totRange_ns[1]) &&
555 (__p->getToT() <= totRange_ns[2]) &&
556 (__q->getToT() <= totRange_ns[2])) {
557
558 h1b->Fill((double) combinatorics.getIndex(p->getPMT(),q->getPMT()), 1.0);
559 }
560 }
561 }
562
563 t1 = p->getT();
564 }
565 }
566
567 p = q;
568 }
569
570 } else if (background == rates_t ||
571 background == count_t) {
572
573 double Rs = 0.0; // total rate [Hz]
574
575 for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
576 if (!veto[i]) {
577 Rs += rate_Hz[i]; // [Hz]
578 }
579 }
580
581 for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
582 for (int j = i; ++j != NUMBER_OF_PMTS; ) {
583
584 if (!veto[i] && !veto[j]) {
585
586 const double R1 = rate_Hz[i]; // [Hz]
587 const double R2 = rate_Hz[j]; // [Hz]
588
589 // evaluate expected counts within a time window of 2 * TMax_ns for the two PMTs and the other PMTs inside the optical module;
590 // expected rate due to random coincidences then is product of the probability for the specific observation and the number of trials
591 // the observation is one hit in PMT 1, one hit in PMT 2 and a number of hits in the other PMTs inside the same module
592 // corresponding to the given multiplicity range
593
594 const double N = getP((R1) * 2 * TMax_ns * 1e-9,
595 (R2) * 2 * TMax_ns * 1e-9,
596 (Rs - R1 - R2) * 2 * TMax_ns * 1e-9,
597 multiplicity.getLowerLimit(),
598 multiplicity.getUpperLimit()) * getFrameTime() / (2*TMax_ns);
599
600 h1b->Fill((double) combinatorics.getIndex(i,j), N);
601 }
602 }
603 }
604 }
605
606 //-----------------------------------------------------------
607 // fill the livetime by PMT pairs
608 //-----------------------------------------------------------
609
610 const double livetime_s = getFrameTime()*1e-9 * exp(-deadTime_ns * numberOfSignalEvents/getFrameTime());
611
612 for (size_t i = 0; i < combinatorics.getNumberOfPairs(); ++i) {
613
614 const int pmt1 = combinatorics.getPair(i).first;
615 const int pmt2 = combinatorics.getPair(i).second;
616
617 if (!veto[pmt1] && !veto[pmt2]) {
618 h1L->Fill(i, livetime_s);
619 }
620 }
621 }
622 }
623 }
624 }
625 STATUS(endl);
626
627 if (background == tails_t) {
628
629 const double R = (2.0*TMax_ns) / ((ymax - ymin)/ny);
630
631 for (vector<JHistogram>::iterator hr = zmap.begin() ; hr != zmap.end() ; ++hr) {
632
633 if (hr->is_valid()) {
634
635 TH2D* h2s = hr->h2s;
636 TH1D* h1b = hr->h1b;
637
638 for (int ix = 1; ix <= h1b->GetXaxis()->GetNbins(); ++ix) {
639
640 double Y = 0.0; // sum of counts in tail regions
641 double W = 0.0; // square of error of Y
642 int N = 0; // number of bins in tail regions
643
644 for (int iy = 1; iy <= h2s->GetYaxis()->GetNbins(); ++iy) {
645
646 const double y = h2s->GetYaxis()->GetBinCenter(iy) ;
647
648 if (Tail_ns(fabs(y))) {
649 Y += h2s->GetBinContent(ix,iy);
650 W += h2s->GetBinError (ix,iy) * h2s->GetBinError(ix,iy);
651 N += 1;
652 }
653 }
654
655 h1b->SetBinContent(ix, Y/N * R);
656 h1b->SetBinError (ix, sqrt(W/N) * R);
657 }
658 }
659 }
660 }
661
662 //---------------------------------------------
663 // save normalisation constants and store histograms
664 //---------------------------------------------
665
666 TH1D weights_hist(weights_hist_t, NULL, 3, 0.0, 3.0);
667
668 weights_hist.GetXaxis()->SetBinLabel(1, W1_overall_t); // [s]
669 weights_hist.GetXaxis()->SetBinLabel(2, WS_t); // [ns]
670 weights_hist.GetXaxis()->SetBinLabel(3, WB_t); // [ns]
671
672 weights_hist.Fill(W1_overall_t, counter*getFrameTime()*1e-9); // [s]
673 weights_hist.Fill(WS_t, (ymax - ymin)/ny); // [ns]
674 weights_hist.Fill(WB_t, 2.0*TMax_ns); // [ns]
675
676 out << weights_hist;
677
678 for (vector<JHistogram>::iterator i = zmap.begin(); i != zmap.end(); ++i) {
679 if (i->is_valid()) {
680 out << *(i->h2s);
681 out << *(i->h1b);
682 out << *(i->h1L);
683 }
684 }
685
686 out.Write();
687 out.Close();
688}
int main(int argc, char **argv)
string outputFile
KM3NeT DAQ constants, bit handling, etc.
Data structure for detector geometry and calibration.
Basic data structure for L0 hit.
Tools for handling different hit types.
Match operator for consecutive hits.
General purpose messaging.
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#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:74
ROOT I/O of application specific meta data.
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...
PMT analogue signal processor.
Utility class to parse command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Auxiliaries for pre-processing of hits.
I/O formatting auxiliaries.
#define MAKE_CSTRING(A)
Make C-string.
Definition JPrint.hh:72
Auxiliary class to define a range between two values.
Scanning of objects from a single file according a format that follows from the extension of each fil...
ROOT TTree parameter settings of various packages.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Detector data structure.
Definition JDetector.hh:96
Address of module in detector data structure.
int first
index 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:76
General exception.
Definition JException.hh:24
virtual const char * what() const override
Get error message.
Definition JException.hh:64
Auxiliary class for object identification.
Definition JObjectID.hh:25
int getID() const
Get identifier.
Definition JObjectID.hh:50
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 reading from a list of files.
File router for fast addressing of summary data.
void update(const JDAQHeader &header)
Update router.
bool hasSummaryFrame(const JDAQModuleIdentifier &module) const
Has summary frame.
const JDAQSummaryFrame & getSummaryFrame(const JDAQModuleIdentifier &module) const
Get summary frame.
Auxiliary class to convert pair of indices to unique index and back.
const pair_type & getPair(const int index) const
Get pair of indices for given index.
int getIndex(const int first, const int second) const
Get index of pair of indices.
static int getSign(const int first, const int second)
Sign of pair of indices.
size_t getNumberOfPairs() const
Get number of pairs.
Range of values.
Definition JRange.hh:42
bool is_valid() const
Check validity of range.
Definition JRange.hh:311
T getLowerLimit() const
Get lower limit.
Definition JRange.hh:202
T getUpperLimit() const
Get upper limit.
Definition JRange.hh:213
Reduced data structure for L0 hit.
Definition JHitR0.hh:27
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 getDetectorID() const
Get detector identifier.
int getRunNumber() const
Get run number.
int getFrameIndex() const
Get frame index.
const JDAQFrameStatus & getDAQFrameStatus() const
Get DAQ frame status.
const JDAQHeader & getDAQHeader() const
Get DAQ header.
Definition JDAQHeader.hh:49
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.
#define R1(x)
static const char *const _2S
Name extension for 2D counts.
static const char *const WB_t
Named bin for value of TMax_ns in JCalibrateK40.
double getP(const double E1, const double E2, const double ED, const int M_min, const int M_max)
Get coincidence probability of two PMTs within one module due to random background.
static const char *const weights_hist_t
Histogram naming.
static const char *const _1B
Name extension for 1D background.
static const char *const WS_t
Named bin for time residual bin width.
static const char *const W1_overall_t
Named bin for duration of the run.
static const char *const _1L
Name extension for 1D live time.
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).
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.
bool is_valid(const json &js)
Check validity of JSon data.
JTriggerParameters getTriggerParameters(const JMultipleFileScanner_t &file_list)
Get trigger parameters.
void configure(const T &value, const JAbstractCollection< JAbscissa_t > &bounds, JBool< false > option)
Configuration of value.
int j
Definition JPolint.hh:801
bool getDAQStatus(const JDAQFrameStatus &frame, const JStatus &status)
Test status of DAQ.
bool getPMTStatus(const JStatus &status)
Test status of PMT.
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
Acoustics hit.
Transmission with position.
Definition JBillabong.cc:70
Auxiliary class to sort pairs of PMT addresses within optical module.
virtual double getNPE(const double tot_ns) const override
Get number of photo-electrons.
double getIntegralOfChargeProbability(const double xmin, const double xmax, const int NPE) const
Get integral of probability.
Auxiliary data structure for streaming of STL containers.
Auxiliary class to select ROOT class based on class name.
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
Auxiliary base class for list of file names.
Simple data structure for histogram binning.
int getNumberOfBins() const
Get number of bins.
Auxiliary class for specifying the way of pre-processing of hits.
static std::vector< JPreprocessor > getOptions()
Get options.
@ filter_t
filter consecutive hits according match criterion