Jpp - the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JMonitorMultiplicity.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <string>
3 
4 #include "TFile.h"
5 #include "TH1D.h"
6 #include "TH2D.h"
7 #include "TMath.h"
8 #include "TROOT.h"
9 #include "TVectorD.h"
10 
11 #include "JAAnet/JHead.hh"
14 #include "JDAQ/JDAQEvaluator.hh"
16 #include "JDAQ/JDAQTimesliceIO.hh"
17 #include "JDetector/JDetector.hh"
23 #include "Jeep/JParser.hh"
24 #include "Jeep/JMessage.hh"
25 #include "JMath/JMathToolkit.hh"
26 #include "Jeep/JTimer.hh"
28 #include "JLang/JString.hh"
29 #include "JROOT/JManager.hh"
30 #include "JGizmo/JGizmoToolkit.hh"
32 #include "JSupport/JMeta.hh"
35 #include "JSupport/JSupport.hh"
36 #include "JSupport/JTreeScanner.hh"
38 #include "JTrigger/JHit.hh"
39 #include "JTrigger/JHitR0.hh"
40 #include "JTrigger/JMatchL0.hh"
44 #include "JTools/JRange.hh"
45 
46 /**
47  * \file
48  *
49  * Auxiliary program to analyse coincidence multiplicity;
50  * \author mlincetto, mkarel
51  */
52 int main(int argc, char **argv)
53 {
54  using namespace std;
55  using namespace KM3NETDAQ;
56  using namespace JPP;
57 
58  //-----------------------------------------------------------
59  // parameter interface
60  //-----------------------------------------------------------
61 
63  JLimit_t& numberOfEvents = inputFile.getLimit();
64  string outputFile;
65  string detectorFile;
66  Double_t TMax_ns;
67  JRange<double> rateVeto_Hz;
68  JRange<double> totVeto_ns;
69  JROOTClassSelector selector;
70  int debug;
71  int filterLevel;
72  int muteChannel;
73  bool twoDim;
74  bool monitorOccupancy;
75  bool notJoin;
76 
77  try {
78 
79  JParser<> zap("Auxiliary program to estimate PMT and hit multiplicities.");
80 
81  zap['f'] = make_field(inputFile);
82  zap['o'] = make_field(outputFile) = "monitormultiplicity.root";
83  zap['a'] = make_field(detectorFile);
84  zap['n'] = make_field(numberOfEvents) = JLimit::max();
85  zap['T'] = make_field(TMax_ns) = 10.0; // [ns]
86  zap['V'] = make_field(rateVeto_Hz) = JRange<double>(0, 1e5); // [Hz]
87  zap['t'] = make_field(totVeto_ns) = JRange<double>(0, 1e4); // min and max value of ToT of hits to be monitored [ns]
88  zap['C'] = make_field(selector) = getROOTClassSelection<JDAQTimesliceTypes_t>();
89  zap['d'] = make_field(debug) = 1;
90  zap['M'] = make_field(muteChannel) = -1;
91  zap['F'] = make_field(filterLevel, "0 = raw data; 1 = filtered data; 2+ = only clean frames") = 1;
92  zap['D'] = make_field(twoDim, "2D mode for background subtraction");
93  zap['O'] = make_field(monitorOccupancy, "produces 2D PMT vs multiplicity plots");
94  zap['j'] = make_field(notJoin, "do not join consecutive hits on same PMT (diagnostic purpose)");
95 
96  zap(argc, argv);
97  }
98  catch(const exception &error) {
99  FATAL(error.what() << endl);
100  }
101 
102  cout.tie(&cerr);
103 
104  //-----------------------------------------------------------
105  // check the input parameters
106  //-----------------------------------------------------------
107 
109  bool hasTriggerParameters = true;
110 
111  try {
112  parameters = getTriggerParameters(inputFile);
113  NOTICE("Get trigger parameters from input." << endl);
114  }
115  catch(const JException& error) {
116  ERROR("No trigger parameters from input." << endl);
117  hasTriggerParameters = false;
118 
119  }
120 
121  // keep track of prescale factor for proper MC normalization.
122  int prescale = 1;
123 
124  if (hasTriggerParameters) {
125 
126  // for L1 (most common for sea analysis) check that time window is not larger than in trigger
127 
128  if (parameters.writeL1.prescale > 0 && (selector == "JDAQTimeslice" || selector == "JDAQTimesliceL1") &&
129  parameters.TMaxLocal_ns < TMax_ns) {
130  FATAL("TMax_ns (-T) " << TMax_ns << " ns is larger than in the trigger " << parameters.TMaxLocal_ns << " ns." << endl);
131  }
132 
133  if (selector == "JDAQTimeslice") {
134  // Jpp v8 legacy compatibility
135  prescale = (parameters.writeL1.prescale != 0) ? parameters.writeL1.prescale : parameters.writeL0.prescale;
136  }
137 
138  if (selector == "JDAQTimesliceL0") {
139  prescale = parameters.writeL0.prescale;
140  }
141 
142  if (selector == "JDAQTimesliceL1") {
143  prescale = parameters.writeL1.prescale;
144  }
145 
146  if (selector == "JDAQTimesliceSN") {
147  prescale = parameters.writeSN.prescale;
148  }
149 
150  }
151 
152  if (prescale == 0) {
153  FATAL("[R] According to trigger parameters, the " << selector << " stream should be empty.");
154  } else {
155  NOTICE("[R] Prescale factor is " << prescale << endl);
156  }
157 
158  if (!totVeto_ns.is_valid()) {
159  FATAL("Invalid time over threshold " << totVeto_ns << endl);
160  }
161 
162  //-----------------------------------------------------------
163  // load the detector file
164  //-----------------------------------------------------------
165 
167 
168  try {
169  load(detectorFile, detector);
170  }
171  catch(const JException& error) {
172  FATAL(error);
173  }
174 
175  if (detector.empty())
176  FATAL("Empty detector." << endl);
177 
178  const JModuleRouter router(detector);
179  JSummaryRouter summaryRouter;
180 
181  int detSize = getNumberOfModules(detector);
182 
183  int detID = detector.getID();
184 
185  JDetectorAddressMap& detectorAddressMap = getDetectorAddressMap(detID);
186 
187 
188  //-----------------------------------------------------------
189  // summary of input parameters
190  //-----------------------------------------------------------
191 
192  NOTICE("[R] Frame time [ms] " << getFrameTime() * 1e-6 << endl);
193  NOTICE("[R] Tmax [ns] " << TMax_ns << endl);
194  NOTICE("[R] Rate range [Hz] " << rateVeto_Hz << endl);
195  NOTICE("[R] ToT range [ns] " << totVeto_ns << endl);
196  NOTICE("[R] Timeslice stream " << selector << endl);
197  NOTICE("[R] Detector file has " << detSize << " DOMs." << endl);
198 
199  //-----------------------------------------------------------
200  // initialise histograms
201  //-----------------------------------------------------------
202 
203  // nominal livetime + effective livetime for each DOM
204  TH1D* h_livetime = new TH1D("LiveTime", "Livetime", 1 + detSize, 0.0, 1.0 + detSize);
205 
206  int ibin = 1;
207  h_livetime->GetXaxis()->SetBinLabel(ibin++, "Nominal");
208  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
209  h_livetime->GetXaxis()->SetBinLabel(ibin++, TString(getLabel(*module)));
210  }
211 
212  h_livetime->GetYaxis()->SetTitle("Livetime (s)");
213 
214 
215  // multiplicity histograms
216  const int nx = 1 + NUMBER_OF_PMTS;
217  const double xmin = -0.5;
218  const double xmax = nx - 0.5;
219 
220  typedef JManager<TString, TH1D> JManager_t;
221  typedef JManager<TString, TH2D> JManager2D_t;
222  JManager_t H(new TH1D("%_H", NULL, nx, xmin, xmax)); H->Sumw2();
223  JManager_t P(new TH1D("%_P", NULL, nx, xmin, xmax)); P->Sumw2();
224 
225  JManager2D_t HT(new TH2D("%_HT", NULL, nx, xmin, xmax, 100, -TMax_ns, +TMax_ns));
226 
227  TH2D* CO_proto = new TH2D("%_CO", NULL, NUMBER_OF_PMTS, 0.5, 0.5 + NUMBER_OF_PMTS, NUMBER_OF_PMTS, -0.5, -0.5 + NUMBER_OF_PMTS);
228 
229  if (monitorOccupancy) {
230  JModuleAddressMap memo = getDetectorAddressMap<JKM3NeT_t>().getDefaultModuleAddressMap();
231  setAxisLabels(CO_proto->GetYaxis(), memo);
232  }
233 
234  JManager2D_t CO(CO_proto);
235 
236  // benchmark histogram for time distribution of high-multiplicities
237  const int sn_mul_min = 6;
238  TH1D* time_distr = new TH1D("T_distr", "Coincidence (M > 6) distribution in timeslice", 100, 0, 1e8);
239 
240  //-----------------------------------------------------------
241  // output file and auxiliary data structures
242  //-----------------------------------------------------------
243 
244  // output file
245  TFile out(outputFile.c_str(), "recreate");
246 
247  // auxiliary counters
248 
249  int count_HRV = 0; // high-rate-veto counter
250  int count_FAF = 0; // fifo-almost-full counter
251  int count_URV = 0; // user-rate-veto counter
252 
253  int treeMismatchCount = 0;
254 
255  //-----------------------------------------------------------
256  // file pre-processing
257  //-----------------------------------------------------------
258 
259  JTreeScanner<JDAQSummaryslice, JDAQEvaluator> summaryFile(inputFile);
260 
262 
263  vector<JHitR0> filteredData;
264 
265  typedef JSuperFrame2D<JHitR0> JSuperFrame2D_t;
266  typedef JSuperFrame1D<JHitR0> JSuperFrame1D_t;
267  const JMatchL0<JHitR0> match(TMax_ns);
268 
269  counter_type counter = 0;
270 
271  // process Monte Carlo aanet header
272  bool hasMCHeader = true ;
273  bool isMupage = false;
274  double liveTimeMC = 0 ;
275  double liveTimeMCErr = 0 ;
276 
277  JHead mc_header;
278 
279  try {
280  mc_header = (JHead)getHeader(inputFile);
281  }
282  catch (const JException& error) {
283  hasMCHeader = false;
284  }
285 
286  NOTICE("[R] File source is " << (hasMCHeader ? "MC" : "DAQ") << "." << endl);
287  if (hasMCHeader) {
288  liveTimeMC = mc_header.livetime.numberOfSeconds;
289  liveTimeMCErr = mc_header.livetime.errorOfSeconds;
290  if (liveTimeMC > 0) {
291  isMupage = true;
292  NOTICE("[R] mupage live time is (" << liveTimeMC << " +/- " << liveTimeMCErr << ") s." << endl);
293  } else {
294  NOTICE("[R] aanet header found but mupage live time is zero. Likely not a mupage file.");
295  }
296  }
297 
298  // retrieving and exporting run number to tree metadata;
299 
300  int runNumber = 0;
301 
302  if (summaryFile.hasNext()) {
303  runNumber = summaryFile.getEntry(0)->getRunNumber();
304  NOTICE("[R] Processing run #" << runNumber << endl);
305  summaryFile.rewind();
306  }
307 
308  //-----------------------------------------------------------
309  // file processing: loop over timeslices
310  //-----------------------------------------------------------
311 
312  for ( ; in.hasNext() && counter != inputFile.getLimit(); ++counter) {
313 
314  STATUS("Entry: " << setw(10) << counter << '\r');
315 
316  const JDAQTimeslice* timeslice = in.next();
317  const Long64_t index = summaryFile.find(*timeslice);
318  JDAQSummaryslice* summary = (index != -1 ? summaryFile.getEntry(index) : NULL);
319 
320  summaryRouter.update(summary);
321 
322  // check correspondence between summary data and timeslice data
323  if (summary != NULL) {
324  if (timeslice->getFrameIndex() != summary->getFrameIndex()) {
325  ++treeMismatchCount;
326  ERROR("[R>T] Frame indices do not match [counter=" << counter << ", timeslice=" << timeslice->getFrameIndex() << ", summaryslice=" << summary->getFrameIndex() << "]" << endl);
327  if (treeMismatchCount > 1) {
328  FATAL("ROOT TTrees not aligned. Run #" << runNumber << endl);
329  }
330  continue;
331  }
332  }
333 
334  //-----------------------------------------------------------
335  // timeslice processing: loop over superframes
336  //-----------------------------------------------------------
337  for (JDAQTimeslice::const_iterator super_frame = timeslice->begin(); super_frame != timeslice->end(); ++super_frame) {
338 
339  // empty superframes are not skipped because SN timeslices could have empty superframes even if the module was active;
340  // if DAQ error or FIFO full -> frame counts as dead time;
341  // if frame is missing but no errors -> frame counts as live time (SN stream has many empty frames);
342 
343  int moduleID = super_frame->getModuleID();
344 
345  if (!router.hasModule(moduleID)) {
346  // module is not mapped in detector file, nothing can be done
347  continue;
348  }
349 
350  JModuleAddressMap moduleAddressMap = detectorAddressMap.get(moduleID);
351 
352  //----------------------------------------------------------
353  // retrieval of summary data from superframe or summaryframe
354  //----------------------------------------------------------
355 
356  vector<double> rate_Hz(NUMBER_OF_PMTS, 0.0);
357  JDAQFrameStatus status;
358 
359  if (summary == NULL) {
360  WARNING("[R>T>F] Missing summary for timeslice." << endl);
361  // PMT rates estimated from timeslice data
362  for (JDAQSuperFrame::const_iterator i = super_frame->begin(); i != super_frame->end(); ++i) {
363  rate_Hz[i->getPMT()] += 1e9 / getFrameTime();
364  }
365 
366  // frame status retrieved from superframe
367  status = super_frame->getDAQFrameStatus();
368 
369  } else {
370 
371  // summary frame lookup
372 
373  JDAQSummaryFrame summary_frame;
374 
375  if (!summaryRouter.hasSummaryFrame(super_frame->getModuleIdentifier())) {
376  summaryRouter.print(cout);
377  FATAL("[R>D>F] Summary frame not found, module = " << super_frame->getModuleID() << endl);
378  } else {
379  summary_frame = summaryRouter.getSummaryFrame(super_frame->getModuleID());
380  }
381 
382  // PMT rates recovered from summaryslice data
383  for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
384  rate_Hz[i] = summary_frame.getRate(i);
385  }
386 
387  // frame status retrieved from summaryframe
388  status = summary_frame.getDAQFrameStatus();
389 
390  }
391 
392  //---------------------------
393  // processing of summary data
394  //---------------------------
395 
396  // DAQ error (e.g. missing UDP trailer) -> discard frame
397  if (!status.testDAQStatus()) {
398  WARNING("[R>D>F] DAQ status not okay in timeslice #" << timeslice->getFrameIndex() << ", DOM=" << super_frame->getModuleID() << endl);
399  continue;
400  }
401 
402  vector<bool> veto(NUMBER_OF_PMTS, false);
403  bool moduleLiveTime = 0;
404 
405  // default status of module is alive if at least one PMT is active
406  for (vector<double>::const_iterator r = rate_Hz.begin(); !moduleLiveTime && r != rate_Hz.end(); ++r) {
407  if ((*r) > 0) {
408  moduleLiveTime = 1;
409  }
410  }
411 
412  int frame_HRV_count = status.countHighRateVeto();
413  int frame_FAF_count = status.countFIFOStatus();
414  int frame_URV_count = 0; // user rate veto
415 
416  for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
417  if (veto[i] = (veto[i] || !rateVeto_Hz(rate_Hz[i]))) {
418  frame_URV_count++;
419  }
420  }
421 
422  count_HRV += frame_HRV_count;
423  count_FAF += frame_FAF_count;
424  count_URV += frame_URV_count;
425 
426  int frame_VETO_count = frame_HRV_count + frame_FAF_count + frame_URV_count;
427 
428 
429  if (filterLevel >= 2 && frame_VETO_count > 0) {
430  moduleLiveTime = 0;
431  fill(veto.begin(), veto.end(), true);
432  } else if (filterLevel >= 1 && frame_HRV_count + frame_FAF_count > 0) {
433  for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
434  veto[i] = (veto[i] || (status.testFIFOStatus(i) || status.testHighRateVeto(i)));
435  }
436  }
437 
438  // recover DOM in case the dead channel corresponds to the muted channel
439 
440  if (muteChannel != -1) {
441 
442  // veto conditions are not exclusive
443  int mute_VETO_count =
444  status.testHighRateVeto(muteChannel) +
445  status.testFIFOStatus(muteChannel ) +
446  !rateVeto_Hz(rate_Hz[muteChannel] );
447 
448  if (mute_VETO_count == frame_VETO_count) {
449  moduleLiveTime = 1;
450  fill(veto.begin(), veto.end(), false);
451  }
452 
453  veto[muteChannel] = true;
454 
455  }
456 
457  // get address and module for current superframe DOM;
458  const JModuleAddress& address = router.getAddress(super_frame->getModuleID());
459  const JModule& module = detector.getModule(address);
460  const TString moduleLabel = getLabel(module);
461 
462  // acknowledge DOM is live in this timeslice for analysis purpose;
463  if (moduleLiveTime) {
464  h_livetime->Fill(moduleLabel, getFrameTime() * 1e-9 * moduleLiveTime);
465  }
466 
467  //--------------------------------------------------------------------
468  // superframe processing: hit joining, calibration, filter, processing
469  //--------------------------------------------------------------------
470 
471  JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*super_frame, router.getModule(super_frame->getModuleID()));
472 
473  if (!notJoin) {
474  for (JSuperFrame2D_t::iterator i = buffer.begin(); i != buffer.end(); ++i) {
475  i->join(match);
476  }
477  }
478 
479  JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer);
480 
481  data.pop_back(); // remove end marker
482 
483  filteredData.clear();
484 
485  for (vector<JHitR0>::const_iterator h = data.begin(); h != data.end(); ++h) {
486 
487  const int pmt = h->getPMT();
488 
489  // veto and calibration;
490  if (!veto[pmt] &&
491  rateVeto_Hz(rate_Hz[pmt]) &&
492  totVeto_ns(h->getToT()) ) {
493  filteredData.push_back(*h);
494  }
495  }
496 
497  // processing signal;
498 
499  if (filteredData.size() > 1) {
500 
501  TH1D* h_hit = H[moduleLabel];
502  TH1D* h_pmt = P[moduleLabel];
503 
504  TH2D* h2_hit = HT[moduleLabel];
505  TH2D* h2_co = CO[moduleLabel];
506 
507  sort(filteredData.begin(), filteredData.end());
508 
509  //-----------------------------------------------------------
510  // superframe processing: loop on calibrated hits
511  //-----------------------------------------------------------
512  for (vector<JHitR0>::const_iterator p = filteredData.begin(); p != filteredData.end(); ) {
513 
515 
516  // coincidence building using a fixed time window
517  while (++q != filteredData.end() && q->getT() - p->getT() <= TMax_ns ) {}
518 
519  // multiplicity estimation
520  int hit_multiplicity = distance(p,q);
521  set<int> pmthit;
522  for (vector<JHitR0>::const_iterator __p = p ; __p != q; ++__p ) { pmthit.insert(__p->getPMT()) ; }
523  int pmt_multiplicity = pmthit.size() ;
524 
525  // histogram filling
526  h_pmt->Fill(pmt_multiplicity);
527  h_hit->Fill(hit_multiplicity);
528 
529  if (twoDim && hit_multiplicity > 1) {
530  const double W = factorial(hit_multiplicity, 2);
531  for (vector<JHitR0>::const_iterator __p = p; __p != q; ++__p) {
532  for (vector<JHitR0>::const_iterator __q = __p; ++__q != q; ) {
533  double dt = JCombinatorics::getSign(__p->getPMT(), __q->getPMT()) * (__q->getT() - __p->getT());
534  h2_hit->Fill(hit_multiplicity, dt, 1.0/W);
535  }
536  }
537  }
538 
539  if (monitorOccupancy) {
540  for (set<int>::const_iterator r = pmthit.begin(); r != pmthit.end(); ++r) {
541  JPMTPhysicalAddress pmtAddress = moduleAddressMap.getAddressTranslator(*r);
542  TString pmtLabel(pmtAddress.toString());
543  h2_co->Fill(pmt_multiplicity, pmtLabel, 1.0/pmt_multiplicity);
544  }
545  }
546 
547  // benchmarking
548  if (hit_multiplicity >= sn_mul_min) {
549  time_distr->Fill(p->getT());
550  }
551 
552  p = q;
553 
554  }
555 
556  }
557 
558  }
559 
560  }
561 
562  NOTICE("[R] " << counter << " timeslices processed." << endl);
563 
564  if (isMupage) {
565  h_livetime->Scale(liveTimeMC / (getFrameTime() * 1e-9 * counter * prescale));
566  }
567 
568  // nominal live time;
569  double nominalLiveTime = (isMupage ? liveTimeMC : counter * getFrameTime() * 1e-9) / prescale;
570  h_livetime->Fill("Nominal", nominalLiveTime);
571 
572  out.cd();
573 
574  int nLiveDOMs = H.size();
575 
576  h_livetime->Write();
577 
578  // writing hit multiplicities is redundant given new options
579  // H.Write(out);
580 
581  P.Write(out);
582 
583  if (twoDim) {
584  HT.Write(out);
585  }
586 
587  if (monitorOccupancy) {
588  CO.Write(out);
589  }
590 
591  time_distr->Write();
592 
593  double rate_HRV = count_HRV / (1.0 * counter);
594  double rate_FAF = count_FAF / (1.0 * counter);
595  double rate_URV = count_URV / (1.0 * counter);
596  double biolumIndex = (rate_HRV + rate_FAF) / (nLiveDOMs * NUMBER_OF_PMTS);
597  NOTICE("[R] Average HRV count per timeslice: " << rate_HRV << endl);
598  NOTICE("[R] Average FIFO-almost-full count per timeslice: " << rate_FAF << endl);
599  NOTICE("[R] Average user rate veto count per timeslice: " << rate_URV << endl);
600  NOTICE("[R] " << H.size() << " DOMs were active in the run." << endl);
601  NOTICE("[R] Bioluminescence proxy index: " << biolumIndex << endl);
602 
603  putObject(&out, JMeta(argc,argv));
604 
605  out.Close();
606 
607  NOTICE(endl);
608 
609  return (counter ? 0 : 1);
610 
611 }
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
General exception.
Definition: JException.hh:23
#define WARNING(A)
Definition: JMessage.hh:65
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.
int countFIFOStatus() const
Count FIFO status.
Auxiliary methods for geometrical methods.
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:57
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
L0 match criterion.
Definition: JMatchL0.hh:27
std::string getLabel(const JLocation &location)
Get module label for monitoring and other applications.
Definition: JLocation.hh:246
#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.
static const double H
Planck constant [eV s].
*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
Lookup table for PMT addresses in detector.
Long64_t counter_type
Type definition for counter.
data_type r[M+1]
Definition: JPolint.hh:742
Dynamic ROOT object management.
Auxiliary class for multiplexing object iterators.
static const JModuleCounter getNumberOfModules
Function object to count unique modules.
Basic data structure for time and time over threshold information of hit.
string outputFile
Data structure for detector geometry and calibration.
long long int factorial(const long long int n)
Determine factorial.
Definition: JMathToolkit.hh:42
int getFrameIndex() const
Get frame index.
1-dimensional frame with time calibrated data from one optical module.
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Lookup table for PMT addresses in optical module.
Detector specific mapping between logical positions and readout channels of PMTs in optical modules...
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
Detector file.
Definition: JHead.hh:196
Hit data structure.
Definition: JDAQHit.hh:34
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
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.
void setAxisLabels(TAxis *axis, const JModuleAddressMap &memo)
Set axis with PMT address labels.
#define NOTICE(A)
Definition: JMessage.hh:64
#define ERROR(A)
Definition: JMessage.hh:66
Data time slice.
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.
Router for fast addressing of summary data in JDAQSummaryslice data structure as a function of the op...
JDetectorAddressMap & getDetectorAddressMap()
Get detector address map.
General purpose messaging.
Monte Carlo run header.
Definition: JHead.hh:1113
#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.
bool testHighRateVeto() const
Test high-rate veto status.
const JPMTAddressTranslator & getAddressTranslator(const int tdc) const
Get PMT address translator.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
Auxiliary class to define a range between two values.
int countHighRateVeto() const
Count high-rate veto status.
General purpose class for object reading from a list of file names.
Utility class to parse command line options.
std::string toString() const
Convert PMT physical address to string.
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.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
Data structure for PMT physical address.
do set_variable DETECTOR_TXT $WORKDIR detector
KM3NeT DAQ constants, bit handling, etc.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
bool testDAQStatus() const
Test DAQ status of packets.
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:38
const JDAQFrameStatus & getDAQFrameStatus() const
Get DAQ frame status.
JTriggerParameters getTriggerParameters(const JMultipleFileScanner_t &file_list)
Get trigger parameters.
then $DIR JPlotNPE PDG P
Definition: JPlotNPE-PDG.sh:60
bool testFIFOStatus() const
Test FIFO status.
bool putObject(TDirectory *dir, const T &object)
Write object to ROOT directory.