Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JTurbot.cc File Reference

Example program to search for correlations between triggered events and timeslice 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 "JDAQ/JDAQ.hh"
#include "JDAQ/JDAQEvent.hh"
#include "JDAQ/JDAQTimeslice.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.

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.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 48 of file JTurbot.cc.

49 {
50  using namespace std;
51  using namespace JPP;
52  using namespace KM3NETDAQ;
53 
54  JSingleFileScanner<> inputFile;
55  JLimit_t& numberOfEvents = inputFile.getLimit();
56  string outputFile;
57  string detectorFile;
58  double TMaxLocal_ns;
59  string pmtFile;
60  int numberOfTimeslices;
61  double binWidth_ns;
62  double Pmin;
63  JROOTClassSelector selector;
64  int debug;
65 
66  try {
67 
68  JParser<> zap("Example program to search for correlations between triggered events and timeslice data.");
69 
70  zap['f'] = make_field(inputFile);
71  zap['o'] = make_field(outputFile) = "";
72  zap['a'] = make_field(detectorFile);
73  zap['n'] = make_field(numberOfEvents) = JLimit::max();
74  zap['T'] = make_field(TMaxLocal_ns) = 20.0;
75  zap['P'] = make_field(pmtFile) = "";
76  zap['N'] = make_field(numberOfTimeslices) = 100;
77  zap['W'] = make_field(binWidth_ns) = 10e3;
78  zap['p'] = make_field(Pmin) = 1.0e-6;
79  zap['C'] = make_field(selector) = "", getROOTClassSelection<JDAQTimesliceTypes_t>();
80  zap['d'] = make_field(debug) = 2;
81 
82  zap(argc, argv);
83  }
84  catch(const exception& error) {
85  FATAL(error.what() << endl);
86  }
87 
88 
89  cout.tie(&cerr);
90 
91  JDetector detector;
92 
93  try {
94  load(detectorFile, detector);
95  }
96  catch(const JException& error) {
97  FATAL(error);
98  }
99 
100  const JDAQHitRouter router(detector);
101 
102  const double TMax_ns = max(binWidth_ns, getMaximalTime(detector));
103 
104  typedef double hit_type;
105 
106  JBuildL0<hit_type> buildL0;
107  JBuildL1<hit_type> buildL1(TMaxLocal_ns, true);
108  vector <hit_type> data;
109 
110  typedef JManager<int, TH1D> JManager_t;
111 
112  const Double_t xmin = -(numberOfTimeslices + 1) * getFrameTime() - 0.5*binWidth_ns;
113  const Double_t xmax = +(numberOfTimeslices + 1) * getFrameTime() + 0.5*binWidth_ns;
114  const Int_t nx = (Int_t) ((xmax - xmin) / binWidth_ns);
115 
116  JManager_t manager(new TH1D("M_%", NULL, nx, xmin, xmax));
117 
119 
120  JSinglePointer< JTreeScannerInterface<JDAQTimeslice> > ps;
121 
122  JAutoTreeScanner<JDAQTimeslice> zmap = JType<JDAQTimesliceTypes_t>();
123 
124  if (selector == "") {
125 
126  if ((ps = new JTreeScanner< JAssertConversion<JDAQTimesliceSN, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
127  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL2, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
128  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL1, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
129  (ps = new JTreeScanner< JAssertConversion<JDAQTimeslice, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
130  (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL0, JDAQTimeslice> >(inputFile))->getEntries() != 0) {
131  } else {
132  FATAL("No timeslice data." << endl);
133  }
134 
135  NOTICE("Selected class " << ps->getClassName() << endl);
136 
137  } else {
138 
139  ps = zmap[selector];
140 
141  ps->configure(inputFile);
142  }
143 
144 
145  typedef map<int, vector<double> > map_type; // frame index -> event times
146 
147  map_type buffer;
148 
149  for (JTreeScanner<JDAQEvent> in(inputFile.getFilename(), inputFile.getLimit()); in.hasNext(); ) {
150 
151  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
152 
153  JDAQEvent* event = in.next();
154 
155  // Average time of triggered hits
156 
157  double t0 = 0.0;
158 
159  for (JDAQEvent::const_iterator<JDAQTriggeredHit> hit = event->begin<JDAQTriggeredHit>(); hit != event->end<JDAQTriggeredHit>(); ++hit) {
160  t0 += getTime(*hit, router.getPMT(*hit));
161  }
162 
163  t0 /= event->size<JDAQTriggeredHit>();
164  t0 += event->getFrameIndex() * getFrameTime();
165 
166  buffer[event->getFrameIndex()].push_back(t0);
167  }
168  STATUS(endl);
169 
170 
171  while (ps->hasNext()) {
172 
173  STATUS("event: " << setw(10) << ps->getCounter() << '\r'); DEBUG(endl);
174 
175  JDAQTimeslice* timeslice = ps->next();
176 
177  map_type::const_iterator p = buffer.lower_bound(timeslice->getFrameIndex() - numberOfTimeslices);
178  map_type::const_iterator q = buffer.upper_bound(timeslice->getFrameIndex() + numberOfTimeslices);
179 
180  int number_of_events = 0;
181 
182  for (map_type::const_iterator i = p; i != q; ++i) {
183  number_of_events += i->second.size();
184  }
185 
186  if (number_of_events == 0) {
187  continue;
188  }
189 
190  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
191 
192  data.clear();
193 
194  buildL1(*frame, router.getModule(frame->getModuleID()), back_inserter(data));
195 
196  TH1D* h1 = manager[frame->getModuleID()];
197 
198  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
199 
200  const double t1 = *hit + frame->getFrameIndex() * getFrameTime();
201 
202  for (map_type::const_iterator i = p; i != q; ++i) {
203  for (map_type::mapped_type::const_iterator j = i->second.begin(); j != i->second.end(); ++j) {
204 
205  const double t0 = *j;
206 
207  h1->Fill(t1 - t0);
208  }
209  }
210  }
211  }
212  }
213  STATUS(endl);
214 
215  if (pmtFile != "") {
216 
217  JPMTParametersMap parameters;
218 
219  try {
220  parameters.load(pmtFile.c_str());
221  }
222  catch(const JException& error) {}
223 
224  TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2])) + [3]");
225 
226  string option = "L";
227 
228  if (debug < debug_t) {
229  option += "Q";
230  }
231 
232  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
233 
234  bool status = false;
235 
236  JManager_t::iterator p = manager.find(module->getID());
237 
238  if (p != manager.end()) {
239 
240  TH1D* h1 = p->second;
241 
242  if (h1->GetEntries() != 0) {
243 
244  // start values
245 
246  Double_t ymax = 0.0;
247  Double_t x0 = 0.0; // [ns]
248  Double_t sigma = 250.0; // [ns]
249  Double_t ymin = 0.0;
250 
251  for (int i = 1; i != h1->GetNbinsX(); ++i) {
252 
253  const Double_t x = h1->GetBinCenter (i);
254  const Double_t y = h1->GetBinContent(i);
255 
256  if (y > ymax) {
257  ymax = y;
258  x0 = x;
259  }
260 
261  ymin += y;
262  }
263 
264  ymin /= h1->GetNbinsX();
265 
266  f1.SetParameter(0, ymax);
267  f1.SetParameter(1, x0);
268  if (binWidth_ns < sigma)
269  f1.SetParameter(2, sigma);
270  else
271  f1.FixParameter(2, binWidth_ns/sqrt(12.0));
272  f1.SetParameter(3, ymin);
273 
274  // fit
275 
276  h1->Fit(&f1, option.c_str(), "same", x0 - 5 * sigma, x0 + 5 * sigma);
277 
278  status = (fabs(f1.GetParameter(1)) <= 0.5*getFrameTime() &&
279  f1.GetParameter(0) >= f1.GetParameter(3));
280  }
281 
282  if (h1->GetEntries() != 0) {
283 
284  // look for peaks at offsets equal to a multiple of frame time
285 
286  JWeight bg; // background
287  map<int, JWeight> sn; // signal(s)
288 
289  for (Int_t i = 1, ns = -(numberOfTimeslices + 1); i <= h1->GetNbinsX(); ++i) {
290 
291  const Double_t x = h1->GetBinCenter (i);
292  const Double_t y = h1->GetBinContent(i);
293 
294  while (x > (ns + 1) * getFrameTime() - TMax_ns) { ++ns; }
295 
296  if (fabs(x - ns * getFrameTime()) < TMax_ns)
297  sn[ns].put(y);
298  else
299  bg .put(y);
300  }
301 
302 
303  for (map<int, JWeight>::iterator i = sn.begin(); i != sn.end(); ) {
304 
305  const Double_t y = bg.getTotal() * i->second.getCount() / bg.getCount();
306  const Double_t P = TMath::PoissonI(i->second.getTotal(), y);
307 
308  if (i->second.getTotal() > y && P < Pmin)
309  ++i; // keep
310  else
311  sn.erase(i++); // remove
312  }
313 
314  if (sn.size() != 1) {
315 
316  ERROR("Module " << setw(8) << module->getID() << " number of peaks " << sn.size() << " != 1" << endl);
317 
318  for (map<int, JWeight>::const_iterator i = sn.begin(); i != sn.end(); ++i) {
319  WARNING("Peak at " << setw(4) << i->first << " [frame time] S/N = "
320  << FIXED(7,1) << i->second.getTotal() << " / "
321  << FIXED(7,1) << bg.getTotal() * i->second.getCount() / bg.getCount() << endl);
322  }
323  }
324 
325  status &= (sn.size() == 1 &&
326  sn.begin()->first == 0);
327  }
328  }
329 
330  if (!status) {
331 
332  NOTICE("Module " << module->getID() << " set QEs of all PMTs to zero." << endl);
333 
334  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
335  parameters[JPMTIdentifier(module->getID(), pmt)].QE = 0.0;
336  }
337  }
338  }
339 
340  ofstream out(pmtFile.c_str());
341 
342  parameters.comment.add(JMeta(argc, argv));
343 
344  out << parameters << endl;
345 
346  out.close();
347  }
348 
349  if (outputFile != "") {
350  manager.Write(outputFile.c_str());
351  }
352 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:71
Utility class to parse command line options.
Definition: JParser.hh:1410
void clear()
Clear data.
#define WARNING(A)
Definition: JMessage.hh:63
debug
Definition: JMessage.hh:27
Auxiliary class to manage set of histograms.
#define STATUS(A)
Definition: JMessage.hh:61
Template const_iterator.
Definition: JDAQEvent.hh:69
Structure to store the ToT mean and standard deviation of the hits produced by a nanobeacon in a sour...
Auxiliary class for a type holder.
Definition: JType.hh:19
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:461
double getTime(const Hit &hit)
Get true time of hit.
string outputFile
int getFrameIndex() const
Get frame index.
JLimit JLimit_t
Type definition of limit.
Definition: JLimit.hh:214
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1836
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
#define NOTICE(A)
Definition: JMessage.hh:62
#define ERROR(A)
Definition: JMessage.hh:64
Data time slice.
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
int debug
debug level
Definition: JSirene.cc:59
double getMaximalTime(const JDetector &detector)
Get maximal time between modules in detector following causality.
#define FATAL(A)
Definition: JMessage.hh:65
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:72
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
JTriggerCounter_t next()
Increment trigger counter.
double hit_type
Definition: JDataFilter.cc:89
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:60