Jpp  15.0.1-rc.1-highQE
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  // status of module
45 
46  enum JStatus_t {
47  DEFAULT = 0, //!< module exists in detector
48  NODATA = -2, //!< module not present in data or no entries
49  OUT_SYNC = -3, //!< module is out-of-sync
50  ERROR = -4, //!< module is out-of-sync, possibly multiple time offsets
51  IN_SYNC = +1 //!< module is in-sync
52  };
53 }
54 
55 /**
56  * \file
57  *
58  * Example program to search for correlations between triggered events and time slice data.\n
59  *
60  * This application can be used to identify optical modules that are out-of-sync with respect to the central White Rabbit clock.\n
61  * To this end, time correlations are searched for between triggered events and local coincidences, referred to as L1 hits.
62  * The local coincidences are independently obtained from time slice data and can thereby cross the time slice border of a triggered event.\n
63  * The QE of the PMTs of the modules that are out-of-sync are set to 0.\n
64  * 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
65  * The larger this value, the bigger time offsets can be covered.
66  * To identify the modules that are out-of-sync as fast as possible (without determining the time offset),
67  * this value can be set to 0 with the caveat that a out-of-sync occurring during the data taking run maybe missed.
68  *
69  * \author mdejong
70  */
71 int main(int argc, char **argv)
72 {
73  using namespace std;
74  using namespace JPP;
75  using namespace KM3NETDAQ;
76 
77  JSingleFileScanner<> inputFile;
78  JLimit_t& numberOfEvents = inputFile.getLimit();
79  string outputFile;
80  string detectorFile;
81  bool overwriteDetector;
82  double TMaxLocal_ns;
83  int numberOfTimeslices;
84  double binWidth_ns;
85  double Pmin;
86  JROOTClassSelector selector;
87  int qaqc;
88  int debug;
89 
90  try {
91 
92  JParser<> zap("Example program to search for correlations between triggered events and timeslice data.");
93 
94  zap['f'] = make_field(inputFile, "input file.");
95  zap['o'] = make_field(outputFile, "optional output file containing ROOT histograms.") = "";
96  zap['a'] = make_field(detectorFile, "detector file.");
97  zap['A'] = make_field(overwriteDetector, "overwrite detector file provided through '-a' with module (PMT) status.");
98  zap['n'] = make_field(numberOfEvents) = JLimit::max();
99  zap['T'] = make_field(TMaxLocal_ns, "time window for local coincidences (L1),") = 20.0;
100  zap['N'] = make_field(numberOfTimeslices, "time slice difference between triggered event and L1.") = 100;
101  zap['W'] = make_field(binWidth_ns, "bin width for output histograms." ) = 10e3;
102  zap['p'] = make_field(Pmin, "minimal probability for background to be signal.") = 1.0e-7;
103  zap['C'] = make_field(selector, "data type selection (default is all).") = "", getROOTClassSelection<JDAQTimesliceTypes_t>();
104  zap['Q'] = make_field(qaqc) = 0;
105  zap['d'] = make_field(debug) = 2;
106 
107  zap(argc, argv);
108  }
109  catch(const exception& error) {
110  FATAL(error.what() << endl);
111  }
112 
113 
114  cout.tie(&cerr);
115 
117 
118  try {
119  load(detectorFile, detector);
120  }
121  catch(const JException& error) {
122  FATAL(error);
123  }
124 
125  if (detector.empty()) {
126  FATAL("Empty detector " << detectorFile << endl);
127  }
128 
129  const JDAQHitRouter router(detector);
130 
131  const double TMax_ns = max(binWidth_ns, getMaximalTime(detector));
132 
133  typedef double hit_type;
134 
135  JBuildL0<hit_type> buildL0;
137  vector <hit_type> data;
138 
139  typedef JManager<int, TH1D> JManager_t;
140 
141  const Double_t xmin = -(numberOfTimeslices + 1) * getFrameTime() - 0.5*binWidth_ns;
142  const Double_t xmax = +(numberOfTimeslices + 1) * getFrameTime() + 0.5*binWidth_ns;
143  const Int_t nx = (Int_t) ((xmax - xmin) / binWidth_ns);
144 
145  JManager_t manager(new TH1D("M_%", NULL, nx, xmin, xmax));
146 
148 
150 
152 
153  if (selector == "") {
154 
155  if ((ps = new JTreeScanner< JAssertConversion<JDAQTimesliceSN, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
156  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL2, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
157  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL1, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
158  (ps = new JTreeScanner< JAssertConversion<JDAQTimeslice, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
159  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL0, JDAQTimeslice> >(inputFile))->getEntries() != 0) {
160  } else {
161  FATAL("No timeslice data." << endl);
162  }
163 
164  NOTICE("Selected class " << ps->getClassName() << endl);
165 
166  } else {
167 
168  ps = zmap[selector];
169 
170  ps->configure(inputFile);
171  }
172 
173  ps->setLimit(inputFile.getLimit());
174 
175  typedef map<int, vector<double> > map_type; // frame index -> event times
176 
177  map_type buffer;
178 
179  for (JTreeScanner<JDAQEvent> in(inputFile.getFilename()); in.hasNext(); ) {
180 
181  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
182 
183  JDAQEvent* event = in.next();
184 
185  // Average time of triggered hits
186 
187  double t0 = 0.0;
188 
189  for (JDAQEvent::const_iterator<JDAQTriggeredHit> hit = event->begin<JDAQTriggeredHit>(); hit != event->end<JDAQTriggeredHit>(); ++hit) {
190  t0 += getTime(*hit, router.getPMT(*hit));
191  }
192 
193  t0 /= event->size<JDAQTriggeredHit>();
194  t0 += event->getFrameIndex() * getFrameTime();
195 
196  buffer[event->getFrameIndex()].push_back(t0);
197  }
198  STATUS(endl);
199 
200 
201  while (ps->hasNext()) {
202 
203  STATUS("event: " << setw(10) << ps->getCounter() << '\r'); DEBUG(endl);
204 
205  JDAQTimeslice* timeslice = ps->next();
206 
207  map_type::const_iterator p = buffer.lower_bound(timeslice->getFrameIndex() - numberOfTimeslices);
208  map_type::const_iterator q = buffer.upper_bound(timeslice->getFrameIndex() + numberOfTimeslices);
209 
210  int number_of_events = 0;
211 
212  for (map_type::const_iterator i = p; i != q; ++i) {
213  number_of_events += i->second.size();
214  }
215 
216  if (number_of_events == 0) {
217  continue;
218  }
219 
220  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
221 
222  data.clear();
223 
224  buildL1(*frame, router.getModule(frame->getModuleID()), back_inserter(data));
225 
226  TH1D* h1 = manager[frame->getModuleID()];
227 
228  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
229 
230  const double t1 = *hit + frame->getFrameIndex() * getFrameTime();
231 
232  for (map_type::const_iterator i = p; i != q; ++i) {
233  for (map_type::mapped_type::const_iterator j = i->second.begin(); j != i->second.end(); ++j) {
234 
235  const double t0 = *j;
236 
237  h1->Fill(t1 - t0);
238  }
239  }
240  }
241  }
242  }
243  STATUS(endl);
244 
245 
246  // fit function
247 
248  TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2])) + [3]");
249 
250  string option = "L";
251 
252  if (debug < debug_t) {
253  option += "Q";
254  }
255 
256 
257  map<int, JStatus_t> status;
258 
259  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
260 
261  if (module->getFloor() != 0) {
262 
263  status[module->getID()] = DEFAULT;
264 
265  JManager_t::iterator p = manager.find(module->getID());
266 
267  if (p == manager.end() || p->second->GetEntries() == 0) {
268 
269  status[module->getID()] = NODATA;
270 
271  continue;
272  }
273 
274  TH1D* h1 = p->second;
275 
276  // start values
277 
278  Double_t ymax = 0.0;
279  Double_t x0 = 0.0; // [ns]
280  Double_t sigma = 250.0; // [ns]
281  Double_t ymin = 0.0;
282 
283  for (int i = 1; i != h1->GetNbinsX(); ++i) {
284 
285  const Double_t x = h1->GetBinCenter (i);
286  const Double_t y = h1->GetBinContent(i);
287 
288  if (y > ymax) {
289  ymax = y;
290  x0 = x;
291  }
292 
293  ymin += y;
294  }
295 
296  ymin /= h1->GetNbinsX();
297 
298  f1.SetParameter(0, ymax);
299  f1.SetParameter(1, x0);
300  if (binWidth_ns < sigma)
301  f1.SetParameter(2, sigma);
302  else
303  f1.FixParameter(2, binWidth_ns/sqrt(12.0));
304  f1.SetParameter(3, ymin);
305 
306  // fit
307 
308  h1->Fit(&f1, option.c_str(), "same", x0 - 5 * sigma, x0 + 5 * sigma);
309 
310 
311  if (fabs(f1.GetParameter(1)) <= 0.5*getFrameTime() && // position
312  f1.GetParameter(0) >= f1.GetParameter(3)) { // S/N
313 
314  status[module->getID()] = IN_SYNC;
315  }
316 
317 
318  // look for peaks at time offsets equal to a multiple of frame time
319 
320  map<int, JWeight> bg; // background
321  map<int, JWeight> sn; // signal(s)
322 
323  for (Int_t i = 1, ns = -(numberOfTimeslices + 1); i <= h1->GetNbinsX(); ++i) {
324 
325  const Double_t x = h1->GetBinCenter (i);
326  const Double_t y = h1->GetBinContent(i);
327 
328  while (x > (ns + 1) * getFrameTime() - TMax_ns) { ++ns; }
329 
330  if (fabs(x - ns * getFrameTime()) < 1.0 * TMax_ns)
331  sn[ns].put(y);
332  else if (fabs(x - ns * getFrameTime()) < 10.0 * TMax_ns)
333  bg[ns].put(y);
334  }
335 
336  DEBUG("Module " << setw(8) << module->getID() << endl);
337 
338  for (map<int, JWeight>::iterator i = sn.begin(); i != sn.end(); ) {
339 
340  const Double_t y = bg[i->first].getTotal() * i->second.getCount() / bg[i->first].getCount();
341  const Double_t P = TMath::PoissonI(i->second.getTotal(), y);
342 
343  DEBUG("Peak at " << setw(4) << i->first << " [frame time] "
344  << "S/N = "
345  << FIXED(7,1) << i->second.getTotal() << " / "
346  << FIXED(7,1) << y << ' '
347  << SCIENTIFIC(12,3) << P << endl);
348 
349  if (i->second.getTotal() > y && P < Pmin)
350  ++i; // keep
351  else
352  sn.erase(i++); // remove
353  }
354 
355  if (!(sn.size() == 1 &&
356  sn.begin()->first == 0)) {
357 
358  ERROR("Module " << setw(8) << module->getID() << " number of peaks " << sn.size() << endl);
359 
360  for (map<int, JWeight>::const_iterator i = sn.begin(); i != sn.end(); ++i) {
361  ERROR("Peak at " << setw(4) << i->first << " [frame time] "
362  << "S/N = "
363  << FIXED(7,1) << i->second.getTotal() << " / "
364  << FIXED(7,1) << bg[i->first].getTotal() * i->second.getCount() / bg[i->first].getCount() << endl);
365  }
366 
367  status[module->getID()] = (sn.size() == 1 ? OUT_SYNC : ERROR);
368  }
369  }
370  }
371 
372  if (overwriteDetector) {
373 
374  detector.comment.add(JMeta(argc, argv));
375 
376  for (map<int, JStatus_t>::const_iterator i = status.begin(); i != status.end(); ++i) {
377 
378  if (i->second != IN_SYNC) {
379 
380  NOTICE("Module " << setw(8) << i->first << " set out-of-sync." << endl);
381 
382  JModule& module = detector.getModule(router.getAddress(i->first));
383 
384  module.getStatus().set(MODULE_OUT_OF_SYNC);
385 
386  for (JModule::iterator pmt = module.begin(); pmt != module.end(); ++pmt) {
387  pmt->getStatus().set(OUT_OF_SYNC);
388  }
389  }
390  }
391 
392  store(detectorFile.c_str(), detector);
393  }
394 
395  if (outputFile != "") {
396  manager.Write(outputFile.c_str());
397  }
398 
399  int nin = 0;
400  int nout = 0;
401 
402  for (map<int, JStatus_t>::const_iterator i = status.begin(); i != status.end(); ++i) {
403  if (i->second == IN_SYNC) {
404  ++nin;
405  }
406  if (i->second != IN_SYNC &&
407  i->second != NODATA) {
408  ++nout;
409  }
410  }
411 
412  NOTICE("Number of modules out/in sync " << nout << '/' << nin << endl);
413 
414  QAQC(nin << ' ' << nout << endl);
415 
416  return 0;
417 }
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:46
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
#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:196
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.
#define NOTICE(A)
Definition: JMessage.hh:64
#define ERROR(A)
Definition: JMessage.hh:66
Data time slice.
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:666
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:41
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