Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JTurbot.cc
Go to the documentation of this file.
1 
2 
3 #include <string>
4 #include <iostream>
5 #include <iomanip>
6 #include <limits>
7 #include <map>
8 #include <vector>
9 
10 #include "TROOT.h"
11 #include "TFile.h"
12 #include "TH1D.h"
13 #include "TF1.h"
14 #include "TMath.h"
15 
17 #include "JDAQ/JDAQEventIO.hh"
18 #include "JDAQ/JDAQTimesliceIO.hh"
19 #include "JDAQ/JDAQEvaluator.hh"
20 #include "JDetector/JDetector.hh"
24 #include "JTrigger/JBuildL0.hh"
25 #include "JTrigger/JBuildL1.hh"
26 #include "JLang/JSinglePointer.hh"
27 #include "JGizmo/JManager.hh"
28 #include "JTools/JWeight.hh"
32 #include "JSupport/JTreeScanner.hh"
34 #include "JSupport/JSupport.hh"
35 #include "JSupport/JMeta.hh"
36 
37 #include "Jeep/JPrint.hh"
38 #include "Jeep/JParser.hh"
39 #include "Jeep/JMessage.hh"
40 
41 namespace {
42 
43  // status of module
44 
45  enum JStatus_t {
46  DEFAULT = 0, //!< module exists in detector
47  NODATA = -2, //!< module not present in data or no entries
48  OUT_OF_SYNC = -3, //!< module is out-of-sync
49  ERROR = -4, //!< module is out-of-sync, possibly multiple time offsets
50  IN_SYNC = +1 //!< module is in-sync
51  };
52 }
53 
54 /**
55  * \file
56  *
57  * Example program to search for correlations between triggered events and time slice data.\n
58  *
59  * This application can be used to identify optical modules that are out-of-sync with respect to the central White Rabbit clock.\n
60  * To this end, time correlations are searched for between triggered events and local coincidences, referred to as L1 hits.
61  * The local coincidences are independently obtained from time slice data and can thereby cross the time slice border of a triggered event.\n
62  * The QE of the PMTs of the modules that are out-of-sync are set to 0.\n
63  * The number of time slices to look for L1 hits away from the triggered event can be set via command line option <tt>-N <number of time slices></tt>.\n
64  * The larger this value, the bigger time offsets can be covered.
65  * To identify the modules that are out-of-sync as fast as possible (without determining the time offset),
66  * this value can be set to 0 with the caveat that a out-of-sync occurring during the data taking run maybe missed.
67  *
68  * \author mdejong
69  */
70 int main(int argc, char **argv)
71 {
72  using namespace std;
73  using namespace JPP;
74  using namespace KM3NETDAQ;
75 
76  JSingleFileScanner<> inputFile;
77  JLimit_t& numberOfEvents = inputFile.getLimit();
78  string outputFile;
79  string detectorFile;
80  double TMaxLocal_ns;
81  string pmtFile;
82  int numberOfTimeslices;
83  double binWidth_ns;
84  double Pmin;
85  JROOTClassSelector selector;
86  int qaqc;
87  int debug;
88 
89  try {
90 
91  JParser<> zap("Example program to search for correlations between triggered events and timeslice data.");
92 
93  zap['f'] = make_field(inputFile, "input file.");
94  zap['o'] = make_field(outputFile, "optional output file containing ROOT histograms.") = "";
95  zap['a'] = make_field(detectorFile, "detector file.");
96  zap['n'] = make_field(numberOfEvents) = JLimit::max();
97  zap['T'] = make_field(TMaxLocal_ns, "time window for local coincidences (L1),") = 20.0;
98  zap['P'] = make_field(pmtFile, "optional PMT file (to set QE of PMTs).") = "";
99  zap['N'] = make_field(numberOfTimeslices, "time slice difference between triggered event and L1.") = 100;
100  zap['W'] = make_field(binWidth_ns, "bin width for output histograms." ) = 10e3;
101  zap['p'] = make_field(Pmin, "minimal probability for background to be signal.") = 1.0e-7;
102  zap['C'] = make_field(selector, "data type selection (default is all).") = "", getROOTClassSelection<JDAQTimesliceTypes_t>();
103  zap['Q'] = make_field(qaqc) = 0;
104  zap['d'] = make_field(debug) = 2;
105 
106  zap(argc, argv);
107  }
108  catch(const exception& error) {
109  FATAL(error.what() << endl);
110  }
111 
112 
113  cout.tie(&cerr);
114 
116 
117  try {
118  load(detectorFile, detector);
119  }
120  catch(const JException& error) {
121  FATAL(error);
122  }
123 
124  if (detector.empty()) {
125  FATAL("Empty detector " << detectorFile << endl);
126  }
127 
128  const JDAQHitRouter router(detector);
129 
130  const double TMax_ns = max(binWidth_ns, getMaximalTime(detector));
131 
132  typedef double hit_type;
133 
134  JBuildL0<hit_type> buildL0;
136  vector <hit_type> data;
137 
138  typedef JManager<int, TH1D> JManager_t;
139 
140  const Double_t xmin = -(numberOfTimeslices + 1) * getFrameTime() - 0.5*binWidth_ns;
141  const Double_t xmax = +(numberOfTimeslices + 1) * getFrameTime() + 0.5*binWidth_ns;
142  const Int_t nx = (Int_t) ((xmax - xmin) / binWidth_ns);
143 
144  JManager_t manager(new TH1D("M_%", NULL, nx, xmin, xmax));
145 
147 
149 
151 
152  if (selector == "") {
153 
154  if ((ps = new JTreeScanner< JAssertConversion<JDAQTimesliceSN, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
155  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL2, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
156  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL1, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
157  (ps = new JTreeScanner< JAssertConversion<JDAQTimeslice, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
158  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL0, JDAQTimeslice> >(inputFile))->getEntries() != 0) {
159  } else {
160  FATAL("No timeslice data." << endl);
161  }
162 
163  NOTICE("Selected class " << ps->getClassName() << endl);
164 
165  } else {
166 
167  ps = zmap[selector];
168 
169  ps->configure(inputFile);
170  }
171 
172  ps->setLimit(inputFile.getLimit());
173 
174  typedef map<int, vector<double> > map_type; // frame index -> event times
175 
176  map_type buffer;
177 
178  for (JTreeScanner<JDAQEvent> in(inputFile.getFilename()); in.hasNext(); ) {
179 
180  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
181 
182  JDAQEvent* event = in.next();
183 
184  // Average time of triggered hits
185 
186  double t0 = 0.0;
187 
188  for (JDAQEvent::const_iterator<JDAQTriggeredHit> hit = event->begin<JDAQTriggeredHit>(); hit != event->end<JDAQTriggeredHit>(); ++hit) {
189  t0 += getTime(*hit, router.getPMT(*hit));
190  }
191 
192  t0 /= event->size<JDAQTriggeredHit>();
193  t0 += event->getFrameIndex() * getFrameTime();
194 
195  buffer[event->getFrameIndex()].push_back(t0);
196  }
197  STATUS(endl);
198 
199 
200  while (ps->hasNext()) {
201 
202  STATUS("event: " << setw(10) << ps->getCounter() << '\r'); DEBUG(endl);
203 
204  JDAQTimeslice* timeslice = ps->next();
205 
206  map_type::const_iterator p = buffer.lower_bound(timeslice->getFrameIndex() - numberOfTimeslices);
207  map_type::const_iterator q = buffer.upper_bound(timeslice->getFrameIndex() + numberOfTimeslices);
208 
209  int number_of_events = 0;
210 
211  for (map_type::const_iterator i = p; i != q; ++i) {
212  number_of_events += i->second.size();
213  }
214 
215  if (number_of_events == 0) {
216  continue;
217  }
218 
219  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
220 
221  data.clear();
222 
223  buildL1(*frame, router.getModule(frame->getModuleID()), back_inserter(data));
224 
225  TH1D* h1 = manager[frame->getModuleID()];
226 
227  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
228 
229  const double t1 = *hit + frame->getFrameIndex() * getFrameTime();
230 
231  for (map_type::const_iterator i = p; i != q; ++i) {
232  for (map_type::mapped_type::const_iterator j = i->second.begin(); j != i->second.end(); ++j) {
233 
234  const double t0 = *j;
235 
236  h1->Fill(t1 - t0);
237  }
238  }
239  }
240  }
241  }
242  STATUS(endl);
243 
244 
245  // fit function
246 
247  TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2])) + [3]");
248 
249  string option = "L";
250 
251  if (debug < debug_t) {
252  option += "Q";
253  }
254 
255 
256  map<int, JStatus_t> status;
257 
258  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
259 
260  status[module->getID()] = DEFAULT;
261 
262  JManager_t::iterator p = manager.find(module->getID());
263 
264  if (p == manager.end() || p->second->GetEntries() == 0) {
265 
266  status[module->getID()] = NODATA;
267 
268  continue;
269  }
270 
271  TH1D* h1 = p->second;
272 
273  // start values
274 
275  Double_t ymax = 0.0;
276  Double_t x0 = 0.0; // [ns]
277  Double_t sigma = 250.0; // [ns]
278  Double_t ymin = 0.0;
279 
280  for (int i = 1; i != h1->GetNbinsX(); ++i) {
281 
282  const Double_t x = h1->GetBinCenter (i);
283  const Double_t y = h1->GetBinContent(i);
284 
285  if (y > ymax) {
286  ymax = y;
287  x0 = x;
288  }
289 
290  ymin += y;
291  }
292 
293  ymin /= h1->GetNbinsX();
294 
295  f1.SetParameter(0, ymax);
296  f1.SetParameter(1, x0);
297  if (binWidth_ns < sigma)
298  f1.SetParameter(2, sigma);
299  else
300  f1.FixParameter(2, binWidth_ns/sqrt(12.0));
301  f1.SetParameter(3, ymin);
302 
303  // fit
304 
305  h1->Fit(&f1, option.c_str(), "same", x0 - 5 * sigma, x0 + 5 * sigma);
306 
307 
308  if (fabs(f1.GetParameter(1)) <= 0.5*getFrameTime() && // position
309  f1.GetParameter(0) >= f1.GetParameter(3)) { // S/N
310 
311  status[module->getID()] = IN_SYNC;
312  }
313 
314 
315  // look for peaks at time offsets equal to a multiple of frame time
316 
317  map<int, JWeight> bg; // background
318  map<int, JWeight> sn; // signal(s)
319 
320  for (Int_t i = 1, ns = -(numberOfTimeslices + 1); i <= h1->GetNbinsX(); ++i) {
321 
322  const Double_t x = h1->GetBinCenter (i);
323  const Double_t y = h1->GetBinContent(i);
324 
325  while (x > (ns + 1) * getFrameTime() - TMax_ns) { ++ns; }
326 
327  if (fabs(x - ns * getFrameTime()) < 1.0 * TMax_ns)
328  sn[ns].put(y);
329  else if (fabs(x - ns * getFrameTime()) < 10.0 * TMax_ns)
330  bg[ns].put(y);
331  }
332 
333  DEBUG("Module " << setw(8) << module->getID() << endl);
334 
335  for (map<int, JWeight>::iterator i = sn.begin(); i != sn.end(); ) {
336 
337  const Double_t y = bg[i->first].getTotal() * i->second.getCount() / bg[i->first].getCount();
338  const Double_t P = TMath::PoissonI(i->second.getTotal(), y);
339 
340  DEBUG("Peak at " << setw(4) << i->first << " [frame time] "
341  << "S/N = "
342  << FIXED(7,1) << i->second.getTotal() << " / "
343  << FIXED(7,1) << y << ' '
344  << SCIENTIFIC(12,3) << P << endl);
345 
346  if (i->second.getTotal() > y && P < Pmin)
347  ++i; // keep
348  else
349  sn.erase(i++); // remove
350  }
351 
352  if (!(sn.size() == 1 &&
353  sn.begin()->first == 0)) {
354 
355  ERROR("Module " << setw(8) << module->getID() << " number of peaks " << sn.size() << endl);
356 
357  for (map<int, JWeight>::const_iterator i = sn.begin(); i != sn.end(); ++i) {
358  ERROR("Peak at " << setw(4) << i->first << " [frame time] "
359  << "S/N = "
360  << FIXED(7,1) << i->second.getTotal() << " / "
361  << FIXED(7,1) << bg[i->first].getTotal() * i->second.getCount() / bg[i->first].getCount() << endl);
362  }
363 
364  status[module->getID()] = (sn.size() == 1 ? OUT_OF_SYNC : ERROR);
365  }
366  }
367 
368 
369  if (pmtFile != "") {
370 
372 
373  try {
374  parameters.load(pmtFile.c_str());
375  }
376  catch(const JException& error) {}
377 
378  for (map<int, JStatus_t>::const_iterator i = status.begin(); i != status.end(); ++i) {
379 
380  if (i->second != IN_SYNC) {
381 
382  NOTICE("Module " << setw(8) << i->first << " set QEs of all PMTs to zero." << endl);
383 
384  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
385  parameters[JPMTIdentifier(i->first, pmt)].QE = 0.0;
386  }
387  }
388  }
389 
390  ofstream out(pmtFile.c_str());
391 
392  parameters.comment.add(JMeta(argc, argv));
393 
394  out << parameters << endl;
395 
396  out.close();
397  }
398 
399  if (outputFile != "") {
400  manager.Write(outputFile.c_str());
401  }
402 
403  int nin = 0;
404  int nout = 0;
405 
406  for (map<int, JStatus_t>::const_iterator i = status.begin(); i != status.end(); ++i) {
407  if (i->second == IN_SYNC) {
408  ++nin;
409  }
410  if (i->second != IN_SYNC &&
411  i->second != NODATA) {
412  ++nout;
413  }
414  }
415 
416  NOTICE("Number of modules out/in sync " << nout << '/' << nin << endl);
417 
418  QAQC(nin << ' ' << nout << endl);
419 
420  return 0;
421 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:71
Utility class to parse command line options.
Definition: JParser.hh:1493
General exception.
Definition: JException.hh:23
Direct access to PMT data in detector data structure for DAQ hits.
void clear()
Clear data.
debug
Definition: JMessage.hh:29
JStatus_t
Definition: JTurbot.cc:45
ROOT TTree parameter settings.
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:80
Template const_iterator.
Definition: JDAQEvent.hh:68
Auxiliary class to select ROOT class based on class name.
esac print_variable DETECTOR INPUT_FILE OUTPUT_FILE CDF for TYPE in
Definition: JSirene.sh:45
*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
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:69
Auxialiary class to assert type conversion.
Definition: JConversion.hh:66
const JModule & getModule(const JDAQKeyHit &hit) const
Get module parameters.
Dynamic ROOT object management.
Auxiliary class for a type holder.
Definition: JType.hh:19
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:481
double getTime(const Hit &hit)
Get true time of hit.
string outputFile
Template definition for direct access of elements in ROOT TChain.
Definition: JTreeScanner.hh:91
Data structure for detector geometry and calibration.
Simple wrapper around JModuleRouter class for direct addressing of PMT data in detector data structur...
int getFrameIndex() const
Get frame index.
esac $JPP_DIR examples JDetector JTransitTime o $OUTPUT_FILE n N $NPE T $TTS_NS d $DEBUG for HISTOGRAM in tts tt2 pmt
Definition: JTransitTime.sh:36
The template JSinglePointer class can be used to hold a pointer to an object.
Auxiliary class to manage set of compatible ROOT objects (e.g.
Definition: JManager.hh:42
Scanning of objects from a single file according a format that follows from the extension of each fil...
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
I/O formatting auxiliaries.
Detector file.
Definition: JHead.hh:130
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
ROOT I/O of application specific meta data.
#define NOTICE(A)
Definition: JMessage.hh:64
#define ERROR(A)
Definition: JMessage.hh:66
Data time slice.
Auxiliary class for map of PMT parameters.
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
int debug
debug level
Definition: JSirene.cc:61
double getMaximalTime(const double R_Hz)
Get maximal time for given rate.
General purpose messaging.
Auxiliary class to select JTreeScanner based on ROOT class name.
Template L1 hit builder.
Definition: JBuildL1.hh:85
#define FATAL(A)
Definition: JMessage.hh:67
void load(const char *file_name)
Load from input file.
Utility class to parse command line options.
#define QAQC(A)
QA/QC output macro.
Definition: JMessage.hh:100
Auxiliary data structure for L1 build parameters.
Definition: JBuildL1.hh:37
int j
Definition: JPolint.hh:634
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
JComment & add(const std::string &comment)
Add comment.
Definition: JComment.hh:51
Template L0 hit builder.
Definition: JBuildL0.hh:35
KM3NeT DAQ constants, bit handling, etc.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
const JPMT & getPMT(const JDAQKeyHit &hit) const
Get PMT parameters.
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:518
int qaqc
QA/QC file descriptor.
then $DIR JPlotNPE PDG P
Definition: JPlotNPE-PDG.sh:60
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
int main(int argc, char *argv[])
Definition: Main.cpp:15