Jpp  16.0.0-rc.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Enumerations | Functions
JTurbot.cc File Reference

Example program to search for correlations between triggered events and time slice data. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <limits>
#include <map>
#include <vector>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TF1.h"
#include "TMath.h"
#include "km3net-dataformat/online/JDAQ.hh"
#include "km3net-dataformat/definitions/pmt_status.hh"
#include "km3net-dataformat/definitions/module_status.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JDAQHitRouter.hh"
#include "JTrigger/JBuildL0.hh"
#include "JTrigger/JBuildL1.hh"
#include "JLang/JSinglePointer.hh"
#include "JROOT/JManager.hh"
#include "JTools/JWeight.hh"
#include "JROOT/JROOTClassSelector.hh"
#include "JSupport/JSingleFileScanner.hh"
#include "JSupport/JTreeScannerInterface.hh"
#include "JSupport/JTreeScanner.hh"
#include "JSupport/JAutoTreeScanner.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JMeta.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Enumerations

enum  JStatus_t
 

Functions

int main (int argc, char **argv)
 

Detailed Description

Example program to search for correlations between triggered events and time slice data.


This application can be used to identify optical modules that are out-of-sync with respect to the central White Rabbit clock.
To this end, time correlations are searched for between triggered events and local coincidences, referred to as L1 hits. The local coincidences are independently obtained from time slice data and can thereby cross the time slice border of a triggered event.
The MODULE_OUT_OF_SYNC of status the modules that are found to be out-of-sync are set.

The number of time slices to look for L1 hits away from the triggered event can be set via command line option -N <number of time slices>.
The larger this value, the bigger time offsets can be covered. To identify the modules that are out-of-sync as fast as possible (without determining the time offset), 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.

Author
mdejong

Definition in file JTurbot.cc.

Enumeration Type Documentation

enum JStatus_t

Definition at line 46 of file JTurbot.cc.

46  {
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  };
#define ERROR(A)
Definition: JMessage.hh:66

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 72 of file JTurbot.cc.

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