Jpp  pmt_effective_area_update_2
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JRipple.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <fstream>
4 #include <iomanip>
5 #include <vector>
6 
7 #include "TROOT.h"
8 #include "TFile.h"
9 #include "TH1D.h"
10 #include "TH2F.h"
11 #include "TF1.h"
12 
14 #include "JDAQ/JDAQEvaluator.hh"
16 
17 #include "JDetector/JDetector.hh"
19 
20 #include "Jeep/JPrint.hh"
21 #include "Jeep/JParser.hh"
22 #include "Jeep/JMessage.hh"
23 
25 
26 #include "JTrigger/JHitR0.hh"
27 
28 #include "JTrigger/JMatchL0.hh"
29 
30 #include "JTrigger/JBuildL0.hh"
31 
33 
37 #include "JTools/JCombinatorics.hh"
38 
39 #include "JROOT/JManager.hh"
41 #include "JMath/JMathToolkit.hh"
42 
44 
45 #include "JSupernova/JSupernova.hh"
47 #include "JSupport/JTreeScanner.hh"
49 #include "JSupport/JSupport.hh"
51 
52 
53 #include "JRipple.hh"
54 
55 
56 using namespace std;
57 using namespace JPP;
58 using namespace KM3NETDAQ;
59 using namespace JSUPERNOVA;
60 
61 /**
62  * \author mlincett
63  * \file
64  * Example application to examine hit rates (and active channels) as a function of time.
65  */
66 int main(int argc, char **argv)
67 {
68 
69  JMultipleFileScanner<> inputFile;
70  JLimit_t& numberOfEvents = inputFile.getLimit();
71  string outputFile;
72  string detectorFile;
73  double TMax_ns;
74  JROOTClassSelector selector;
75  int debug;
76  int binWidth_ms;
77  bool backVeto;
78  bool globalOutputOnly;
79  bool mode2D;
80  JRange<int> multiplicityRange;
81  JDAQHitToTSelector totSelector_ns;
82  // Check input parameters
83 
84  try {
85 
86  JParser<> zap("Example program to examine rates as a function of time on ms-level timescales.");
87 
88  zap['f'] = make_field(inputFile);
89  zap['o'] = make_field(outputFile) = "ripple.root";
90  zap['n'] = make_field(numberOfEvents) = JLimit::max();
91  zap['a'] = make_field(detectorFile);
92  zap['C'] = make_field(selector) = getROOTClassSelection<JDAQTimesliceTypes_t>();
93  zap['T'] = make_field(TMax_ns, "Time window for local coincidences (if 0 run in L0 mode)") = 0.0;
94  zap['B'] = make_field(binWidth_ms, "Bin width (experimental)") = 1;
95  zap['V'] = make_field(backVeto, "Apply retroactive veto.");
96  zap['D'] = make_field(mode2D, "L1 mode: create 2D histogram with time differences of coincidences (heavy memory usage, ignored if TMax_ns = 0)");
97  zap['O'] = make_field(globalOutputOnly, "Write only aggregate histograms");
98  zap['M'] = make_field(multiplicityRange, "L1 mode: multiplicity range (ignored if TMax_ns = 0)") = JRange<int>(2,31);
99  zap['t'] = make_field(totSelector_ns, "ToT acceptance range") = JDAQHitToTSelector(3, 255);
100  zap['d'] = make_field(debug) = 1;
101 
102 
103  zap(argc, argv);
104  }
105  catch(const exception &error) {
106  FATAL(error.what() << endl);
107  }
108 
109  cout.tie(&cerr);
110 
111  if (((int)(getFrameTime() / 1e6)) % binWidth_ms) {
112  FATAL("Frame time must be an integer multiple of bin width");
113  }
114 
116 
117  try {
118  load(detectorFile, detector);
119  }
120  catch(const JException& error) {
121  FATAL(error);
122  }
123 
124  // Configure input streams
125 
127 
129 
130  JTreeScannerInterface<JDAQTimeslice>* pts = zmap[selector];
131 
132  pts->configure(inputFile);
133 
134  int fEnd = pts->rbegin()->getFrameIndex();
135  int fStart = pts->begin( )->getFrameIndex();
136 
137  if (fEnd > inputFile.getUpperLimit()) {
138  fEnd = fStart + inputFile.getUpperLimit();
139  }
140 
141  double tEnd_ms = fEnd * getFrameTime() / 1.0e6;
142  double tStart_ms = (fStart - 1) * getFrameTime() / 1.0e6;
143 
144  int runNumber = pts->begin()->getRunNumber();
145 
146  TString runTag = Form("%d" , runNumber);
147 
148  double tRun_ms = tEnd_ms - tStart_ms;
149 
150  NOTICE("START/END/DURATION [s]: " << tStart_ms / 1000 << " " << tEnd_ms / 1000 << " " << tRun_ms / 1000 << endl);
151 
152  // hit containers and data routers
153 
154  typedef JHitR0 hit_type;
155  typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
156  typedef JSuperFrame1D<hit_type> JSuperFrame1D_t;
157 
158  const JMatchL0<hit_type> match(TMax_ns); // time window self-coincidences
159 
160  const JModuleRouter moduleRouter(detector);
161  const JDAQHitRouter hitRouter(detector);
162 
163  JSummaryRouter summaryRouter;
164 
165  // output histograms
166 
167  typedef JManager<string, TH1F> JManager_t;
168 
169  int nx = (tEnd_ms - tStart_ms) / binWidth_ms;
170 
171  JManager_t RT_DOM(new TH1F("RT_%", NULL, nx , tStart_ms, tEnd_ms));
172  JManager_t NC_DOM(new TH1F("NC_%", NULL, nx / 100, tStart_ms, tEnd_ms));
173 
174  TString rt_tag = runTag + TString(".RT_DET_%");
175  TString nc_tag = runTag + TString(".NC_DET_%");
176 
177  JManager_t RT_DET(new TH1F(rt_tag, NULL, nx, tStart_ms, tEnd_ms));
178  JManager_t NC_DET(new TH1F(nc_tag, NULL, nx / 100, tStart_ms, tEnd_ms));
179 
180  h2d_t* hT = NULL;
181 
182 
183  if (mode2D) {
184 
185  const int ny = 50;
186 
187  const int max_size = floor(h2d_limit / (sizeof(h2d_bintype) * ny));
188 
189  if (nx >= max_size) {
190  FATAL("2D histogram size not supported by ROOT file output; limit input size (-n) below " << floor(max_size / 100.0) << endl);
191  }
192 
193  TString rt2d_tag = runTag + TString(".RT2D_DET");
194 
195  hT = new h2d_t(rt2d_tag, NULL, nx, tStart_ms, tEnd_ms, ny, -TMax_ns, +TMax_ns);
196  hT->SetDirectory(0);
197 
198  }
199 
200  //---------------------
201  // Timeslice processing
202  //---------------------
203 
204  counter_type counter = 0;
205 
206  // preload first timeslice
207 
208  JDAQTimeslice curr;
209 
210  if (pts->hasNext()) {
211  curr = *(pts->next());
212  ++counter;
213  } else {
214  FATAL("Input file is too short.");
215  }
216 
217  // loop over timeslices
218 
219  for ( ; pts->hasNext() && counter != inputFile.getLimit(); ++counter) {
220 
221  STATUS("timeslice: " << setw(10) << counter << '\r'); DEBUG(endl);
222 
223  JDAQTimeslice next = *(pts->next());
224 
225  const int ic = curr.getFrameIndex();
226 
227  const int in = next.getFrameIndex();
228 
229  const double tTimeslice_ms = getTimeOfRTS(ic) / 1.0e6;
230 
231  // look up next summaryslice
232 
233  JDAQSummaryslice* nextSummary = NULL;
234 
235  if (backVeto && ((in - ic) == 1)) {
236  nextSummary = new JDAQSummaryslice(next);
237  summaryRouter.update(nextSummary);
238  }
239 
240  // loop over frames
241 
242  for (JDAQTimeslice::const_iterator frame = curr.begin(); frame != curr.end(); ++frame) {
243 
244  const int moduleID = frame->getModuleID();
245  const JModule& module = moduleRouter.getModule(moduleID);
246  const string moduleLabel = getLabel(module);
247 
248  TH1F* RD = RT_DOM[moduleLabel];
249 
250  vector<bool> veto(NUMBER_OF_PMTS, false);
251 
252  // current veto
253 
254  if (frame->testHighRateVeto() || frame->testFIFOStatus()) {
255  for (int pmt = 0; pmt < NUMBER_OF_PMTS; pmt++) {
256  veto[pmt] = ( frame->testHighRateVeto(pmt) || frame->testFIFOStatus(pmt) );
257  }
258  }
259 
260  // retroactive veto
261 
262  if (nextSummary != NULL) {
263 
264  JDAQSummaryFrame nextFrame = summaryRouter.getSummaryFrame(moduleID);
265 
266  if (nextFrame.testHighRateVeto() || nextFrame.testFIFOStatus()) {
267  for (int pmt = 0; pmt < NUMBER_OF_PMTS; pmt++) {
268  veto[pmt] = veto[pmt] || ( nextFrame.testHighRateVeto(pmt) || nextFrame.testFIFOStatus(pmt) );
269  }
270  }
271 
272  }
273 
274  // track number of active channel
275 
276  NC_DOM[moduleLabel]->Fill(tTimeslice_ms, count(veto.begin(), veto.end(), false));
277 
278  // veto
279 
280  JSuperFrame2D_t& buffer2D = JSuperFrame2D_t::demultiplex(*frame, module, totSelector_ns);
281 
282  for (JSuperFrame2D_t::iterator i = buffer2D.begin(); i != buffer2D.end(); ++i) {
283  if (veto[i->getPMTAddress()]) {
284  i->reset();
285  }
286  }
287 
288  buffer2D.preprocess(JPreprocessor::join_t, match);
289 
290  JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer2D); data.pop_back();
291 
292 
293  // data process
294 
295  if (data.size() > 1) {
296 
297  // sort for L1 processing
298 
299  if (TMax_ns > 0) {
300  sort(data.begin(), data.end());
301  }
302 
303  // loop over hits in frame
304 
305  for (vector<hit_type>::const_iterator p = data.begin(); p != data.end(); p++) {
306 
307  const double tHit_ms = tTimeslice_ms + (p->getT() / 1.0e6);
308 
309  if (TMax_ns == 0.0) { // L0 mode
310 
311  RD->Fill(tHit_ms);
312 
313  } else { // L1 mode
314 
316 
317  while (++q != data.end() && q->getT() - p->getT() <= TMax_ns) {}
318 
319  int M = distance(p, q);
320 
321  if (multiplicityRange(M)) {
322 
323  RD->Fill(tHit_ms);
324 
325  if (mode2D) { // L1 2D mode
326 
327  const double W = factorial(M, 2);
328 
329  for (vector<hit_type>::const_iterator __p = p; __p != q; ++__p) {
330 
331  for (vector<hit_type>::const_iterator __q = __p; ++__q != q; ) {
332 
333  double dt = JCombinatorics::getSign(__p->getPMT(), __q->getPMT()) * (__q->getT() - __p->getT());
334 
335  hT->Fill(tHit_ms, dt, 1.0 / W);
336 
337  }
338  }
339  }
340 
341  }
342 
343  p = (q - 1);
344 
345  }
346 
347  }
348  }
349  }
350 
351  // cleanup
352 
353  if (nextSummary != NULL) { delete nextSummary; }
354 
355  // update pointers
356 
357  curr = next;
358 
359  }
360 
361 
362  //--------------------------------------
363  // Histogram processing
364  //--------------------------------------
365 
366  NOTICE("Processing histograms." << endl);
367 
368  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
369 
370  string moduleLabel = getLabel(*module);
371 
372  if (RT_DOM.count(moduleLabel) && NC_DOM.count(moduleLabel)) {
373 
374  TH1F* rt = RT_DOM.at(moduleLabel);
375  TH1F* nc = NC_DOM.at(moduleLabel);
376 
377  for (int b = 1; b <= rt->GetXaxis()->GetNbins(); b++) {
378 
379  double r = rt->GetBinContent(b);
380  double t = rt->GetBinCenter( b);
381 
382  RT_DET["SUM"]->Fill(t, r);
383 
384  }
385 
386 
387  for (int b = 1; b <= nc->GetXaxis()->GetNbins(); b++) {
388 
389  double n = nc->GetBinContent(b);
390  double t = nc->GetBinCenter( b);
391 
392  NC_DET["SUM"]->Fill(t, n);
393 
394  }
395 
396  } else {
397 
398  DEBUG(moduleLabel << " not active." << endl);
399 
400  }
401 
402  }
403 
404  NOTICE("Writing output file" << endl);
405 
406  if (outputFile != "") {
407 
408  TFile out(outputFile.c_str(), "RECREATE");
409 
410  NOTICE("Writing 1D histograms" << endl);
411 
412  RT_DET.Write(out);
413  NC_DET.Write(out);
414 
415  NOTICE("Writing 2D histogram" << endl);
416 
417  if (hT != NULL) {
418  hT->Write();
419  }
420 
421  if (!globalOutputOnly) {
422 
423  TString dir_tag = runTag + TString(".Modules");
424 
425  NOTICE("Writing individual modules histograms" << endl);
426 
427  TDirectory* dir = out.mkdir(dir_tag);
428  RT_DOM.Write(*dir);
429  NC_DOM.Write(*dir);
430  }
431 
432  NOTICE("Closing file" << endl);
433 
434  out.Close();
435  }
436 
437 }
438 
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
void update(const JDAQSummaryslice *ps)
Update router.
int main(int argc, char *argv[])
Definition: Main.cc:15
TH2F h2d_t
Definition: JRipple.hh:7
ROOT TTree parameter settings of various packages.
do $JPP JMEstimator M
Definition: JMEstimator.sh:37
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:66
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:81
Auxiliary class to select ROOT class based on class name.
Router for direct addressing of module data in detector data structure.
const JDAQSummaryFrame & getSummaryFrame(const JDAQModuleIdentifier &module) const
Get summary frame.
Long64_t counter_type
Type definition for counter.
data_type r[M+1]
Definition: JPolint.hh:742
Dynamic ROOT object management.
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:91
Data structure for detector geometry and calibration.
long long int factorial(const long long int n)
Determine factorial.
Definition: JMathToolkit.hh:42
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.
Float_t h2d_bintype
Definition: JRipple.hh:6
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
I/O formatting auxiliaries.
Reduced data structure for L0 hit.
Definition: JHitR0.hh:25
#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.
#define NOTICE(A)
Definition: JMessage.hh:64
Support methods.
Data time slice.
int debug
debug level
Definition: JSirene.cc:63
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.
Auxiliary class to select JTreeScanner based on ROOT class name.
#define FATAL(A)
Definition: JMessage.hh:67
Scanning of objects from multiple files according a format that follows from the extension of each fi...
bool testHighRateVeto() const
Test high-rate veto status.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
std::vector< int > count
Definition: JAlgorithm.hh:180
Utility class to parse command line options.
alias put_queue eval echo n
Definition: qlib.csh:19
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
2-dimensional frame with time calibrated data from one optical module.
Auxiliary class to select DAQ hits based on time-over-treshold value.
do set_variable DETECTOR_TXT $WORKDIR detector
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
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:40
const unsigned int h2d_limit
Definition: JRipple.hh:11
bool testFIFOStatus() const
Test FIFO status.