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