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