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