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