Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 "JDAQ/JDAQClock.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JDAQ/JDAQSummaryslice.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/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

int main ( int  argc,
char **  argv 
)

Definition at line 61 of file JRipple.cc.

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