Jpp
Enumerations | Functions
JTurbot.cc File Reference
#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 "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JDAQHitRouter.hh"
#include "JDetector/JPMTParametersMap.hh"
#include "JTrigger/JBuildL0.hh"
#include "JTrigger/JBuildL1.hh"
#include "JLang/JSinglePointer.hh"
#include "JGizmo/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 timeslice data.

Author
mdejong

Definition in file JTurbot.cc.

Enumeration Type Documentation

◆ JStatus_t

enum JStatus_t

Definition at line 45 of file JTurbot.cc.

45  {
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  };

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 60 of file JTurbot.cc.

61 {
62  using namespace std;
63  using namespace JPP;
64  using namespace KM3NETDAQ;
65 
66  JSingleFileScanner<> inputFile;
67  JLimit_t& numberOfEvents = inputFile.getLimit();
68  string outputFile;
69  string detectorFile;
70  double TMaxLocal_ns;
71  string pmtFile;
72  int numberOfTimeslices;
73  double binWidth_ns;
74  double Pmin;
75  JROOTClassSelector selector;
76  int qaqc;
77  int debug;
78 
79  try {
80 
81  JParser<> zap("Example program to search for correlations between triggered events and timeslice data.");
82 
83  zap['f'] = make_field(inputFile);
84  zap['o'] = make_field(outputFile) = "";
85  zap['a'] = make_field(detectorFile);
86  zap['n'] = make_field(numberOfEvents) = JLimit::max();
87  zap['T'] = make_field(TMaxLocal_ns) = 20.0;
88  zap['P'] = make_field(pmtFile) = "";
89  zap['N'] = make_field(numberOfTimeslices) = 100;
90  zap['W'] = make_field(binWidth_ns) = 10e3;
91  zap['p'] = make_field(Pmin) = 1.0e-6;
92  zap['C'] = make_field(selector) = "", getROOTClassSelection<JDAQTimesliceTypes_t>();
93  zap['Q'] = make_field(qaqc) = 0;
94  zap['d'] = make_field(debug) = 2;
95 
96  zap(argc, argv);
97  }
98  catch(const exception& error) {
99  FATAL(error.what() << endl);
100  }
101 
102 
103  cout.tie(&cerr);
104 
106 
107  try {
108  load(detectorFile, detector);
109  }
110  catch(const JException& error) {
111  FATAL(error);
112  }
113 
114  if (detector.empty()) {
115  FATAL("Empty detector " << detectorFile << endl);
116  }
117 
118  const JDAQHitRouter router(detector);
119 
120  const double TMax_ns = max(binWidth_ns, getMaximalTime(detector));
121 
122  typedef double hit_type;
123 
124  JBuildL0<hit_type> buildL0;
125  JBuildL1<hit_type> buildL1(JBuildL1Parameters(TMaxLocal_ns, true));
126  vector <hit_type> data;
127 
128  typedef JManager<int, TH1D> JManager_t;
129 
130  const Double_t xmin = -(numberOfTimeslices + 1) * getFrameTime() - 0.5*binWidth_ns;
131  const Double_t xmax = +(numberOfTimeslices + 1) * getFrameTime() + 0.5*binWidth_ns;
132  const Int_t nx = (Int_t) ((xmax - xmin) / binWidth_ns);
133 
134  JManager_t manager(new TH1D("M_%", NULL, nx, xmin, xmax));
135 
137 
139 
141 
142  if (selector == "") {
143 
144  if ((ps = new JTreeScanner< JAssertConversion<JDAQTimesliceSN, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
145  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL2, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
146  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL1, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
147  (ps = new JTreeScanner< JAssertConversion<JDAQTimeslice, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
148  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL0, JDAQTimeslice> >(inputFile))->getEntries() != 0) {
149  } else {
150  FATAL("No timeslice data." << endl);
151  }
152 
153  NOTICE("Selected class " << ps->getClassName() << endl);
154 
155  } else {
156 
157  ps = zmap[selector];
158 
159  ps->configure(inputFile);
160  }
161 
162  ps->setLimit(inputFile.getLimit());
163 
164  typedef map<int, vector<double> > map_type; // frame index -> event times
165 
166  map_type buffer;
167 
168  for (JTreeScanner<JDAQEvent> in(inputFile.getFilename()); in.hasNext(); ) {
169 
170  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
171 
172  JDAQEvent* event = in.next();
173 
174  // Average time of triggered hits
175 
176  double t0 = 0.0;
177 
178  for (JDAQEvent::const_iterator<JDAQTriggeredHit> hit = event->begin<JDAQTriggeredHit>(); hit != event->end<JDAQTriggeredHit>(); ++hit) {
179  t0 += getTime(*hit, router.getPMT(*hit));
180  }
181 
182  t0 /= event->size<JDAQTriggeredHit>();
183  t0 += event->getFrameIndex() * getFrameTime();
184 
185  buffer[event->getFrameIndex()].push_back(t0);
186  }
187  STATUS(endl);
188 
189 
190  while (ps->hasNext()) {
191 
192  STATUS("event: " << setw(10) << ps->getCounter() << '\r'); DEBUG(endl);
193 
194  JDAQTimeslice* timeslice = ps->next();
195 
196  map_type::const_iterator p = buffer.lower_bound(timeslice->getFrameIndex() - numberOfTimeslices);
197  map_type::const_iterator q = buffer.upper_bound(timeslice->getFrameIndex() + numberOfTimeslices);
198 
199  int number_of_events = 0;
200 
201  for (map_type::const_iterator i = p; i != q; ++i) {
202  number_of_events += i->second.size();
203  }
204 
205  if (number_of_events == 0) {
206  continue;
207  }
208 
209  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
210 
211  data.clear();
212 
213  buildL1(*frame, router.getModule(frame->getModuleID()), back_inserter(data));
214 
215  TH1D* h1 = manager[frame->getModuleID()];
216 
217  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
218 
219  const double t1 = *hit + frame->getFrameIndex() * getFrameTime();
220 
221  for (map_type::const_iterator i = p; i != q; ++i) {
222  for (map_type::mapped_type::const_iterator j = i->second.begin(); j != i->second.end(); ++j) {
223 
224  const double t0 = *j;
225 
226  h1->Fill(t1 - t0);
227  }
228  }
229  }
230  }
231  }
232  STATUS(endl);
233 
234 
235  // fit function
236 
237  TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2])) + [3]");
238 
239  string option = "L";
240 
241  if (debug < debug_t) {
242  option += "Q";
243  }
244 
245 
246  map<int, JStatus_t> status;
247 
248  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
249 
250  status[module->getID()] = DEFAULT;
251 
252  JManager_t::iterator p = manager.find(module->getID());
253 
254  if (p == manager.end() || p->second->GetEntries() == 0) {
255 
256  status[module->getID()] = NODATA;
257 
258  continue;
259  }
260 
261  TH1D* h1 = p->second;
262 
263  // start values
264 
265  Double_t ymax = 0.0;
266  Double_t x0 = 0.0; // [ns]
267  Double_t sigma = 250.0; // [ns]
268  Double_t ymin = 0.0;
269 
270  for (int i = 1; i != h1->GetNbinsX(); ++i) {
271 
272  const Double_t x = h1->GetBinCenter (i);
273  const Double_t y = h1->GetBinContent(i);
274 
275  if (y > ymax) {
276  ymax = y;
277  x0 = x;
278  }
279 
280  ymin += y;
281  }
282 
283  ymin /= h1->GetNbinsX();
284 
285  f1.SetParameter(0, ymax);
286  f1.SetParameter(1, x0);
287  if (binWidth_ns < sigma)
288  f1.SetParameter(2, sigma);
289  else
290  f1.FixParameter(2, binWidth_ns/sqrt(12.0));
291  f1.SetParameter(3, ymin);
292 
293  // fit
294 
295  h1->Fit(&f1, option.c_str(), "same", x0 - 5 * sigma, x0 + 5 * sigma);
296 
297 
298  if (fabs(f1.GetParameter(1)) <= 0.5*getFrameTime() && // position
299  f1.GetParameter(0) >= f1.GetParameter(3)) { // S/N
300 
301  status[module->getID()] = IN_SYNC;
302  }
303 
304 
305  // look for peaks at time offsets equal to a multiple of frame time
306 
307  JWeight bg; // background
308  map<int, JWeight> sn; // signal(s)
309 
310  for (Int_t i = 1, ns = -(numberOfTimeslices + 1); i <= h1->GetNbinsX(); ++i) {
311 
312  const Double_t x = h1->GetBinCenter (i);
313  const Double_t y = h1->GetBinContent(i);
314 
315  while (x > (ns + 1) * getFrameTime() - TMax_ns) { ++ns; }
316 
317  if (fabs(x - ns * getFrameTime()) < TMax_ns)
318  sn[ns].put(y);
319  else
320  bg .put(y);
321  }
322 
323  DEBUG("Module " << setw(8) << module->getID() << endl);
324 
325  for (map<int, JWeight>::iterator i = sn.begin(); i != sn.end(); ) {
326 
327  const Double_t y = bg.getTotal() * i->second.getCount() / bg.getCount();
328  const Double_t P = TMath::PoissonI(i->second.getTotal(), y);
329 
330  DEBUG("Peak at " << setw(4) << i->first << " [frame time] "
331  << "S/N = "
332  << FIXED(7,1) << i->second.getTotal() << " / "
333  << FIXED(7,1) << y << ' '
334  << SCIENTIFIC(12,3) << P << endl);
335 
336  if (i->second.getTotal() > y && P < Pmin)
337  ++i; // keep
338  else
339  sn.erase(i++); // remove
340  }
341 
342  if (!(sn.size() == 1 &&
343  sn.begin()->first == 0)) {
344 
345  WARNING("Module " << setw(8) << module->getID() << " number of peaks " << sn.size() << endl);
346 
347  for (map<int, JWeight>::const_iterator i = sn.begin(); i != sn.end(); ++i) {
348  WARNING("Peak at " << setw(4) << i->first << " [frame time] "
349  << "S/N = "
350  << FIXED(7,1) << i->second.getTotal() << " / "
351  << FIXED(7,1) << bg.getTotal() * i->second.getCount() / bg.getCount() << endl);
352  }
353 
354  status[module->getID()] = (sn.size() == 1 ? OUT_OF_SYNC : ERROR);
355  }
356  }
357 
358 
359  if (pmtFile != "") {
360 
361  JPMTParametersMap parameters;
362 
363  try {
364  parameters.load(pmtFile.c_str());
365  }
366  catch(const JException& error) {}
367 
368  for (map<int, JStatus_t>::const_iterator i = status.begin(); i != status.end(); ++i) {
369 
370  if (i->second != IN_SYNC) {
371 
372  NOTICE("Module " << setw(8) << i->first << " set QEs of all PMTs to zero." << endl);
373 
374  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
375  parameters[JPMTIdentifier(i->first, pmt)].QE = 0.0;
376  }
377  }
378  }
379 
380  ofstream out(pmtFile.c_str());
381 
382  parameters.comment.add(JMeta(argc, argv));
383 
384  out << parameters << endl;
385 
386  out.close();
387  }
388 
389  if (outputFile != "") {
390  manager.Write(outputFile.c_str());
391  }
392 
393  int nin = 0;
394  int nout = 0;
395 
396  for (map<int, JStatus_t>::const_iterator i = status.begin(); i != status.end(); ++i) {
397  if (i->second == IN_SYNC) {
398  ++nin;
399  }
400  if (i->second != IN_SYNC &&
401  i->second != NODATA) {
402  ++nout;
403  }
404  }
405 
406  NOTICE("Number of modules out/in sync " << nout << '/' << nin << endl);
407 
408  QAQC(nin << ' ' << nout << endl);
409 
410  return 0;
411 }
JTOOLS::JWeight::getCount
long long int getCount() const
Get total count.
Definition: JWeight.hh:68
JSUPPORT::JAutoTreeScanner
Auxiliary class to select JTreeScanner based on ROOT class name.
Definition: JAutoTreeScanner.hh:34
KM3NETDAQ::JDAQEvent
DAQ Event.
Definition: JDAQEvent.hh:30
JLANG::JSinglePointer
The template JSinglePointer class can be used to hold a pointer to an object.
Definition: JSinglePointer.hh:24
JSUPPORT::JLimit
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
JLANG::JObjectStreamIO::load
void load(const char *file_name)
Load from input file.
Definition: JObjectStreamIO.hh:30
FIXED
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:481
JLANG::JType
Auxiliary class for a type holder.
Definition: JType.hh:19
JTRIGGER::JBuildL1
Template L1 hit builder.
Definition: JBuildL1.hh:85
JDETECTOR::load
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
Definition: JDetectorToolkit.hh:476
JTOOLS::JWeight
Weight calculator.
Definition: JWeight.hh:23
std::vector
Definition: JSTDTypes.hh:12
KM3NETDAQ::JDAQChronometer::getFrameIndex
int getFrameIndex() const
Get frame index.
Definition: JDAQChronometer.hh:132
JTOOLS::j
int j
Definition: JPolint.hh:634
KM3NETDAQ::NUMBER_OF_PMTS
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
KM3NETDAQ::JDAQTimeslice
Data time slice.
Definition: JDAQTimeslice.hh:30
JPARSER::JParser
Utility class to parse command line options.
Definition: JParser.hh:1493
JGIZMO::JManager
Auxiliary class to manage set of compatible ROOT objects (e.g.
Definition: JManager.hh:40
NOTICE
#define NOTICE(A)
Definition: JMessage.hh:64
QAQC
#define QAQC(A)
QA/QC output macro.
Definition: JMessage.hh:100
KM3NETDAQ::getFrameTime
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
JPP
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Definition: JAAnetToolkit.hh:37
KM3NETDAQ::JDAQTriggeredHit
DAQ triggered hit.
Definition: JDAQTriggeredHit.hh:20
ERROR
#define ERROR(A)
Definition: JMessage.hh:66
JTRIGGER::JBuildL0
Template L0 hit builder.
Definition: JBuildL0.hh:34
JDETECTOR::getMaximalTime
double getMaximalTime(const JDetector &detector)
Get maximal time between modules in detector following causality.
Definition: JDetectorToolkit.hh:306
JTRIGGER::JBuildL1Parameters
Auxiliary data structure for L1 build parameters.
Definition: JBuildL1.hh:37
WARNING
#define WARNING(A)
Definition: JMessage.hh:65
JROOT::JROOTClassSelector
Auxiliary class to select ROOT class based on class name.
Definition: JROOTClassSelector.hh:32
debug
int debug
debug level
Definition: JSirene.cc:59
JSUPPORT::JTreeScanner
Template definition for direct access of elements in ROOT TChain.
Definition: JTreeScanner.hh:91
KM3NETDAQ::JDAQEvent::const_iterator
Template const_iterator.
Definition: JDAQEvent.hh:68
KM3NETDAQ::JDAQTriggerCounter::next
JTriggerCounter_t next()
Increment trigger counter.
Definition: JDAQTriggerCounter.hh:121
JLANG::JAssertConversion
Auxialiary class to assert type conversion.
Definition: JConversion.hh:66
JDETECTOR::JPMTParametersMap
Auxiliary class for map of PMT parameters.
Definition: JPMTParametersMap.hh:88
JAANET::getTime
double getTime(const Hit &hit)
Get true time of hit.
Definition: JAAnetToolkit.hh:87
SCIENTIFIC
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:518
std::map
Definition: JSTDTypes.hh:16
STATUS
#define STATUS(A)
Definition: JMessage.hh:63
JTOOLS::JWeight::getTotal
double getTotal() const
Get total weight.
Definition: JWeight.hh:79
make_field
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
JDETECTOR::JPMTIdentifier
PMT identifier.
Definition: JPMTIdentifier.hh:30
JDETECTOR::JDetector
Detector data structure.
Definition: JDetector.hh:80
JDETECTOR::JPMTParametersMap::comment
JComment comment
Definition: JPMTParametersMap.hh:406
JAANET::detector
Detector file.
Definition: JHead.hh:130
JSUPPORT::JMeta
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:71
DEBUG
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
JDETECTOR::JPMTParametersMap::QE
double QE
Global QE.
Definition: JPMTParametersMap.hh:417
std
Definition: jaanetDictionary.h:36
KM3NETDAQ
KM3NeT DAQ data structures and auxiliaries.
Definition: DataQueue.cc:39
JTOOLS::JWeight::put
void put(const double w)
Put weight.
Definition: JWeight.hh:55
qaqc
int qaqc
QA/QC file descriptor.
FATAL
#define FATAL(A)
Definition: JMessage.hh:67
outputFile
string outputFile
Definition: JDAQTimesliceSelector.cc:37
JEEP::JComment::add
JComment & add(const std::string &comment)
Add comment.
Definition: JComment.hh:55
JEEP::debug_t
debug
Definition: JMessage.hh:29
JLANG::JException
General exception.
Definition: JException.hh:23
JSUPPORT::JSingleFileScanner<>
JDETECTOR::JDAQHitRouter
Simple wrapper around JModuleRouter class for direct addressing of PMT data in detector data structur...
Definition: JDAQHitRouter.hh:24