1 #ifndef PhysicsTools_FWLite_Scanner_h
2 #define PhysicsTools_FWLite_Scanner_h
9 #include <TObjString.h>
10 #include <TObjArray.h>
11 #include <TDirectory.h>
15 #if !defined(__CINT__) && !defined(__MAKECINT__)
17 #include <RooArgList.h>
18 #include <RooDataSet.h>
19 #include <RooRealVar.h>
20 #include <RooCategory.h>
45 template <
typename Collection>
80 void scan(
const char *exprs,
const char *
cut =
"",
int nmax = -1) {
84 TObjArray *exprArray = TString(exprs).Tokenize(
exprSep_);
87 printf(
" : %9s : %4s : %9s : %3s",
"RUN",
"LUMI",
"EVENT",
"#IT");
88 rowline += 3 * 4 + 9 + 4 + 9 + 3 - 1;
90 printf(
" : %5s : %3s",
"EVENT",
"#IT");
91 rowline += 3 + 6 + 3 + 3 - 1;
93 for (
int i = 0;
i < exprArray->GetEntries(); ++
i) {
94 TString
str = ((TObjString *)(*exprArray)[
i])->GetString();
97 if ((ex[0] ==
'@') && (ex.find(
'=') != std::string::npos)) {
98 lb = lb.substr(1, ex.find(
'=') - 1);
99 ex = ex.substr(ex.find(
'=') + 1);
103 (lb.size() > 8 ? lb.substr(lb.size() - 8) : lb)
111 TString rule(
'-', rowline);
112 std::cout <<
" " << rule <<
" " << std::endl;
126 const Collection &vals = *
handle_;
127 for (
size_t j = 0,
n = vals.size();
j <
n; ++
j) {
128 if (!scanner.
test(&vals[
j]))
132 printf(
" : %9u : %4u : %9llu : %3lu",
id.
run(),
id.luminosityBlock(),
id.
event(), (
unsigned long)j);
134 printf(
" : %5d : %3lu", iev, (
unsigned long)j);
136 scanner.
print(&vals[j]);
168 const Collection &vals = *
handle_;
169 for (
size_t j = 0,
n = vals.size();
j <
n; ++
j) {
170 if (scanner.
test(&vals[
j]))
199 TH1 *
draw(
const char *expr,
const char *
cut, TString drawopt, TH1 *
hist) {
209 if (hist ==
nullptr) {
210 std::cerr <<
"Method draw(expr, cut, drawopt, hist) cannot be called with null 'hist'. Use the other draw "
224 const Collection &vals = *
handle_;
225 for (
size_t j = 0,
n = vals.size();
j <
n; ++
j) {
226 scanner.
fill1D(&vals[
j], hist);
230 if (drawopt.Contains(
"NORM", TString::kIgnoreCase) && (hist->Integral() != 0)) {
232 hist->Scale(1.0 / hist->Integral());
234 drawopt(TRegexp(
"[Nn][Oo][Rr][Mm]")) =
"";
237 if (!drawopt.Contains(
"GOFF", TString::kIgnoreCase))
249 const char *
cut =
"",
250 TString drawopt =
"",
251 const char *hname =
"htemp",
252 const TH1 *htemplate =
nullptr) {
254 if (htemplate !=
nullptr) {
255 if ((strcmp(hname,
"htemp") == 0) && (strcmp(hname, htemplate->GetName()) != 0))
257 hist = (TH1 *)htemplate->Clone(hname);
258 }
else if (drawopt.Contains(
"SAME", TString::kIgnoreCase)) {
263 if (hist ==
nullptr) {
264 if (strcmp(hname,
"htemp") == 0)
266 hist =
new TH1F(hname,
"", gEnv->GetValue(
"Hist.Binning.1D.x", 100), 0, 0);
267 hist->SetCanExtend(TH1::kAllAxes);
269 hist->SetTitle((strlen(
cut) ? TString(expr) +
"{" +
cut +
"}" : TString(expr)));
270 hist->GetXaxis()->SetTitle(expr);
271 return draw(expr,
cut, drawopt, hist);
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;
286 return draw(expr,
cut, drawopt, hsame);
288 if (strcmp(hname,
"htemp") == 0)
290 TH1 *htemp =
new TH1F(hname, expr, nbins, xlow, xhigh);
292 htemp->SetTitle(TString(expr) +
"{" +
cut +
"}");
293 htemp->GetXaxis()->SetTitle(expr);
294 return draw(expr,
cut, drawopt, htemp);
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;
307 return draw(expr,
cut, drawopt, hsame);
309 if (strcmp(hname,
"htemp") == 0)
311 TH1 *htemp =
new TH1F(hname, expr, nbins, xbins);
313 htemp->SetTitle(TString(expr) +
"{" +
cut +
"}");
314 htemp->GetXaxis()->SetTitle(expr);
315 return draw(expr,
cut, drawopt, htemp);
320 TProfile *
drawProf(TString xexpr, TString yexpr,
const char *
cut, TString drawopt, TProfile *
hist) {
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."
347 const Collection &vals = *
handle_;
348 for (
size_t j = 0,
n = vals.size();
j <
n; ++
j) {
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))
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))
374 hist = (TProfile *)hist->Clone(hname);
375 }
else if (drawopt.Contains(
"SAME", TString::kIgnoreCase)) {
380 if (hist ==
nullptr) {
381 if (strcmp(hname,
"htemp") == 0)
383 hist =
new TProfile(hname,
"", gEnv->GetValue(
"Hist.Binning.1D.x", 100), 0., 0.);
384 hist->SetCanExtend(TH1::kAllAxes);
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;
401 return drawProf(xexpr, yexpr,
cut, drawopt, hsame);
403 if (strcmp(hname,
"htemp") == 0)
405 TProfile *htemp =
new TProfile(hname,
"", bins, xlow, xhigh);
406 return drawProf(xexpr, yexpr,
cut, drawopt, htemp);
411 TH2 *
draw2D(TString xexpr, TString yexpr,
const char *
cut, TString drawopt, TH2 *
hist) {
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."
438 const Collection &vals = *
handle_;
439 for (
size_t j = 0,
n = vals.size();
j <
n; ++
j) {
440 scanner.
fill2D(&vals[
j], hist);
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))
458 const char *
cut =
"",
459 TString drawopt =
"",
460 const char *hname =
"htemp",
461 TH2 *htemplate =
nullptr) {
463 if (htemplate !=
nullptr) {
464 if ((strcmp(hname,
"htemp") == 0) && (strcmp(hname, htemplate->GetName()) != 0))
466 hist = (TH2 *)hist->Clone(hname);
467 }
else if (drawopt.Contains(
"SAME", TString::kIgnoreCase)) {
472 if (hist ==
nullptr) {
483 if (strcmp(hname,
"htemp") == 0)
494 const Collection &vals = *
handle_;
495 for (
size_t j = 0,
n = vals.size();
j <
n; ++
j) {
496 if (!scanner.
test(&vals[
j]))
498 double x = scanner.
eval(&vals[j], 0);
499 double y = scanner.
eval(&vals[j], 1);
502 if ((xmin == 0) || (x <= xmin))
510 hist =
new TH2F(hname,
512 gEnv->GetValue(
"Hist.Binning.2D.x", 20),
515 gEnv->GetValue(
"Hist.Binning.2D.y", 20),
519 return draw2D(xexpr, yexpr,
cut, drawopt, hist);
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;
537 return draw2D(xexpr, yexpr,
cut, drawopt, hsame);
539 if (strcmp(hname,
"htemp") == 0)
541 TH2 *htemp =
new TH2F(
"htemp",
"", xbins, xlow, xhigh, ybins, ylow, yhigh);
542 return draw2D(xexpr, yexpr,
cut, drawopt, htemp);
546 TGraph *
drawGraph(TString xexpr, TString yexpr,
const char *
cut, TString drawopt, TGraph *graph) {
558 if (graph ==
nullptr) {
559 graph =
new TGraph();
560 graph->SetNameTitle(
"htemp", (strlen(cut) ? yexpr +
":" + xexpr +
"{" + cut +
"}" : yexpr +
":" + xexpr));
571 const Collection &vals = *
handle_;
572 for (
size_t j = 0,
n = vals.size();
j <
n; ++
j) {
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);
590 TString xexpr, TString yexpr,
const char *
cut =
"", TString drawopt =
"AP",
const char *gname =
"htemp") {
591 if (strcmp(gname,
"htemp") == 0)
593 TGraph *graph =
new TGraph();
594 graph->SetNameTitle(gname, (strlen(
cut) ? yexpr +
":" + xexpr +
"{" +
cut +
"}" : yexpr +
":" + xexpr));
606 const char *boolvars,
607 const char *
cut =
"",
608 const char *
name =
"data") {
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();
621 if ((ex[0] ==
'@') && (ex.find(
'=') != std::string::npos)) {
622 lb = lb.substr(1, ex.find(
'=') - 1);
623 ex = ex.substr(ex.find(
'=') + 1);
626 std::cerr <<
"Filed to define real variable '" << lb <<
"', expr = '" << ex <<
"'" << std::endl;
630 RooRealVar *
var =
new RooRealVar(lb.c_str(), lb.c_str(), 0.0);
633 for (
int i = 0;
i < nbools; ++
i) {
634 TString
str = ((TObjString *)(*catArray)[
i])->GetString();
637 if ((ex[0] ==
'@') && (ex.find(
'=') != std::string::npos)) {
638 lb = lb.substr(1, ex.find(
'=') - 1);
639 ex = ex.substr(ex.find(
'=') + 1);
642 std::cerr <<
"Filed to define bool variable '" << lb <<
"', cut = '" << ex <<
"'" << std::endl;
645 RooCategory *
cat =
new RooCategory(lb.c_str(), lb.c_str());
646 cat->defineType(
"fail", 0);
647 cat->defineType(
"pass", 1);
651 RooDataSet *ds =
new RooDataSet(
name,
name, vars);
666 const Collection &vals = *
handle_;
667 for (
size_t j = 0,
n = vals.size();
j <
n; ++
j) {
668 if (!scanner.
test(&vals[
j]))
670 for (
int i = 0;
i < nreals; ++
i) {
671 RooRealVar *
var = (RooRealVar *)vars.at(
i);
672 var->setVal(scanner.
eval(&vals[j],
i));
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)));
722 fprintf(stderr,
"Type <CR> to continue or q to quit ==> ");
724 int readch = getchar(), answer = readch;
726 while (readch !=
'\n' && readch != EOF)
729 return !(answer ==
'q' || answer ==
'Q');
734 TObject *
obj = gDirectory->Get(
"htemp");
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);
747 hist->SetLineColor(kBlack);
748 hist->SetMarkerColor(kBlack);
751 std::cerr <<
"There is no 'htemp' histogram from which to 'SAME'." << std::endl;
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);
763 hist->SetLineColor(kBlack);
764 hist->SetMarkerColor(kBlack);
767 std::cerr <<
"There is no 'htemp' histogram from which to 'SAME'." << std::endl;
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);
779 hist->SetLineColor(kBlack);
780 hist->SetMarkerColor(kBlack);
783 std::cerr <<
"There is no 'htemp' histogram from which to 'SAME'." << std::endl;
789 #endif // PhysicsTools_FWLite_Scanner_h
void clearEventSelector()
void print(const void *obj) const
virtual edm::EventAuxiliary const & eventAuxiliary() const =0
void addEventSelector(fwlite::EventSelector *selector)
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's usuall...
void setPrintFullEventId(bool printIt=true)
void fill2D(const void *obj, TH2 *hist2d) const
static PFTauRenderPlugin instance
TObjArray eventSelectors_
TH2 * draw2D(TString xexpr, TString yexpr, const char *cut="", TString drawopt="", const char *hname="htemp", TH2 *htemplate=nullptr)
size_t count(const char *cut)
void setIgnoreExceptions(bool ignoreThem)
void setMaxLinesToPrint(int lines)
TH2 * getSameH2(const char *hname)
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's usually (y...
void setExpressionSeparator(TString separator)
TH1 * draw(const char *expr, int nbins, double *xbins, const char *cut="", const char *drawopt="", const char *hname="htemp")
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's usuall...
TH1 * getSameH1(const char *hname)
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")
fwlite::Handle< Collection > HandleT
The type of the Handle to read the Ts from the event. Needed to resolve its Type. ...
TH1 * draw(const char *expr, int nbins, double xlow, double xhigh, const char *cut="", const char *drawopt="", const char *hname="htemp")
fwlite::Scanner<C>, a way to inspect or plots elements of a collection C by using the StringParser...
edm::TypeWithDict objType
fwlite::EventBase * event_
list var
if using global norm cols_to_minmax = ['t_delta', 't_hmaxNearP','t_emaxNearP', 't_hAnnular', 't_eAnnular','t_pt','t_nVtx','t_ieta','t_eHcal10', 't_eHcal30','t_rhoh','t_eHcal'] df[cols_to_minmax] = df[cols_to_minmax].apply(lambda x: (x - x.min()) / (x.max() - x.min()) if (x.max() - x.min() > 0) else 1.0/200.0)
printf("params %d %f %f %f\n", minT, eps, errmax, chi2max)
bool test(const void *obj, size_t icut=0) const
double eval(const void *obj, size_t iexpr=0) const
bool failedToGet() const
Returns true only if Handle was used in a 'get' call and the data could not be found.
void fill1D(const void *obj, TH1 *hist) const
TGraph * drawGraph(TString xexpr, TString yexpr, const char *cut, TString drawopt, TGraph *graph)
Scanner(fwlite::EventBase *ev, const char *label, const char *instance="", const char *process="")
bool addExpression(const char *expr)
void scan(const char *exprs, const char *cut="", int nmax=-1)
void setIgnoreExceptions(bool ignoreThem)
RooDataSet * fillDataSet(const char *realvars, const char *boolvars, const char *cut="", const char *name="data")
TProfile * getSameProf(const char *hname)
TH1 * draw(const char *expr, const char *cut="", TString drawopt="", const char *hname="htemp", const TH1 *htemplate=nullptr)
static edm::TypeWithDict elementType(const edm::TypeWithDict &wrapperType)
Perform the type deduction form edm::Wrapper<C> to C::value_type and resolves typedefs.
TH1 * draw(const char *expr, const char *cut, TString drawopt, TH1 *hist)
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=nullptr, const char *iProcessLabel=nullptr)
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's usuall...
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's usually (y...
void fillProf(const void *obj, TProfile *prof) const
void setMaxEvents(int max)
virtual bool atEnd() const =0
static std::string const separator(":")
HitContainer const *__restrict__ TkSoA const *__restrict__ Quality const *__restrict__ CAHitNtupletGeneratorKernelsGPU::HitToTuple const *__restrict__ int32_t int iev
bool selectEvent(const fwlite::EventBase &ev) const
virtual EventBase const & toBegin()=0
void fillGraph(const void *obj, TGraph *graph) const
TObjArray & eventSelectors()
static std::type_info const & typeInfo()
__host__ __device__ void printIt(C *m, const char *prefix="")