Jpp  18.2.1-ARCA-DF-PATCH
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 
15 #include "JDAQ/JDAQTimesliceIO.hh"
17 #include "JDAQ/JDAQEvaluator.hh"
20 #include "JTrigger/JHitR0.hh"
21 #include "JTrigger/JMatchL0.hh"
22 #include "JTrigger/JHitToolkit.hh"
27 #include "JDetector/JDetector.hh"
34 #include "JSupport/JSupport.hh"
35 #include "JSupport/JMeta.hh"
37 #include "JTools/JRange.hh"
39 
41 #include "JROOT/JRootToolkit.hh"
42 #include "JROOT/JRootFileWriter.hh"
45 
46 #include "Jeep/JPrint.hh"
47 #include "Jeep/JParser.hh"
48 #include "Jeep/JMessage.hh"
49 
50 namespace {
51 
52  using JLANG::JObjectID;
55  using JDETECTOR::JModule;
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  * Get coincidence probability of two PMTs within one module due to random background.
136  *
137  * \param E1 expected counts of 1st PMT
138  * \param E2 expected counts of 2nd PMT
139  * \param ED expected counts of module, excluding PMT 1 and PMT 2
140  * \param M_min multiplicity lower limit
141  * \param M_max multiplicity upper limit
142  * \return probability
143  */
144  inline double getP(const double E1, const double E2, const double ED, const int M_min, const int M_max)
145  {
146  // 1st, always calculate the probability for twofold multiplicity
147 
148  double P = E1 * E2 * exp(-(E1 + E2 + ED));
149 
150  // 2nd, calculate the probability of multiplicity at lower limit
151 
152  int MD = 1;
153 
154  for ( ; MD + 2 <= M_min; ++MD) {
155  P *= ED/MD;
156  }
157 
158  // 3rd, add the probabilities of multiplicities up to upper limit
159 
160  double P_sum = P;
161 
162  for ( ; MD + 2 <= M_max; ++MD) {
163  P *= ED/MD;
164  P_sum += P;
165  }
166 
167  return P_sum;
168  }
169 }
170 
171 
172 /**
173  * \file
174  *
175  * Auxiliary program to determine PMT parameters from K40 data.
176  *
177  * By default, the combinatorial background is estimated from the singles rate.\n
178  * In case L0 data are taken (and no summary data are available),
179  * the singles rate can be determined from the amount of L0 data.\n
180  * One can also estimate the background from the tails in the time residual distribution.\n
181  * In that case, option -B can be used to specify the minimal and maximal time offset in ns.\n
182  * For example -B "10 20".\n
183  * Note that the minimal and maximal time should be larger that the width of
184  * the time residual distribution and less than the time window of the L1 coincidence, respectively.\n
185  * The time window is applied symmetrically around a time offset of zero.\n
186  * Finally, if L0 data are available, one can estimate the background using
187  * the same procedure but with a random (read wrong) time calibration.
188  *
189  * The option -t can be used to specify three time-over-threshold values.\n
190  * 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
191  * Furthermore, the time-over-threshold of one of the two hits should be greater or equals the second value.
192  *
193  * The option -M can be used to specify a range of multiplicities.\n
194  * For example -M "2 2" will strictly select two-fold coincidences.
195  * Note that such a selection can bias the result.
196  *
197  * Also consult <a href="https://common.pages.km3net.de/jpp/JMonitor.PDF">documentation</a>.
198  * \author mdejong
199  */
200 int main(int argc, char **argv)
201 {
202  using namespace std;
203  using namespace JPP;
204  using namespace KM3NETDAQ;
205 
206  JMultipleFileScanner_t inputFile;
207  counter_type numberOfEvents;
208  string outputFile;
209  string detectorFile;
210  Double_t TMax_ns;
211  JRange<double> rateRange_Hz;
212  array <double, 3> totRange_ns;
213  JRange<double> Tail_ns;
214  JRange<int> multiplicity;
215  double deadTime_us;
216  JROOTClassSelector selector;
217  string background;
218  JPreprocessor option;
219  int debug;
220 
221  try {
222 
223  JParser<> zap("Auxiliary program to determine PMT parameters from K40 data.");
224 
225  zap['f'] = make_field(inputFile, "input file.");
226  zap['o'] = make_field(outputFile, "output file.") = "calibrate_k40.root";
227  zap['n'] = make_field(numberOfEvents) = JLimit::max();
228  zap['a'] = make_field(detectorFile, "detector file.");
229  zap['T'] = make_field(TMax_ns, "time window [ns].") = 20.0;
230  zap['V'] = make_field(rateRange_Hz, "PMT rate range [Hz].") = JRange<double>(0.0, 20.0e3);
231  zap['t'] = make_field(totRange_ns, "PMT time-over-threshold range [ns].") = { 4.0, 7.0, ToTmax_ns };
232  zap['b'] = make_field(background, "background estimation method.") = rates_t, count_t, tails_t, rndms_t;
233  zap['B'] = make_field(Tail_ns, "time window used for background estimation.") = JRange<double>(15.0, 20.0);
234  zap['M'] = make_field(multiplicity, "multiplicity range of hits on DOM.") = JRange<int>(2, 31);
235  zap['D'] = make_field(deadTime_us, "L1 dead time (us)") = 0.0;
236  zap['C'] = make_field(selector, "timeslice selector, e.g. JDAQTimesliceL1.") = getROOTClassSelection<JDAQTimesliceTypes_t>();
237  zap['O'] = make_field(option, "hit pre-processing option.") = JPreprocessor::getOptions(JPreprocessor::filter_t);
238  zap['d'] = make_field(debug, "debug flag.") = 1;
239 
240  zap(argc, argv);
241  }
242  catch(const exception &error) {
243  FATAL(error.what() << endl);
244  }
245 
246  //-----------------------------------------------------------
247  // check the input parameters
248  //-----------------------------------------------------------
249 
250 
251  if (selector == JDAQTimesliceL2::Class_Name() ||
252  selector == JDAQTimesliceSN::Class_Name()) {
253  FATAL("Option -C <selector> " << selector << " not compatible with calibration method." << endl);
254  }
255 
256  if (selector == JDAQTimeslice ::Class_Name() ||
257  selector == JDAQTimesliceL1::Class_Name()) {
258 
260 
261  try {
262  parameters = getTriggerParameters(inputFile);
263  }
264  catch(const JException& error) {
265  FATAL("No trigger parameters from input:" << error.what() << endl);
266  }
267 
268  if ((selector == JDAQTimeslice ::Class_Name() && parameters.writeL1.prescale > 0) ||
269  (selector == JDAQTimesliceL1::Class_Name())) {
270 
271  if (parameters.TMaxLocal_ns < TMax_ns) {
272  FATAL("Option -T <TMax_ns> = " << TMax_ns << " is larger than in the trigger " << parameters.TMaxLocal_ns << endl);
273  }
274 
275  if (background == rndms_t ||
276  background == count_t) {
277  FATAL("Option -C <selector> " << selector << " incompatible with option -b <background> " << background << endl);
278  }
279  }
280  }
281 
282  if (!multiplicity.is_valid() || multiplicity.getLowerLimit() < 2 || multiplicity.getUpperLimit() > NUMBER_OF_PMTS) {
283  FATAL("Invalid option -M <multiplicity> " << multiplicity << endl);
284  }
285 
286  if (totRange_ns[0] > totRange_ns[1] ||
287  totRange_ns[1] > totRange_ns[2]) {
288  FATAL("Invalid option -t <totRange_ns> " << JEEPZ() << totRange_ns << endl);
289  }
290 
291  if (background == tails_t) {
292 
293  if (!Tail_ns.is_valid()) {
294  FATAL("Invalid option -B <Tail_ns> " << Tail_ns << endl);
295  }
296 
297  if (Tail_ns.getUpperLimit() > TMax_ns) {
298  FATAL("Invalid option -B <Tail_ns> " << Tail_ns << "; upper limit larger than option -T <TMax_ns> " << TMax_ns << endl);
299  }
300  }
301 
302  //-----------------------------------------------------------
303  // load detector file
304  //-----------------------------------------------------------
305 
307 
308  try {
309  load(detectorFile, detector);
310  }
311  catch(const JException& error) {
312  FATAL(error);
313  }
314 
315  if (detector.empty()) {
316  FATAL("Empty detector." << endl);
317  }
318 
319  const JModuleRouter router(detector);
320 
321  //-----------------------------------------------------------
322  // determine time offsets for background method rndms_t
323  //-----------------------------------------------------------
324 
325  vector<double> t0; // time offsets for random background evaluation [ns]
326 
327  for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
328  t0.push_back(i * 2 * TMax_ns);
329  }
330 
331  //-----------------------------------------------------------
332  // correction factor for rate due to specified time-over-threshold range
333  //-----------------------------------------------------------
334 
335  const JPMTAnalogueSignalProcessor cpu;
336 
337  const int NPE = 1;
338  const double RTU = (cpu.getIntegralOfChargeProbability(cpu.getNPE(totRange_ns[0]), cpu.getNPE(totRange_ns[2]), NPE)
339  /
340  cpu.getIntegralOfChargeProbability(cpu.getNPE(0.0), cpu.getNPE(ToTmax_ns), NPE));
341 
342 
343  //-----------------------------------------------------------
344  // initialise histograms
345  //-----------------------------------------------------------
346 
347  vector<JHistogram> zmap(detector.size());
348 
349  const Double_t ymin = -floor(TMax_ns) + 0.5;
350  const Double_t ymax = +floor(TMax_ns) - 0.5;
351  const Int_t ny = (Int_t) ((ymax - ymin) / 1.0);
352 
353  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
354 
355  const JModuleAddress& address = router.getAddress(module->getID());
356 
357  if (!module->empty()) {
358 
359  NOTICE("Booking histograms for module " << module->getID() << endl);
360 
361  zmap[address.first] = JHistogram(*module, JAbstractHistogram<Double_t>(ny, ymin, ymax));
362  }
363  }
364 
365 
366  typedef JHitR0 hit_type;
367  typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
368  typedef JSuperFrame1D<hit_type> JSuperFrame1D_t;
369 
370  const JMatchL0<hit_type> match(TMax_ns); // time window self-coincidences [ns]
371 
372  const double deadTime_ns = deadTime_us * 1e3;
373 
374 
375  TFile out(outputFile.c_str(), "recreate");
376 
377  putObject(out, JMeta(argc, argv));
378 
379  counter_type counter = 0;
380 
381  for (JMultipleFileScanner_t::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
382 
383  STATUS("Processing: " << *i << endl);
384 
385  JSummaryFileRouter summary(*i);
386 
389 
390  for (JDAQHeader header; in.hasNext() && counter != numberOfEvents; ++counter) {
391 
392  STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
393 
394  const JDAQTimeslice* timeslice = in.next();
395 
396  if (header.getDetectorID() != timeslice->getDetectorID() ||
397  header.getRunNumber () != timeslice->getRunNumber ()) {
398 
399  header = timeslice->getDAQHeader();
400 
401  putObject(out, header);
402  }
403 
404  if (background == rates_t) {
405 
406  summary.update(timeslice->getDAQHeader());
407 
408  if (timeslice->getFrameIndex() != summary.getSummaryslice()->getFrameIndex()) {
409 
410  ERROR("Frame indices do not match at "
411  << "[counter = " << counter << "]"
412  << "[timeslice = " << timeslice->getFrameIndex() << "]"
413  << "[summaryslice = " << summary.getSummaryslice()->getFrameIndex() << "]" << endl);
414 
415  continue;
416  }
417  }
418 
419  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
420 
421  if (router.hasModule(frame->getModuleID())) {
422 
423  const JModule& module = router.getModule(frame->getModuleID());
424 
425  //-----------------------------------------------------------
426  // determine background rates and veto per PMT
427  //-----------------------------------------------------------
428 
429  vector<double> rate_Hz(NUMBER_OF_PMTS, 0.0);
430  vector<bool> veto (NUMBER_OF_PMTS, false);
431 
432  if (background == count_t) {
433 
434  for (JDAQSuperFrame::const_iterator i = frame->begin(); i != frame->end(); ++i) {
435  rate_Hz[i->getPMT()] += RTU * 1e9 / getFrameTime();
436  }
437 
438  } else if (background == tails_t) {
439 
440  // see below
441 
442  } else if (background == rndms_t) {
443 
444  // see below
445 
446  } else if (background == rates_t) {
447 
448  if (summary.hasSummaryFrame(frame->getModuleID())) {
449 
450  const JDAQSummaryFrame& sum = summary.getSummaryFrame(frame->getModuleID());
451 
452  for (int i = 0; i != NUMBER_OF_PMTS; ++i){
453  rate_Hz[i] = sum.getRate (i) *RTU;
454  veto [i] = sum.getValue(i) == 0;
455  }
456 
457  } else {
458 
459  FATAL("Summary frame of module " << frame->getModuleID() << " not found at frame index " << timeslice->getFrameIndex() << endl);
460  }
461  }
462 
463  const JDAQFrameStatus status = frame->getDAQFrameStatus();
464 
465  // Set veto according DAQ or PMT status and rate;
466  // Hits of PMTs with veto will not be counted in livetime nor coincidences nor background
467 
468  for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
469  if (!getDAQStatus(status, module, i) ||
470  !getPMTStatus(status, module, i) ||
471  !rateRange_Hz(rate_Hz[i])) {
472  veto[i] = true;
473  }
474  }
475 
476  const JHistogram& histogram = zmap[router.getAddress(frame->getModuleID()).first];
477  const JCombinatorics& combinatorics = histogram.getCombinatorics();
478 
479  TH2D* h2s = histogram.h2s;
480  TH1D* h1b = histogram.h1b;
481  TH1D* h1L = histogram.h1L;
482 
483  //-----------------------------------------------------------
484  // process data
485  //-----------------------------------------------------------
486 
487  JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*frame, module);
488 
489  buffer.preprocess(option, match);
490 
491  for (JSuperFrame2D_t::iterator i = buffer.begin(); i != buffer.end(); ++i) {
492  if (veto[i->getPMTAddress()]) {
493  i->reset();
494  }
495  }
496 
497  JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer);
498 
499  DEBUG("Number of hits " << timeslice->getFrameIndex() << ":" << frame->getModuleID() << ' ' << frame->size() << ' ' << data.size() << endl);
500 
501  // Signal;
502  // Hits from PMTs that do not comply with rate cuts have been exluded above
503 
504  size_t numberOfSignalEvents = 0;
505 
506  double t1 = numeric_limits<double>::lowest();
507 
508  for (vector<hit_type>::const_iterator p = data.begin(); p != data.end(); ) {
509 
511 
512  while (++q != data.end() && q->getT() - p->getT() <= TMax_ns) {}
513 
514  const int N = distance(p,q);
515 
516  if (multiplicity(N)) {
517 
518  numberOfSignalEvents += 1;
519 
520  if (p->getT() > t1 + deadTime_ns) {
521 
522  for (vector<hit_type>::const_iterator __p = p; __p != q; ++__p) {
523  for (vector<hit_type>::const_iterator __q = __p; ++__q != q; ) {
524 
525  if ((__p->getToT() >= totRange_ns[0]) &&
526  (__q->getToT() >= totRange_ns[0]) &&
527  (__p->getToT() >= totRange_ns[1] ||
528  __q->getToT() >= totRange_ns[1]) &&
529  (__p->getToT() <= totRange_ns[2]) &&
530  (__q->getToT() <= totRange_ns[2])) {
531 
532  h2s->Fill((double) combinatorics.getIndex(__p->getPMT(),__q->getPMT()),
533  JCombinatorics::getSign(__p->getPMT(),__q->getPMT()) * (__q->getT() - __p->getT()));
534  }
535  }
536  }
537 
538  t1 = p->getT();
539  }
540  }
541 
542  p = q;
543  }
544 
545  // Background;
546  // Note that rate cuts and veto have already been accounted for when data buffer was filled
547 
548  if (background == rndms_t) {
549 
550  for (vector<hit_type>::iterator i = data.begin(); i != data.end(); ++i) {
551  *i = hit_type(i->getPMT(), JHit(i->getT() + t0[i->getPMT()], i->getToT()));
552  }
553 
554  sort(data.begin(), data.end());
555 
556  double t1 = numeric_limits<double>::lowest();
557 
558  for (vector<hit_type>::const_iterator p = data.begin(); p != data.end(); ) {
559 
561 
562  while (++q != data.end() && q->getT() - p->getT() <= TMax_ns) {}
563 
564  const int N = distance(p,q);
565 
566  if (multiplicity(N)) {
567 
568  if (p->getT() > t1 + deadTime_ns) {
569 
570  for (vector<hit_type>::const_iterator __p = p; __p != q; ++__p) {
571  for (vector<hit_type>::const_iterator __q = __p; ++__q != q; ) {
572 
573  if ((__p->getToT() >= totRange_ns[0]) &&
574  (__q->getToT() >= totRange_ns[0]) &&
575  (__p->getToT() >= totRange_ns[1] ||
576  __q->getToT() >= totRange_ns[1]) &&
577  (__p->getToT() <= totRange_ns[2]) &&
578  (__q->getToT() <= totRange_ns[2])) {
579 
580  h1b->Fill((double) combinatorics.getIndex(p->getPMT(),q->getPMT()), 1.0);
581  }
582  }
583  }
584 
585  t1 = p->getT();
586  }
587  }
588 
589  p = q;
590  }
591 
592  } else if (background == rates_t ||
593  background == count_t) {
594 
595  double Rs = 0.0; // total rate [Hz]
596 
597  for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
598  if (!veto[i]) {
599  Rs += rate_Hz[i]; // [Hz]
600  }
601  }
602 
603  for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
604  for (int j = i; ++j != NUMBER_OF_PMTS; ) {
605 
606  if (!veto[i] && !veto[j]) {
607 
608  const double R1 = rate_Hz[i]; // [Hz]
609  const double R2 = rate_Hz[j]; // [Hz]
610 
611  // evaluate expected counts within a time window of 2 * TMax_ns for the two PMTs and the other PMTs inside the optical module;
612  // expected rate due to random coincidences then is product of the probability for the specific observation and the number of trials
613  // 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
614  // corresponding to the given multiplicity range
615 
616  const double N = getP((R1) * 2 * TMax_ns * 1e-9,
617  (R2) * 2 * TMax_ns * 1e-9,
618  (Rs - R1 - R2) * 2 * TMax_ns * 1e-9,
619  multiplicity.getLowerLimit(),
620  multiplicity.getUpperLimit()) * getFrameTime() / (2*TMax_ns);
621 
622  h1b->Fill((double) combinatorics.getIndex(i,j), N);
623  }
624  }
625  }
626  }
627 
628  //-----------------------------------------------------------
629  // fill the livetime by PMT pairs
630  //-----------------------------------------------------------
631 
632  const double livetime_s = getFrameTime()*1e-9 * exp(-deadTime_ns * numberOfSignalEvents/getFrameTime());
633 
634  for (size_t i = 0; i < combinatorics.getNumberOfPairs(); ++i) {
635 
636  const int pmt1 = combinatorics.getPair(i).first;
637  const int pmt2 = combinatorics.getPair(i).second;
638 
639  if (!veto[pmt1] && !veto[pmt2]) {
640  h1L->Fill(i, livetime_s);
641  }
642  }
643  }
644  }
645  }
646  }
647  STATUS(endl);
648 
649  if (background == tails_t ) {
650 
651  const double R = (2.0*TMax_ns) / ((ymax - ymin)/ny);
652 
653  for (vector<JHistogram>::iterator hr = zmap.begin() ; hr != zmap.end() ; ++hr) {
654 
655  if (hr->is_valid()) {
656 
657  TH2D* h2s = hr->h2s;
658  TH1D* h1b = hr->h1b;
659 
660  for (int ix = 1; ix <= h1b->GetXaxis()->GetNbins(); ++ix) {
661 
662  double Y = 0.0; // sum of counts in tail regions
663  double W = 0.0; // square of error of Y
664  int N = 0; // number of bins in tail regions
665 
666  for (int iy = 1; iy <= h2s->GetYaxis()->GetNbins(); ++iy) {
667 
668  const double y = h2s->GetYaxis()->GetBinCenter(iy) ;
669 
670  if (Tail_ns(fabs(y))) {
671  Y += h2s->GetBinContent(ix,iy);
672  W += h2s->GetBinError (ix,iy) * h2s->GetBinError(ix,iy);
673  N += 1;
674  }
675  }
676 
677  h1b->SetBinContent(ix, Y/N * R);
678  h1b->SetBinError (ix, sqrt(W/N) * R);
679  }
680  }
681  }
682  }
683 
684  //---------------------------------------------
685  // save normalisation constants and store histograms
686  //---------------------------------------------
687 
688  TH1D weights_hist(weights_hist_t, NULL, 3, 0.0, 3.0);
689 
690  weights_hist.GetXaxis()->SetBinLabel(1, W1_overall_t); // [s]
691  weights_hist.GetXaxis()->SetBinLabel(2, WS_t); // [ns]
692  weights_hist.GetXaxis()->SetBinLabel(3, WB_t); // [ns]
693 
694  weights_hist.Fill(W1_overall_t, counter*getFrameTime()*1e-9); // [s]
695  weights_hist.Fill(WS_t, (ymax - ymin)/ny); // [ns]
696  weights_hist.Fill(WB_t, 2.0*TMax_ns); // [ns]
697 
698  out << weights_hist;
699 
700  for (vector<JHistogram>::iterator i = zmap.begin(); i != zmap.end(); ++i) {
701  if (i->is_valid()) {
702  out << *(i->h2s);
703  out << *(i->h1b);
704  out << *(i->h1L);
705  }
706  }
707 
708  out.Write();
709  out.Close();
710 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
const double xmax
Definition: JQuadrature.cc:24
Utility class to parse command line options.
Definition: JParser.hh:1514
Auxiliaries for pre-processing of hits.
General exception.
Definition: JException.hh:24
Auxiliary class to convert pair of indices to unique index and back.
JAbscissa_t getLowerLimit() const
Get lower limit.
Definition: JRange.hh:202
int main(int argc, char *argv[])
Definition: Main.cc:15
ROOT TTree parameter settings of various packages.
double getRate(const int tdc, const double factor=1.0) const
Get count rate.
Basic data structure for L0 hit.
const JModule & getModule(const JObjectID &id) const
Get module parameters.
Data structure for a composite optical module.
Definition: JModule.hh:68
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
int getDetectorID() const
Get detector identifier.
static const char *const WB_t
Named bin for value of TMax_ns in JCalibrateK40.
L0 match criterion.
Definition: JMatchL0.hh:27
#define STATUS(A)
Definition: JMessage.hh:63
bool is_valid(const json &js)
Check validity of JSon data.
void update(const JDAQHeader &header)
Update router.
Detector data structure.
Definition: JDetector.hh:89
Auxiliary class to select ROOT class based on class name.
Router for direct addressing of module data in detector data structure.
#define R1(x)
*fatal Wrong number of arguments esac JCookie sh typeset Z DETECTOR typeset Z SOURCE_RUN typeset Z TARGET_RUN set_variable PARAMETERS_FILE $WORKDIR parameters
Definition: diff-Tuna.sh:38
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:136
Long64_t counter_type
Type definition for counter.
then warning Cannot perform comparison test for histogram
then fatal Wrong number of arguments fi set_variable STRING $argv[1] set_variable DETECTORXY_TXT $WORKDIR $DETECTORXY_TXT tail read X Y CHI2 RMS printf optimum n $X $Y $CHI2 $RMS awk v Y
const JDAQSummaryFrame & getSummaryFrame() const
Get default summary frame.
Auxiliary class for multiplexing object iterators.
Simple data structure for histogram binning.
string outputFile
int getNumberOfBins() const
Get number of bins.
int getRunNumber() const
Get run number.
Data structure for detector geometry and calibration.
Tools for handling different hit types.
Acoustics hit.
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
int getFrameIndex() const
Get frame index.
1-dimensional frame with time calibrated data from one optical module.
Q JDAQEvent livetime_s
Definition: JDataQuality.sh:57
int first
index of module in detector data structure
Scanning of objects from a single file according a format that follows from the extension of each fil...
static const char *const WS_t
Named bin for time residual bin width.
I/O formatting auxiliaries.
Auxiliary class to sort pairs of PMT addresses within optical module.
Detector file.
Definition: JHead.hh:226
Reduced data structure for L0 hit.
Definition: JHitR0.hh:25
const JDAQSummaryslice * getSummaryslice() const
Get summary slice.
Hit data structure.
Definition: JDAQHit.hh:34
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1989
static const char *const _1B
Name extension for 1D background.
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
Data storage class for rate measurements of all PMTs in one module.
int getID() const
Get identifier.
Definition: JObjectID.hh:50
ROOT I/O of application specific meta data.
JAbscissa_t getUpperLimit() const
Get upper limit.
Definition: JRange.hh:213
#define NOTICE(A)
Definition: JMessage.hh:64
#define ERROR(A)
Definition: JMessage.hh:66
Data time slice.
File router for fast addressing of summary data.
virtual const pointer_type & next() override
Get next element.
Address of module in detector data structure.
Match operator for consecutive hits.
virtual bool hasNext() override
Check availability of next element.
bool putObject(TDirectory &dir, const TObject &object)
Write object to ROOT directory.
Auxiliary data structure for streaming of STL containers.
Definition: JPrint.hh:65
General purpose messaging.
double getP(const double expval, bool hit)
Get Poisson probability to observe a hit or not for given expectation value for the number of hits...
Definition: JFitToolkit.hh:41
#define FATAL(A)
Definition: JMessage.hh:67
Scanning of objects from multiple files according a format that follows from the extension of each fi...
Direct access to module in detector data structure.
static const char *const weights_hist_t
Histogram naming.
static const char *const _1L
Name extension for 1D live time.
then JCookie sh JDataQuality D $DETECTOR_ID R
Definition: JDataQuality.sh:41
const double xmin
Definition: JQuadrature.cc:23
Auxiliary base class for list of file names.
virtual const char * what() const override
Get error message.
Definition: JException.hh:64
bool getPMTStatus(const JStatus &status)
Test status of PMT.
then usage $script< input file >[option[primary[working directory]]] nWhere option can be N
Definition: JMuonPostfit.sh:40
void configure(const T &value, const JAbstractCollection< JAbscissa_t > &bounds, JBool< false > option)
Configuration of value.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
Auxiliary class to define a range between two values.
Utility class to parse command line options.
bool hasSummaryFrame(const JDAQModuleIdentifier &module) const
Has summary frame.
const JDAQHeader & getDAQHeader() const
Get DAQ header.
Definition: JDAQHeader.hh:49
Auxiliary class for object identification.
Definition: JObjectID.hh:22
static const char *const W1_overall_t
Named bin for duration of the run.
const JModuleAddress & getAddress(const JObjectID &id) const
Get address of module.
bool hasModule(const JObjectID &id) const
Has module.
2-dimensional frame with time calibrated data from one optical module.
PMT analogue signal processor.
Object reading from a list of files.
int j
Definition: JPolint.hh:792
do set_variable DETECTOR_TXT $WORKDIR detector
bool getDAQStatus(const JDAQFrameStatus &frame, const JStatus &status)
Test status of DAQ.
KM3NeT DAQ constants, bit handling, etc.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
Definition: JCanberra.sh:46
Auxiliary class for specifying the way of pre-processing of hits.
const JDAQFrameStatus & getDAQFrameStatus() const
Get DAQ frame status.
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable STRING $argv[2] set_array QUANTILES set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable MODULE `getModule -a $DETECTOR -L "$STRING 0"` source JAcousticsToolkit.sh typeset -A TRIPODS get_tripods $WORKDIR/tripod.txt TRIPODS XMEAN
JTriggerParameters getTriggerParameters(const JMultipleFileScanner_t &file_list)
Get trigger parameters.
int debug
debug level
double getValue(const int tdc) const
Get value.
then $DIR JPlotNPE PDG P
Definition: JPlotNPE-PDG.sh:62
static const char *const _2S
Name extension for 2D counts.
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
Template definition of histogram object interface.
Definition: JHistogram.hh:27