CMS 3D CMS Logo

Scanner.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_FWLite_Scanner_h
2 #define PhysicsTools_FWLite_Scanner_h
3 // these includes are FWLite-safe
6 // these are from ROOT, so they're safe too
7 #include <TString.h>
8 #include <TRegexp.h>
9 #include <TObjString.h>
10 #include <TObjArray.h>
11 #include <TDirectory.h>
12 #include <TEnv.h>
13 #include <TClass.h>
14 
15 #if !defined(__CINT__) && !defined(__MAKECINT__)
16 #include <RooFit.h>
17 #include <RooArgList.h>
18 #include <RooDataSet.h>
19 #include <RooRealVar.h>
20 #include <RooCategory.h>
21 
23 #endif
24 
26 
27 namespace fwlite {
28 
45  template <typename Collection>
46  class Scanner {
47  public:
50 
52  Scanner(fwlite::EventBase *ev, const char *label, const char *instance = "", const char *process = "")
53  : event_(ev),
54  label_(label),
56  printFullEventId_(ev->isRealData()),
58  exprSep_(":"),
59  maxEvents_(-1),
60  maxLinesToPrint_(50) {
62  }
63 
64  //------------------------------------------------------------------------------------------------------------------------------------
80  void scan(const char *exprs, const char *cut = "", int nmax = -1) {
83 
84  TObjArray *exprArray = TString(exprs).Tokenize(exprSep_);
85  int rowline = 0;
86  if (printFullEventId_) {
87  printf(" : %9s : %4s : %9s : %3s", "RUN", "LUMI", "EVENT", "#IT");
88  rowline += 3 * 4 + 9 + 4 + 9 + 3 - 1; // -1 as first char remain blank
89  } else {
90  printf(" : %5s : %3s", "EVENT", "#IT");
91  rowline += 3 + 6 + 3 + 3 - 1; // -1 as first char remain blank
92  }
93  for (int i = 0; i < exprArray->GetEntries(); ++i) {
94  TString str = ((TObjString *)(*exprArray)[i])->GetString();
95  std::string lb = str.Data();
96  std::string ex = str.Data();
97  if ((ex[0] == '@') && (ex.find('=') != std::string::npos)) {
98  lb = lb.substr(1, ex.find('=') - 1);
99  ex = ex.substr(ex.find('=') + 1);
100  }
101  scanner.addExpression(ex.c_str());
102  printf(" : %8s",
103  (lb.size() > 8 ? lb.substr(lb.size() - 8) : lb)
104  .c_str()); // the rightmost part is usually the more interesting one
105  rowline += 3 + 8;
106  }
107  std::cout << " :" << std::endl;
108  rowline += 2;
109  delete exprArray;
110 
111  TString rule('-', rowline);
112  std::cout << " " << rule << " " << std::endl;
113 
114  if (strlen(cut))
115  scanner.setCut(cut);
116 
117  int iev = 0, line = 0;
118  for (event_->toBegin(); (iev != nmax) && !event_->atEnd(); ++iev, ++(*event_)) {
119  if (!selectEvent(*event_))
120  continue;
121  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
122  if (handle_.failedToGet()) {
123  if (ignoreExceptions_)
124  continue;
125  }
126  const Collection &vals = *handle_;
127  for (size_t j = 0, n = vals.size(); j < n; ++j) {
128  if (!scanner.test(&vals[j]))
129  continue;
130  if (printFullEventId_) {
132  printf(" : %9u : %4u : %9llu : %3lu", id.run(), id.luminosityBlock(), id.event(), (unsigned long)j);
133  } else {
134  printf(" : %5d : %3lu", iev, (unsigned long)j);
135  }
136  scanner.print(&vals[j]);
137  std::cout << " :" << std::endl;
138  if (++line == maxLinesToPrint_) {
139  line = 0;
140  if (!wantMore()) {
141  iev = nmax - 1; // this is to exit the outer loop
142  break; // and this to exit the inner one
143  }
144  }
145  }
146  }
147  std::cout << std::endl;
148  }
149 
150  //------------------------------------------------------------------------------------------------------------------------------------
154  size_t count(const char *cut) {
155  helper::ScannerBase scanner(objType);
157 
158  scanner.setCut(cut);
159 
160  size_t npass = 0;
161  int iev = 0;
162  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
163  if (maxEvents_ > -1 && iev > maxEvents_)
164  break;
165  if (!selectEvent(*event_))
166  continue;
167  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
168  const Collection &vals = *handle_;
169  for (size_t j = 0, n = vals.size(); j < n; ++j) {
170  if (scanner.test(&vals[j]))
171  npass++;
172  }
173  }
174  return npass;
175  }
176 
178  size_t countEvents() {
179  size_t npass = 0;
180  int iev = 0;
181  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
182  if (maxEvents_ > -1 && iev > maxEvents_)
183  break;
184  if (selectEvent(*event_))
185  npass++;
186  }
187  return npass;
188  }
189 
190  //------------------------------------------------------------------------------------------------------------------------------------
199  TH1 *draw(const char *expr, const char *cut, TString drawopt, TH1 *hist) {
200  // prep the machinery
201  helper::ScannerBase scanner(objType);
203  if (!scanner.addExpression(expr))
204  return nullptr;
205  if (strlen(cut))
206  scanner.setCut(cut);
207 
208  // check histo
209  if (hist == nullptr) {
210  std::cerr << "Method draw(expr, cut, drawopt, hist) cannot be called with null 'hist'. Use the other draw "
211  "methods instead."
212  << std::endl;
213  return nullptr;
214  }
215 
216  // fill histogram
217  int iev = 0;
218  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
219  if (maxEvents_ > -1 && iev > maxEvents_)
220  break;
221  if (!selectEvent(*event_))
222  continue;
223  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
224  const Collection &vals = *handle_;
225  for (size_t j = 0, n = vals.size(); j < n; ++j) {
226  scanner.fill1D(&vals[j], hist);
227  }
228  }
229 
230  if (drawopt.Contains("NORM", TString::kIgnoreCase) && (hist->Integral() != 0)) {
231  hist->Sumw2();
232  hist->Scale(1.0 / hist->Integral());
233  // remove the "NORM" because THistPainter doesn't understand it
234  drawopt(TRegexp("[Nn][Oo][Rr][Mm]")) = "";
235  }
236 
237  if (!drawopt.Contains("GOFF", TString::kIgnoreCase))
238  hist->Draw(drawopt);
239  return hist;
240  }
241 
248  TH1 *draw(const char *expr,
249  const char *cut = "",
250  TString drawopt = "",
251  const char *hname = "htemp",
252  const TH1 *htemplate = nullptr) {
253  TH1 *hist = nullptr;
254  if (htemplate != nullptr) {
255  if ((strcmp(hname, "htemp") == 0) && (strcmp(hname, htemplate->GetName()) != 0))
256  htempDelete();
257  hist = (TH1 *)htemplate->Clone(hname);
258  } else if (drawopt.Contains("SAME", TString::kIgnoreCase)) {
259  hist = getSameH1(hname);
260  }
261 
262  // if in the end we found no way to make "hist"
263  if (hist == nullptr) {
264  if (strcmp(hname, "htemp") == 0)
265  htempDelete();
266  hist = new TH1F(hname, "", gEnv->GetValue("Hist.Binning.1D.x", 100), 0, 0);
267  hist->SetCanExtend(TH1::kAllAxes);
268  }
269  hist->SetTitle((strlen(cut) ? TString(expr) + "{" + cut + "}" : TString(expr)));
270  hist->GetXaxis()->SetTitle(expr);
271  return draw(expr, cut, drawopt, hist);
272  }
273 
276  TH1 *draw(const char *expr,
277  int nbins,
278  double xlow,
279  double xhigh,
280  const char *cut = "",
281  const char *drawopt = "",
282  const char *hname = "htemp") {
283  if (TString(drawopt).Contains("SAME", TString::kIgnoreCase)) {
284  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
285  TH1 *hsame = getSameH1(hname);
286  return draw(expr, cut, drawopt, hsame);
287  }
288  if (strcmp(hname, "htemp") == 0)
289  htempDelete();
290  TH1 *htemp = new TH1F(hname, expr, nbins, xlow, xhigh);
291  if (strlen(cut))
292  htemp->SetTitle(TString(expr) + "{" + cut + "}");
293  htemp->GetXaxis()->SetTitle(expr);
294  return draw(expr, cut, drawopt, htemp);
295  }
298  TH1 *draw(const char *expr,
299  int nbins,
300  double *xbins,
301  const char *cut = "",
302  const char *drawopt = "",
303  const char *hname = "htemp") {
304  if (TString(drawopt).Contains("SAME", TString::kIgnoreCase)) {
305  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
306  TH1 *hsame = getSameH1(hname);
307  return draw(expr, cut, drawopt, hsame);
308  }
309  if (strcmp(hname, "htemp") == 0)
310  htempDelete();
311  TH1 *htemp = new TH1F(hname, expr, nbins, xbins);
312  if (strlen(cut))
313  htemp->SetTitle(TString(expr) + "{" + cut + "}");
314  htemp->GetXaxis()->SetTitle(expr);
315  return draw(expr, cut, drawopt, htemp);
316  }
317 
318  //------------------------------------------------------------------------------------------------------------------------------------
320  TProfile *drawProf(TString xexpr, TString yexpr, const char *cut, TString drawopt, TProfile *hist) {
321  // prep the machinery
322  helper::ScannerBase scanner(objType);
324  if (!scanner.addExpression(xexpr.Data()))
325  return nullptr;
326  if (!scanner.addExpression(yexpr.Data()))
327  return nullptr;
328  if (strlen(cut))
329  scanner.setCut(cut);
330 
331  // check histo
332  if (hist == nullptr) {
333  std::cerr << "Method drawProf(xexpr, yexpr, cut, drawopt, hist) cannot be called with null 'hist'. Use the "
334  "other draw methods instead."
335  << std::endl;
336  return nullptr;
337  }
338 
339  // fill histogram
340  int iev = 0;
341  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
342  if (maxEvents_ > -1 && iev > maxEvents_)
343  break;
344  if (!selectEvent(*event_))
345  continue;
346  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
347  const Collection &vals = *handle_;
348  for (size_t j = 0, n = vals.size(); j < n; ++j) {
349  scanner.fillProf(&vals[j], hist);
350  }
351  }
352 
353  if (!strlen(hist->GetTitle()))
354  hist->SetTitle((strlen(cut) ? yexpr + ":" + xexpr + "{" + cut + "}" : yexpr + ":" + xexpr));
355  if (!strlen(hist->GetXaxis()->GetTitle()))
356  hist->GetXaxis()->SetTitle(xexpr);
357  if (!strlen(hist->GetYaxis()->GetTitle()))
358  hist->GetYaxis()->SetTitle(yexpr);
359  if (!TString(drawopt).Contains("GOFF", TString::kIgnoreCase))
360  hist->Draw(drawopt);
361  return hist;
362  }
364  TProfile *drawProf(TString xexpr,
365  TString yexpr,
366  const char *cut = "",
367  TString drawopt = "",
368  const char *hname = "htemp",
369  TProfile *htemplate = nullptr) {
370  TProfile *hist = nullptr;
371  if (htemplate != nullptr) {
372  if ((strcmp(hname, "htemp") == 0) && (strcmp(hname, htemplate->GetName()) != 0))
373  htempDelete();
374  hist = (TProfile *)hist->Clone(hname);
375  } else if (drawopt.Contains("SAME", TString::kIgnoreCase)) {
376  hist = getSameProf(hname);
377  }
378 
379  // if in the end we found no way to make "hist"
380  if (hist == nullptr) {
381  if (strcmp(hname, "htemp") == 0)
382  htempDelete();
383  hist = new TProfile(hname, "", gEnv->GetValue("Hist.Binning.1D.x", 100), 0., 0.);
384  hist->SetCanExtend(TH1::kAllAxes);
385  }
386  return drawProf(xexpr, yexpr, cut, drawopt, hist);
387  }
388 
390  TProfile *drawProf(TString xexpr,
391  int bins,
392  double xlow,
393  double xhigh,
394  TString yexpr,
395  const char *cut = "",
396  const char *drawopt = "",
397  const char *hname = "htemp") {
398  if (TString(drawopt).Contains("SAME", TString::kIgnoreCase)) {
399  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
400  TProfile *hsame = getSameProf(hname);
401  return drawProf(xexpr, yexpr, cut, drawopt, hsame);
402  }
403  if (strcmp(hname, "htemp") == 0)
404  htempDelete();
405  TProfile *htemp = new TProfile(hname, "", bins, xlow, xhigh);
406  return drawProf(xexpr, yexpr, cut, drawopt, htemp);
407  }
408 
409  //------------------------------------------------------------------------------------------------------------------------------------
411  TH2 *draw2D(TString xexpr, TString yexpr, const char *cut, TString drawopt, TH2 *hist) {
412  // prep the machinery
413  helper::ScannerBase scanner(objType);
415  if (!scanner.addExpression((const char *)xexpr))
416  return nullptr;
417  if (!scanner.addExpression((const char *)yexpr))
418  return nullptr;
419  if (strlen(cut))
420  scanner.setCut(cut);
421 
422  // check histo
423  if (hist == nullptr) {
424  std::cerr << "Method draw2D(xexpr, yexpr, cut, drawopt, hist) cannot be called with null 'hist'. Use the other "
425  "draw methods instead."
426  << std::endl;
427  return nullptr;
428  }
429 
430  // fill histogram
431  int iev = 0;
432  for (event_->toBegin(), iev = 0; !event_->atEnd(); ++(*event_), ++iev) {
433  if (maxEvents_ > -1 && iev > maxEvents_)
434  break;
435  if (!selectEvent(*event_))
436  continue;
437  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
438  const Collection &vals = *handle_;
439  for (size_t j = 0, n = vals.size(); j < n; ++j) {
440  scanner.fill2D(&vals[j], hist);
441  }
442  }
443 
444  if (!strlen(hist->GetTitle()))
445  hist->SetTitle((strlen(cut) ? yexpr + ":" + xexpr + "{" + cut + "}" : yexpr + ":" + xexpr));
446  if (!strlen(hist->GetXaxis()->GetTitle()))
447  hist->GetXaxis()->SetTitle(xexpr);
448  if (!strlen(hist->GetYaxis()->GetTitle()))
449  hist->GetYaxis()->SetTitle(yexpr);
450  if (!TString(drawopt).Contains("GOFF", TString::kIgnoreCase))
451  hist->Draw(drawopt);
452  return hist;
453  }
456  TH2 *draw2D(TString xexpr,
457  TString yexpr,
458  const char *cut = "",
459  TString drawopt = "",
460  const char *hname = "htemp",
461  TH2 *htemplate = nullptr) {
462  TH2 *hist = nullptr;
463  if (htemplate != nullptr) {
464  if ((strcmp(hname, "htemp") == 0) && (strcmp(hname, htemplate->GetName()) != 0))
465  htempDelete();
466  hist = (TH2 *)hist->Clone(hname);
467  } else if (drawopt.Contains("SAME", TString::kIgnoreCase)) {
468  hist = getSameH2(hname);
469  }
470 
471  // if in the end we found no way to make "hist"
472  if (hist == nullptr) {
473  // prep the machinery
474  helper::ScannerBase scanner(objType);
476  if (!scanner.addExpression((const char *)xexpr))
477  return nullptr;
478  if (!scanner.addExpression((const char *)yexpr))
479  return nullptr;
480  if (strlen(cut))
481  scanner.setCut(cut);
482 
483  if (strcmp(hname, "htemp") == 0)
484  htempDelete();
485  // ok this is much more a hack than for the 1D case
486  double xmin = 0, xmax = -1, ymin = 0, ymax = -1;
487  int iev;
488  for (event_->toBegin(), iev = 0; !event_->atEnd(); ++(*event_), ++iev) {
489  if (maxEvents_ > -1 && iev > maxEvents_)
490  break;
491  if (!selectEvent(*event_))
492  continue;
493  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
494  const Collection &vals = *handle_;
495  for (size_t j = 0, n = vals.size(); j < n; ++j) {
496  if (!scanner.test(&vals[j]))
497  continue;
498  double x = scanner.eval(&vals[j], 0);
499  double y = scanner.eval(&vals[j], 1);
500  if ((xmax == -1) || (x >= xmax))
501  xmax = x;
502  if ((xmin == 0) || (x <= xmin))
503  xmin = x;
504  if ((ymax == -1) || (y >= ymax))
505  ymax = y;
506  if ((ymin == 0) || (y <= ymin))
507  ymin = y;
508  }
509  }
510  hist = new TH2F(hname,
511  "",
512  gEnv->GetValue("Hist.Binning.2D.x", 20),
513  xmin,
514  xmax,
515  gEnv->GetValue("Hist.Binning.2D.y", 20),
516  ymin,
517  ymax);
518  }
519  return draw2D(xexpr, yexpr, cut, drawopt, hist);
520  }
521 
523  TH2 *draw2D(TString xexpr,
524  int xbins,
525  double xlow,
526  double xhigh,
527  TString yexpr,
528  int ybins,
529  double ylow,
530  double yhigh,
531  const char *cut = "",
532  const char *drawopt = "",
533  const char *hname = "htemp") {
534  if (TString(drawopt).Contains("SAME", TString::kIgnoreCase)) {
535  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
536  TH2 *hsame = getSameH2(hname);
537  return draw2D(xexpr, yexpr, cut, drawopt, hsame);
538  }
539  if (strcmp(hname, "htemp") == 0)
540  htempDelete();
541  TH2 *htemp = new TH2F("htemp", "", xbins, xlow, xhigh, ybins, ylow, yhigh);
542  return draw2D(xexpr, yexpr, cut, drawopt, htemp);
543  }
544 
546  TGraph *drawGraph(TString xexpr, TString yexpr, const char *cut, TString drawopt, TGraph *graph) {
547  // prep the machinery
548  helper::ScannerBase scanner(objType);
550  if (!scanner.addExpression((const char *)xexpr))
551  return nullptr;
552  if (!scanner.addExpression((const char *)yexpr))
553  return nullptr;
554  if (strlen(cut))
555  scanner.setCut(cut);
556 
557  // make graph, if needed
558  if (graph == nullptr) {
559  graph = new TGraph();
560  graph->SetNameTitle("htemp", (strlen(cut) ? yexpr + ":" + xexpr + "{" + cut + "}" : yexpr + ":" + xexpr));
561  }
562 
563  // fill graph
564  int iev = 0;
565  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
566  if (maxEvents_ > -1 && iev > maxEvents_)
567  break;
568  if (!selectEvent(*event_))
569  continue;
570  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
571  const Collection &vals = *handle_;
572  for (size_t j = 0, n = vals.size(); j < n; ++j) {
573  scanner.fillGraph(&vals[j], graph);
574  }
575  }
576 
577  if (!strlen(graph->GetTitle()))
578  graph->SetTitle((strlen(cut) ? yexpr + ":" + xexpr + "{" + cut + "}" : yexpr + ":" + xexpr));
579  if (!strlen(graph->GetXaxis()->GetTitle()))
580  graph->GetXaxis()->SetTitle(xexpr);
581  if (!strlen(graph->GetYaxis()->GetTitle()))
582  graph->GetYaxis()->SetTitle(yexpr);
583  if (!TString(drawopt).Contains("GOFF", TString::kIgnoreCase))
584  graph->Draw(drawopt);
585  return graph;
586  }
587 
589  TGraph *drawGraph(
590  TString xexpr, TString yexpr, const char *cut = "", TString drawopt = "AP", const char *gname = "htemp") {
591  if (strcmp(gname, "htemp") == 0)
592  htempDelete();
593  TGraph *graph = new TGraph();
594  graph->SetNameTitle(gname, (strlen(cut) ? yexpr + ":" + xexpr + "{" + cut + "}" : yexpr + ":" + xexpr));
595  return drawGraph(xexpr, yexpr, cut, drawopt, graph);
596  }
597 
598  //------------------------------------------------------------------------------------------------------------------------------------
605  RooDataSet *fillDataSet(const char *realvars,
606  const char *boolvars,
607  const char *cut = "",
608  const char *name = "data") {
609  helper::ScannerBase scanner(objType);
611 
612  RooArgList vars;
613  TObjArray *exprArray = TString(realvars).Tokenize(exprSep_);
614  TObjArray *catArray = TString(boolvars).Tokenize(exprSep_);
615  int nreals = exprArray->GetEntries();
616  int nbools = catArray->GetEntries();
617  for (int i = 0; i < nreals; ++i) {
618  TString str = ((TObjString *)(*exprArray)[i])->GetString();
619  std::string lb = str.Data();
620  std::string ex = str.Data();
621  if ((ex[0] == '@') && (ex.find('=') != std::string::npos)) {
622  lb = lb.substr(1, ex.find('=') - 1);
623  ex = ex.substr(ex.find('=') + 1);
624  }
625  if (!scanner.addExpression(ex.c_str())) {
626  std::cerr << "Filed to define real variable '" << lb << "', expr = '" << ex << "'" << std::endl;
627  return nullptr;
628  }
629  // FIXME: I have to leave it dangling on the HEAP otherwise ROOT segfaults...
630  RooRealVar *var = new RooRealVar(lb.c_str(), lb.c_str(), 0.0);
631  vars.add(*var);
632  }
633  for (int i = 0; i < nbools; ++i) {
634  TString str = ((TObjString *)(*catArray)[i])->GetString();
635  std::string lb = str.Data();
636  std::string ex = str.Data();
637  if ((ex[0] == '@') && (ex.find('=') != std::string::npos)) {
638  lb = lb.substr(1, ex.find('=') - 1);
639  ex = ex.substr(ex.find('=') + 1);
640  }
641  if (!scanner.addExtraCut(ex.c_str())) {
642  std::cerr << "Filed to define bool variable '" << lb << "', cut = '" << ex << "'" << std::endl;
643  return nullptr;
644  }
645  RooCategory *cat = new RooCategory(lb.c_str(), lb.c_str());
646  cat->defineType("fail", 0);
647  cat->defineType("pass", 1);
648  vars.add(*cat);
649  }
650 
651  RooDataSet *ds = new RooDataSet(name, name, vars);
652 
653  if (strlen(cut))
654  scanner.setCut(cut);
655  int iev = 0;
656  for (event_->toBegin(); !event_->atEnd(); ++iev, ++(*event_)) {
657  if (maxEvents_ > -1 && iev > maxEvents_)
658  break;
659  if (!selectEvent(*event_))
660  continue;
661  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
662  if (handle_.failedToGet()) {
663  if (ignoreExceptions_)
664  continue;
665  }
666  const Collection &vals = *handle_;
667  for (size_t j = 0, n = vals.size(); j < n; ++j) {
668  if (!scanner.test(&vals[j]))
669  continue;
670  for (int i = 0; i < nreals; ++i) {
671  RooRealVar *var = (RooRealVar *)vars.at(i);
672  var->setVal(scanner.eval(&vals[j], i));
673  }
674  for (int i = 0; i < nbools; ++i) {
675  RooCategory *cat = (RooCategory *)vars.at(i + nreals);
676  cat->setIndex(int(scanner.test(&vals[j], i + 1))); // 0 is the event selection cut
677  }
678  ds->add(vars);
679  }
680  }
681 
682  delete exprArray;
683  delete catArray;
684 
685  return ds;
686  }
687 
690  void setIgnoreExceptions(bool ignoreThem) { ignoreExceptions_ = ignoreThem; }
691  void setMaxLinesToPrint(int lines) { maxLinesToPrint_ = (lines > 0 ? lines : 2147483647); }
692 
693  void addEventSelector(fwlite::EventSelector *selector) { eventSelectors_.Add(selector); }
695  TObjArray &eventSelectors() { return eventSelectors_; }
696  bool selectEvent(const fwlite::EventBase &ev) const {
697  for (int i = 0, n = eventSelectors_.GetEntries(); i < n; ++i) {
698  if (!((fwlite::EventSelector *)(eventSelectors_[i]))->accept(ev))
699  return false;
700  }
701  return true;
702  }
703 
704  void setMaxEvents(int max) { maxEvents_ = max; }
705 
706  private:
711  TString exprSep_;
714 
715  TObjArray eventSelectors_;
716 
718 
720  bool wantMore() const {
721  // ask if user wants more
722  fprintf(stderr, "Type <CR> to continue or q to quit ==> ");
723  // read first char
724  int readch = getchar(), answer = readch;
725  // poll out remaining chars from buffer
726  while (readch != '\n' && readch != EOF)
727  readch = getchar();
728  // check first char
729  return !(answer == 'q' || answer == 'Q');
730  }
731 
732  void htempDelete() {
733  if (gDirectory) {
734  TObject *obj = gDirectory->Get("htemp");
735  if (obj)
736  obj->Delete();
737  }
738  }
739 
742  TH1 *getSameH1(const char *hname) {
743  if (gDirectory && gDirectory->Get("htemp") != nullptr &&
744  gDirectory->Get("htemp")->IsA()->InheritsFrom(TH1::Class())) {
745  TH1 *hist = (TH1 *)((TH1 *)gDirectory->Get("htemp"))->Clone(hname);
746  hist->Reset();
747  hist->SetLineColor(kBlack);
748  hist->SetMarkerColor(kBlack);
749  return hist;
750  } else {
751  std::cerr << "There is no 'htemp' histogram from which to 'SAME'." << std::endl;
752  return nullptr;
753  }
754  }
755 
758  TH2 *getSameH2(const char *hname) {
759  if (gDirectory && gDirectory->Get("htemp") != nullptr &&
760  gDirectory->Get("htemp")->IsA()->InheritsFrom(TH2::Class())) {
761  TH2 *hist = (TH2 *)((TH2 *)gDirectory->Get("htemp"))->Clone(hname);
762  hist->Reset();
763  hist->SetLineColor(kBlack);
764  hist->SetMarkerColor(kBlack);
765  return hist;
766  } else {
767  std::cerr << "There is no 'htemp' histogram from which to 'SAME'." << std::endl;
768  return nullptr;
769  }
770  }
771 
774  TProfile *getSameProf(const char *hname) {
775  if (gDirectory && gDirectory->Get("htemp") != nullptr &&
776  gDirectory->Get("htemp")->IsA()->InheritsFrom(TProfile::Class())) {
777  TProfile *hist = (TProfile *)((TProfile *)gDirectory->Get("htemp"))->Clone(hname);
778  hist->Reset();
779  hist->SetLineColor(kBlack);
780  hist->SetMarkerColor(kBlack);
781  return hist;
782  } else {
783  std::cerr << "There is no 'htemp' histogram from which to 'SAME'." << std::endl;
784  return nullptr;
785  }
786  }
787  };
788 } // namespace fwlite
789 #endif // PhysicsTools_FWLite_Scanner_h
void clearEventSelector()
Definition: Scanner.h:694
bool ignoreExceptions_
Definition: Scanner.h:710
void fill2D(const void *obj, TH2 *hist2d) const
virtual edm::EventAuxiliary const & eventAuxiliary() const =0
string separator
Definition: mps_merge.py:79
void addEventSelector(fwlite::EventSelector *selector)
Definition: Scanner.h:693
TProfile * drawProf(TString xexpr, TString yexpr, const char *cut="", TString drawopt="", const char *hname="htemp", TProfile *htemplate=nullptr)
Just like draw() except that it uses TProfile. Note that the order is (x,y) while in ROOT it&#39;s usuall...
Definition: Scanner.h:364
const double xbins[]
void setPrintFullEventId(bool printIt=true)
Definition: Scanner.h:688
void fillGraph(const void *obj, TGraph *graph) const
static PFTauRenderPlugin instance
TObjArray eventSelectors_
Definition: Scanner.h:715
TH2 * draw2D(TString xexpr, TString yexpr, const char *cut="", TString drawopt="", const char *hname="htemp", TH2 *htemplate=nullptr)
Definition: Scanner.h:456
size_t count(const char *cut)
Definition: Scanner.h:154
void setIgnoreExceptions(bool ignoreThem)
void setMaxLinesToPrint(int lines)
Definition: Scanner.h:691
TH2 * getSameH2(const char *hname)
Definition: Scanner.h:758
TH2 * draw2D(TString xexpr, int xbins, double xlow, double xhigh, TString yexpr, int ybins, double ylow, double yhigh, const char *cut="", const char *drawopt="", const char *hname="htemp")
Just like draw() except that it uses TH2. Note that the order is (x,y) while in ROOT it&#39;s usually (y...
Definition: Scanner.h:523
void setExpressionSeparator(TString separator)
Definition: Scanner.h:689
bool printFullEventId_
Definition: Scanner.h:709
TH1 * draw(const char *expr, int nbins, double *xbins, const char *cut="", const char *drawopt="", const char *hname="htemp")
Definition: Scanner.h:298
TProfile * drawProf(TString xexpr, TString yexpr, const char *cut, TString drawopt, TProfile *hist)
Just like draw() except that it uses TProfile. Note that the order is (x,y) while in ROOT it&#39;s usuall...
Definition: Scanner.h:320
TString exprSep_
Definition: Scanner.h:711
HandleT handle_
Definition: Scanner.h:712
bool test(const void *obj, size_t icut=0) const
TH1 * getSameH1(const char *hname)
Definition: Scanner.h:742
std::string label_
Definition: Scanner.h:708
bool addExtraCut(const char *cut)
Add one extra cut that can be evaluated separately (as if it was an expression)
TGraph * drawGraph(TString xexpr, TString yexpr, const char *cut="", TString drawopt="AP", const char *gname="htemp")
Definition: Scanner.h:589
char const * label
fwlite::Handle< Collection > HandleT
The type of the Handle to read the Ts from the event. Needed to resolve its Type. ...
Definition: Scanner.h:49
TH1 * draw(const char *expr, int nbins, double xlow, double xhigh, const char *cut="", const char *drawopt="", const char *hname="htemp")
Definition: Scanner.h:276
def cat(path)
Definition: eostools.py:401
double eval(const void *obj, size_t iexpr=0) const
fwlite::Scanner<C>, a way to inspect or plots elements of a collection C by using the StringParser...
Definition: Scanner.h:46
edm::TypeWithDict objType
Definition: Scanner.h:713
size_t countEvents()
Definition: Scanner.h:178
fwlite::EventBase * event_
Definition: Scanner.h:707
void fill1D(const void *obj, TH1 *hist) const
void fillProf(const void *obj, TProfile *prof) const
bool wantMore() const
Definition: Scanner.h:720
TGraph * drawGraph(TString xexpr, TString yexpr, const char *cut, TString drawopt, TGraph *graph)
Definition: Scanner.h:546
Scanner(fwlite::EventBase *ev, const char *label, const char *instance="", const char *process="")
Definition: Scanner.h:52
bool addExpression(const char *expr)
void scan(const char *exprs, const char *cut="", int nmax=-1)
Definition: Scanner.h:80
void setIgnoreExceptions(bool ignoreThem)
Definition: Scanner.h:690
__shared__ Hist hist
std::string instance_
Definition: Scanner.h:708
RooDataSet * fillDataSet(const char *realvars, const char *boolvars, const char *cut="", const char *name="data")
Definition: Scanner.h:605
TProfile * getSameProf(const char *hname)
Definition: Scanner.h:774
TH1 * draw(const char *expr, const char *cut="", TString drawopt="", const char *hname="htemp", const TH1 *htemplate=nullptr)
Definition: Scanner.h:248
int maxLinesToPrint_
Definition: Scanner.h:719
static edm::TypeWithDict elementType(const edm::TypeWithDict &wrapperType)
Perform the type deduction form edm::Wrapper<C> to C::value_type and resolves typedefs.
HitContainer const *__restrict__ TkSoA const *__restrict__ Quality const *__restrict__ CAHitNtupletGeneratorKernelsGPU::HitToTuple const *__restrict__ int32_t int32_t int iev
TH1 * draw(const char *expr, const char *cut, TString drawopt, TH1 *hist)
Definition: Scanner.h:199
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=nullptr, const char *iProcessLabel=nullptr)
Definition: Handle.h:100
bool setCut(const char *cut)
Set the default cut that is applied to the events.
TProfile * drawProf(TString xexpr, int bins, double xlow, double xhigh, TString yexpr, const char *cut="", const char *drawopt="", const char *hname="htemp")
Just like draw() except that it uses TProfile. Note that the order is (x,y) while in ROOT it&#39;s usuall...
Definition: Scanner.h:390
TH2 * draw2D(TString xexpr, TString yexpr, const char *cut, TString drawopt, TH2 *hist)
Just like draw() except that it uses TH2. Note that the order is (x,y) while in ROOT it&#39;s usually (y...
Definition: Scanner.h:411
bool failedToGet() const
Returns true only if Handle was used in a &#39;get&#39; call and the data could not be found.
Definition: Handle.h:63
void setMaxEvents(int max)
Definition: Scanner.h:704
void print(const void *obj) const
float x
virtual bool atEnd() const =0
void htempDelete()
Definition: Scanner.h:732
vars
Definition: DeepTauId.cc:30
#define str(s)
virtual EventBase const & toBegin()=0
TObjArray & eventSelectors()
Definition: Scanner.h:695
std::string process_
Definition: Scanner.h:708
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:696
static std::type_info const & typeInfo()
Definition: Wrapper.h:43
__host__ __device__ void printIt(C *m, const char *prefix="")
Definition: FitUtils.h:53