Jpp
Functions
JRipple.cc File Reference
#include <string>
#include <iostream>
#include <fstream>
#include <iomanip>
#include <vector>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2F.h"
#include "TF1.h"
#include "km3net-dataformat/online/JDAQClock.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"
#include "JTrigger/JHitR0.hh"
#include "JTrigger/JMatchL0.hh"
#include "JTrigger/JBuildL0.hh"
#include "JTrigger/JSummaryRouter.hh"
#include "JTrigger/JSuperFrame1D.hh"
#include "JTrigger/JSuperFrame2D.hh"
#include "JTrigger/JTimesliceRouter.hh"
#include "JTools/JCombinatorics.hh"
#include "JGizmo/JManager.hh"
#include "JLang/JObjectMultiplexer.hh"
#include "JMath/JMathToolkit.hh"
#include "JROOT/JROOTClassSelector.hh"
#include "JSupernova/JSupernova.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JSupport/JAutoTreeScanner.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JSupportToolkit.hh"
#include "JRipple.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Author
mlincett

Example application to examine hit rates (and active channels) as a function of time.

Definition in file JRipple.cc.

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 62 of file JRipple.cc.

63 {
64 
65  JMultipleFileScanner<> inputFile;
66  JLimit_t& numberOfEvents = inputFile.getLimit();
67  string outputFile;
68  string detectorFile;
69  double TMax_ns;
70  double TMaxSelf_ns;
71  JROOTClassSelector selector;
72  int debug;
73  int binWidth_ms;
74  bool backVeto;
75  bool globalOutputOnly;
76  bool mode2D;
77  JRange<int> multiplicityRange;
78  JRange<double> totRange_ns;
79  // Check input parameters
80 
81  try {
82 
83  JParser<> zap("Example program to examine rates as a function of time on ms-level timescales.");
84 
85  zap['f'] = make_field(inputFile);
86  zap['o'] = make_field(outputFile) = "ripple.root";
87  zap['n'] = make_field(numberOfEvents) = JLimit::max();
88  zap['a'] = make_field(detectorFile);
89  zap['C'] = make_field(selector) = getROOTClassSelection<JDAQTimesliceTypes_t>();
90  zap['T'] = make_field(TMax_ns, "Time window for local coincidences (if 0 run in L0 mode)") = 0.0;
91  zap['S'] = make_field(TMaxSelf_ns, "Time window for merge of double pulses") = 10.0;
92  zap['B'] = make_field(binWidth_ms, "Bin width (experimental)") = 1;
93  zap['V'] = make_field(backVeto, "Apply retroactive veto.");
94  zap['D'] = make_field(mode2D, "L1 mode: create 2D histogram with time differences of coincidences (heavy memory usage, ignored if TMax_ns = 0)");
95  zap['O'] = make_field(globalOutputOnly, "Write only aggregate histograms");
96  zap['M'] = make_field(multiplicityRange, "L1 mode: multiplicity range (ignored if TMax_ns = 0)") = JRange<int>(2,31);
97  zap['t'] = make_field(totRange_ns, "ToT acceptance range") = JRange<double>(0.0, 255.0);
98  zap['d'] = make_field(debug) = 1;
99 
100 
101  zap(argc, argv);
102  }
103  catch(const exception &error) {
104  FATAL(error.what() << endl);
105  }
106 
107  cout.tie(&cerr);
108 
109  if (((int)(getFrameTime() / 1e6)) % binWidth_ms) {
110  FATAL("Frame time must be an integer multiple of bin width");
111  }
112 
114 
115  try {
116  load(detectorFile, detector);
117  }
118  catch(const JException& error) {
119  FATAL(error);
120  }
121 
122  // Configure input streams
123 
125 
127 
128  JTreeScannerInterface<JDAQTimeslice>* pts = zmap[selector];
129 
130  pts->configure(inputFile);
131 
132  int fEnd = pts->rbegin()->getFrameIndex();
133  int fStart = pts->begin( )->getFrameIndex();
134 
135  if (fEnd > inputFile.getUpperLimit()) {
136  fEnd = fStart + inputFile.getUpperLimit();
137  }
138 
139  double tEnd_ms = fEnd * getFrameTime() / 1.0e6;
140  double tStart_ms = (fStart - 1) * getFrameTime() / 1.0e6;
141 
142  int runNumber = pts->begin()->getRunNumber();
143 
144  TString runTag = Form("%d" , runNumber);
145 
146  double tRun_ms = tEnd_ms - tStart_ms;
147 
148  NOTICE("START/END/DURATION [s]: " << tStart_ms / 1000 << " " << tEnd_ms / 1000 << " " << tRun_ms / 1000 << endl);
149 
150  // Configure data structure and routers
151 
152  typedef JSuperFrame2D<JHitR0> JSuperFrame2D_t;
153 
154  typedef JSuperFrame1D<JHitR0> JSuperFrame1D_t;
155 
156  const JMatchL0<JHitR0> match(TMaxSelf_ns); // time window self-coincidences [ns]
157 
158  const JModuleRouter moduleRouter(detector);
159 
160  const JDAQHitRouter hitRouter(detector);
161 
162  //JTimesliceRouter timesliceRouter;
163  JSummaryRouter summaryRouter;
164 
165  // Configure output histograms
166 
167  typedef JManager<string, TH1F> JManager_t;
168 
169  int nx = (tEnd_ms - tStart_ms) / binWidth_ms;
170 
171  int x_inf = tStart_ms;
172 
173  JManager_t RT_DOM(new TH1F("RT_%", NULL, nx , x_inf, tEnd_ms));
174  JManager_t NC_DOM(new TH1F("NC_%", NULL, nx / 100, x_inf, tEnd_ms));
175 
176  TString rt_tag = runTag + TString(".RT_DET_%");
177  TString nc_tag = runTag + TString(".NC_DET_%");
178 
179  JManager_t RT_DET(new TH1F(rt_tag, NULL, nx, x_inf, tEnd_ms));
180  JManager_t NC_DET(new TH1F(nc_tag, NULL, nx / 100, x_inf, tEnd_ms));
181 
182  h2d_t* hT = NULL;
183 
184 
185  if (mode2D) {
186 
187  const int ny = 50;
188 
189  const int max_size = floor(h2d_limit / (sizeof(h2d_bintype) * ny));
190 
191  if (nx >= max_size) {
192  FATAL("2D histogram size not supported by ROOT file output; limit input size (-n) below " << floor(max_size / 100.0) << endl);
193  }
194 
195  TString rt2d_tag = runTag + TString(".RT2D_DET");
196 
197  hT = new h2d_t(rt2d_tag, NULL, nx, x_inf, tEnd_ms, ny, -TMax_ns, +TMax_ns);
198  hT->SetDirectory(0);
199 
200  }
201 
202  //---------------------
203  // Timeslice processing
204  //---------------------
205 
206  counter_type counter = 0;
207 
208  // preload first timeslice
209 
210  JDAQTimeslice curr;
211 
212  if (pts->hasNext()) {
213  curr = *(pts->next());
214  ++counter;
215  } else {
216  FATAL("Input file is too short.");
217  }
218 
219  // loop over timeslices
220 
221  for ( ; pts->hasNext() && counter != inputFile.getLimit(); ++counter) {
222 
223  STATUS("timeslice: " << setw(10) << counter << '\r'); DEBUG(endl);
224 
225  JDAQTimeslice next = *(pts->next());
226 
227  const int ic = curr.getFrameIndex();
228 
229  const int in = next.getFrameIndex();
230 
231  const double tTimeslice_ms = getTimeOfRTS(ic) / 1.0e6;
232 
233  // look up next summaryslice
234 
235  JDAQSummaryslice* nextSummary = NULL;
236 
237  if (backVeto && ((in - ic) == 1)) {
238  nextSummary = new JDAQSummaryslice(next);
239  summaryRouter.update(nextSummary);
240  }
241 
242  // loop over frames
243 
244  for (JDAQTimeslice::const_iterator frame = curr.begin(); frame != curr.end(); ++frame) {
245 
246  const int DOMID = frame->getModuleID();
247  const JModule& module = moduleRouter.getModule(DOMID);
248  const string moduleLabel = getModuleLabel(module);
249 
250  TH1F* RD = RT_DOM[moduleLabel];
251 
252  vector<bool> veto(NUMBER_OF_PMTS, false);
253 
254  // current veto
255 
256  if (frame->testHighRateVeto() || frame->testFIFOStatus()) {
257  for (int pmt = 0; pmt < NUMBER_OF_PMTS; pmt++) {
258  veto[pmt] = ( frame->testHighRateVeto(pmt) || frame->testFIFOStatus(pmt) );
259  }
260  }
261 
262 
263  // retroactive veto
264 
265  if (nextSummary != NULL) {
266 
267  JDAQSummaryFrame nextFrame = summaryRouter.getSummaryFrame(DOMID);
268 
269  if (nextFrame.testHighRateVeto() || nextFrame.testFIFOStatus()) {
270  for (int pmt = 0; pmt < NUMBER_OF_PMTS; pmt++) {
271  veto[pmt] = veto[pmt] || ( nextFrame.testHighRateVeto(pmt) || nextFrame.testFIFOStatus(pmt) );
272  }
273  }
274 
275  }
276 
277 
278  // determine number of active channels
279 
280  int nActiveChannels = count(veto.begin(), veto.end(), false);
281 
282  NC_DOM[moduleLabel]->Fill(tTimeslice_ms, nActiveChannels);
283 
284 
285  // merge double pulses
286 
287  JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*frame, module);
288 
289  for (JSuperFrame2D_t::iterator i = buffer.begin(); i != buffer.end(); ++i) {
290  if (veto[i->getPMTAddress()]) {
291  i->reset();
292  } else {
293  i->join(match);
294  }
295  }
296 
297  JSuperFrame1D_t& buffer1D = JSuperFrame1D_t::multiplex(buffer);
298  buffer1D.pop_back();
299 
300  // process hit data
301 
302  vector<JHitR0> data;
303 
304  for (vector<JHitR0>::const_iterator h = buffer1D.begin(); h != buffer1D.end(); h++) {
305  if (totRange_ns(h->getToT())) {
306  data.push_back(*h);
307  }
308  }
309 
310  if (data.size() > 1) {
311 
312  // data need to be sorted for L1 processing
313 
314  if (TMax_ns > 0) {
315  sort(data.begin(), data.end());
316  }
317 
318 
319  // loop over hits in frame
320 
321  for (vector<JHitR0>::const_iterator p = data.begin(); p != data.end(); p++) {
322 
323  const double tHit = p->getT();
324  const double tHit_ms = tTimeslice_ms + (tHit / 1.0e6);
325 
326  if (TMax_ns == 0.0) {
327 
328  // L0 mode
329 
330  RD->Fill(tHit_ms);
331 
332  } else {
333 
334  // L1 mode
335 
337 
338  while (++q != data.end() && q->getT() - p->getT() <= TMax_ns) {}
339 
340  int M = distance(p, q);
341 
342  if (multiplicityRange(M)) {
343 
344  RD->Fill(tHit_ms);
345 
346  if (mode2D) {
347 
348  // L1 2D mode
349 
350  const double W = factorial(M, 2);
351 
352  for (vector<JHitR0>::const_iterator __p = p; __p != q; ++__p) {
353 
354  for (vector<JHitR0>::const_iterator __q = __p; ++__q != q; ) {
355 
356  double dt = JCombinatorics::getSign(__p->getPMT(), __q->getPMT()) * (__q->getT() - __p->getT());
357 
358  hT->Fill(tHit_ms, dt, 1.0 / W);
359 
360  }
361  }
362  }
363 
364  }
365 
366  p = (q - 1);
367 
368  }
369 
370  }
371  }
372  }
373 
374  // update pointers
375 
376  curr = next;
377 
378  if (nextSummary != NULL) { delete nextSummary; }
379 
380  }
381 
382 
383  //--------------------------------------
384  // Histogram processing
385  //--------------------------------------
386 
387  NOTICE("Processing histograms." << endl);
388 
389  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
390 
391  string moduleLabel = getModuleLabel(*module);
392 
393  if (RT_DOM.count(moduleLabel) && NC_DOM.count(moduleLabel)) {
394 
395  TH1F* rt = RT_DOM.at(moduleLabel);
396  TH1F* nc = NC_DOM.at(moduleLabel);
397 
398  for (int b = 1; b <= rt->GetXaxis()->GetNbins(); b++) {
399 
400  double r = rt->GetBinContent(b);
401  double t = rt->GetBinCenter( b);
402 
403  RT_DET["SUM"]->Fill(t, r);
404 
405  }
406 
407 
408  for (int b = 1; b <= nc->GetXaxis()->GetNbins(); b++) {
409 
410  double n = nc->GetBinContent(b);
411  double t = nc->GetBinCenter(b );
412 
413  NC_DET["SUM"]->Fill(t, n);
414 
415  }
416 
417  } else {
418 
419  DEBUG(moduleLabel << " not active." << endl);
420 
421  }
422 
423  }
424 
425  NOTICE("Writing output file" << endl);
426 
427  if (outputFile != "") {
428 
429  TFile out(outputFile.c_str(), "RECREATE");
430 
431  NOTICE("Writing 1D histograms" << endl);
432 
433  RT_DET.Write(out);
434  NC_DET.Write(out);
435 
436  NOTICE("Writing 2D histogram" << endl);
437 
438  if (hT != NULL) {
439  hT->Write();
440  }
441 
442  if (!globalOutputOnly) {
443 
444  TString dir_tag = runTag + TString(".Modules");
445 
446  NOTICE("Writing individual modules histograms" << endl);
447 
448  TDirectory* dir = out.mkdir(dir_tag);
449  RT_DOM.Write(*dir);
450  NC_DOM.Write(*dir);
451  }
452 
453  NOTICE("Closing file" << endl);
454 
455  out.Close();
456  }
457 
458 }
JSUPPORT::JAutoTreeScanner
Auxiliary class to select JTreeScanner based on ROOT class name.
Definition: JAutoTreeScanner.hh:34
KM3NETDAQ::JDAQFrameStatus::testHighRateVeto
bool testHighRateVeto() const
Test high-rate veto status.
Definition: JDAQFrameStatus.hh:201
JDETECTOR::getModuleLabel
std::string getModuleLabel(const JModuleLocation &location)
Get module label for monitoring and other applications.
Definition: JDetectorToolkit.hh:708
JMATH::factorial
long long int factorial(const long long int n)
Determine factorial.
Definition: JMathToolkit.hh:41
JSUPPORT::JLimit
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
JLANG::JType
Auxiliary class for a type holder.
Definition: JType.hh:19
JSUPERNOVA::h2d_limit
const unsigned int h2d_limit
Definition: JRipple.hh:11
JDETECTOR::load
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
Definition: JDetectorToolkit.hh:476
JTRIGGER::JSuperFrame1D
1-dimensional frame with time calibrated data from one optical module.
Definition: JSuperFrame1D.hh:35
JTOOLS::n
const int n
Definition: JPolint.hh:628
std::vector
Definition: JSTDTypes.hh:12
KM3NETDAQ::JDAQChronometer::getFrameIndex
int getFrameIndex() const
Get frame index.
Definition: JDAQChronometer.hh:132
JTOOLS::JRange< int >
JTRIGGER::JSummaryRouter::getSummaryFrame
const JDAQSummaryFrame & getSummaryFrame(const JDAQModuleIdentifier &module) const
Get summary frame.
Definition: JSummaryRouter.hh:90
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
JGIZMO::JManager
Auxiliary class to manage set of compatible ROOT objects (e.g.
Definition: JManager.hh:40
NOTICE
#define NOTICE(A)
Definition: JMessage.hh:64
KM3NETDAQ::getFrameTime
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
KM3NETDAQ::JDAQSummaryFrame
Data storage class for rate measurements of all PMTs in one module.
Definition: JDAQSummaryFrame.hh:320
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
Template definition for direct access of elements in ROOT TChain.
Definition: JTreeScanner.hh:91
JSUPPORT::JTreeScannerInterface
Auxiliary interface for direct access of elements in ROOT TChain.
Definition: JTreeScannerInterface.hh:34
KM3NETDAQ::JDAQSummaryslice
Data summary slice.
Definition: JDAQSummaryslice.hh:25
JDETECTOR::JModule
Data structure for a composite optical module.
Definition: JModule.hh:49
KM3NETDAQ::JDAQFrameStatus::testFIFOStatus
bool testFIFOStatus() const
Test FIFO status.
Definition: JDAQFrameStatus.hh:245
JTRIGGER::JSummaryRouter
Router for fast addressing of summary data in JDAQSummaryslice data structure as a function of the op...
Definition: JSummaryRouter.hh:32
STATUS
#define STATUS(A)
Definition: JMessage.hh:63
JDETECTOR::JModuleRouter
Router for direct addressing of module data in detector data structure.
Definition: JModuleRouter.hh:34
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
KM3NETDAQ::getTimeOfRTS
double getTimeOfRTS(const JDAQChronometer &chronometer)
Get time of last RTS in ns since start of run for a given chronometer.
Definition: JDAQChronometer.hh:223
JSUPPORT::JMultipleFileScanner<>
JDETECTOR::JDetector
Detector data structure.
Definition: JDetector.hh:80
JAANET::detector
Detector file.
Definition: JHead.hh:130
JSUPERNOVA::h2d_bintype
Float_t h2d_bintype
Definition: JRipple.hh:6
DEBUG
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
JTRIGGER::JSummaryRouter::update
void update(JDAQSummaryslice *ps)
Update router.
Definition: JSummaryRouter.hh:49
JTOOLS::r
data_type r[M+1]
Definition: JPolint.hh:709
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
JROOT::counter_type
Long64_t counter_type
Type definition for counter.
Definition: JCounter.hh:24
JLANG::JException
General exception.
Definition: JException.hh:23
JSUPERNOVA::h2d_t
TH2F h2d_t
Definition: JRipple.hh:7
JDETECTOR::JDAQHitRouter
Simple wrapper around JModuleRouter class for direct addressing of PMT data in detector data structur...
Definition: JDAQHitRouter.hh:24