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