FairRoot/PandaRoot
cutfinderx.C
Go to the documentation of this file.
1 #include <algorithm>
2 #include <map>
3 #include <utility>
4 #include <iostream>
5 
6 #include "TFile.h"
7 #include "TTree.h"
8 #include "TLeaf.h"
9 #include "TString.h"
10 #include "TH1F.h"
11 #include "TCanvas.h"
12 #include "TROOT.h"
13 #include "TEventList.h"
14 #include "TDirectory.h"
15 #include "TLine.h"
16 #include "TLatex.h"
17 #include "TStyle.h"
18 #include "TObjArray.h"
19 #include "TPRegexp.h"
20 #include "TRegexp.h"
21 #include "TGraph.h"
22 
23 #define MAX 1000
24 #define BINS 100
25 #define NCAN 3
26 
27 typedef std::vector<pair<double, int> > ValueMap;
28 typedef std::map<int, int> CountMap;
29 
31 std::map<TString, TString> mctvar;
32 
34 
35 double signi[MAX];
36 int idx[MAX];
39 double cut[MAX];
40 TGraph grl[MAX];
41 TGraph grr[MAX];
42 
43 double minh[MAX], maxh[MAX];
44 
45 double N0_sig, N0_bg;
46 double Nsigev, Nbgev;
47 bool dstarmode;
48 
49 Float_t fbranch[MAX];
50 Int_t ibranch[MAX];
52 Int_t ev, run, mode, rec, nbranch;
53 
54 // ---------------------------------------------------------------
55 
56 int gettype(TTree *t, TString varname)
57 {
58  TString leaftype = t->GetLeaf(varname)->GetTypeName();
59 
60  if (leaftype=="Float_t") return 0;
61  else if (leaftype=="Int_t") return 1;
62  else if (leaftype=="Bool_t") return 2;
63 
64  return -1;
65 }
66 
67 // ---------------------------------------------------------------
68 
69 int init(TTree *t)
70 {
71  t->SetBranchAddress("ev",&ev);
72  t->SetBranchAddress("run",&run);
73  t->SetBranchAddress("mode",&mode);
74  t->SetBranchAddress("recmode",&rec);
75 
76  TObjArray* branches = t->GetListOfBranches();
77 
78  nbranch = 0;
79 
80  for(int i=0; i<=branches->GetLast(); ++i)
81  {
82  TBranch* branch = (TBranch*)branches->UncheckedAt(i);
83  TString v=branch->GetName();
84 
85  if ( v=="ev" || v=="mode" || v=="run" || v=="nsig" ) continue;
86  if ( v=="essumpt" || v=="essumptc") continue;
87  if ( v.Contains("pdg") || v.Contains("beam") || v.Contains("mct") ) continue;
88  if ( v.BeginsWith("t") && v!="thr") continue;
89  if ( v.EndsWith("vx") || v.EndsWith("vy") || v.EndsWith("vz") || v.EndsWith("pocmag")) continue;
90 
91  bool ok=false;
92 
93  if ( v.BeginsWith("es") || v=="mmiss" || v.EndsWith("d0m") ) ok=true;
94  if ( v.EndsWith("p") || v.EndsWith("tht") || v.EndsWith("pcm") || v.EndsWith("thtcm") || v.EndsWith("pt") ) ok=true;
95  if ( ok || v.EndsWith("ang") || v.Contains("poc") ) ok=true;
96  if ( ok || v.Contains("pid") || v.Contains("min") || v.Contains("max") || v.Contains("sum") || v.Contains("fw") ) ok=true;
97  if ( ok || v.EndsWith("sph") || v.EndsWith("apl") || v.EndsWith("pla") || v.EndsWith("thr") || v.EndsWith("cir") ) ok=true;
98 
99  if ( v=="xmdif" && dstarmode ) ok = true;
100  if ( v=="xd0m" && dstarmode ) ok = false;
101 
102  if (!ok) continue;
103 
104  Int_t dtype=gettype(t, v);
105  switch (dtype)
106  {
107  case 0: t->SetBranchAddress(v, &(fbranch[nbranch])); break;
108  case 1: t->SetBranchAddress(v, &(ibranch[nbranch])); break;
109  case 2: t->SetBranchAddress(v, &(bbranch[nbranch])); break;
110  }
111 
112  vars[nbranch] = v;
113  nbranch++;
114  }
115 
116  return nbranch;
117 }
118 
119 // ---------------------------------------------------------------
120 
121 bool mycompare(int i, int j)
122 {
123  return signi[i]>signi[j];
124 }
125 
126 // ---------------------------------------------------------------
127 
128 int uid(int lev, int lrun, int lmode)
129 {
130  return lev+10000*lrun+(((lmode/100)%10)*20+lmode%10)*100000;
131 }
132 
133 // ---------------------------------------------------------------
134 
135 int countEvents(TTree *t, TEventList &el)
136 {
137  t->SetBranchStatus("*",0);
138  t->SetBranchStatus("ev",1);
139  t->SetBranchStatus("run",1);
140  t->SetBranchStatus("mode",1);
141  t->SetBranchStatus("recmode",1);
142 
143  evcnt.clear();
144  for (int j=0;j<10;++j) evcntrec[j].clear();
145 
146  for (int i=0;i<el.GetN();++i)
147  {
148  t->GetEntry(el.GetEntry(i));
149  evcnt[uid(ev,run,mode)]+=1;
150  if (rec<10) evcntrec[rec][uid(ev,run,mode)]+=1;
151  }
152  t->SetBranchStatus("*",1);
153 
154  return evcnt.size();
155 }
156 
157 // ---------------------------------------------------------------
158 
159 void makeMaps(TTree *t, TString varname, TEventList &els, TEventList &elb, ValueMap &sig, ValueMap &bg, int id)
160 {
161  int i;
162  t->SetBranchStatus("*",0);
163  t->SetBranchStatus("ev",1);
164  t->SetBranchStatus("mode",1);
165  t->SetBranchStatus("run",1);
166  t->SetBranchStatus(varname,1);
167 
168  Float_t var;
169  Int_t dtype=gettype(t, varname);
170 
171  sig.clear();
172  bg.clear();
173  double min = 1e9, max=-1e9;
174 
175  // prepare the pairs of variable value, eventnumber for signal and count signals
176  for (i=0;i<els.GetN();++i)
177  {
178  t->GetEntry(els.GetEntry(i));
179  switch (dtype) {case 2: var = (Float_t)bbranch[id]; break; case 1: var = (Float_t) ibranch[id]; break; default: var=fbranch[id]; }
180 
181  sig.push_back(std::make_pair(var, uid(ev,run,mode) ));
182  if (var>max) max = var;
183  if (var<min) min = var;
184  }
185 
186  // prepare the pairs of variable value, eventnumber for background and count background
187  for (i=0;i<elb.GetN();++i)
188  {
189  t->GetEntry(elb.GetEntry(i));
190  switch (dtype) {case 2: var = (Float_t)bbranch[id]; break; case 1: var = (Float_t) ibranch[id]; break; default: var=fbranch[id]; }
191 
192  bg.push_back(std::make_pair(var, uid(ev,run,mode) ));
193  if (var>max) max = var;
194  if (var<min) min = var;
195  }
196 
197  if (max>25) max=25;
198  minh[id]=min-(max-min)*0.05;
199  maxh[id]=max+(max-min)*0.05;
200 }
201 
202 // ---------------------------------------------------------------
203 
204 double bestEffEvt(TTree *t, TString varname, TEventList &els, TEventList &elb, double &bestcut, double supr, int id)
205 {
206  CountMap bgcnt, sigcnt, sigcnt2;
207  ValueMap sigvals, bgvals;
208 
209  int Nbgval = elb.GetN();
210  int dtype = gettype(t,varname);
211 
212  makeMaps(t, varname, els, elb, sigvals, bgvals, id);
213 
214 
215  // sort signals by first element = variable value
216  sort(bgvals.begin(), bgvals.end());
217 
218  int i=0;
219  bgcnt.clear();
220  // find cut for current efficiency requirement left cut
221  while ( (bgcnt.size()/(double)Nbgev)<(1-supr) && i<Nbgval ) bgcnt[bgvals[i++].second]+=1;
222  double leftcut = bgvals[i].first;
223 
224  bgcnt.clear();
225  i=bgvals.size()-1;
226  // find cut for current efficiency requirement right cut
227  while ( (bgcnt.size()/(double)Nbgev)<(1-supr) && i>=0) bgcnt[bgvals[i--].second]+=1;
228  double rightcut = bgvals[i].first;
229 
230  sigcnt.clear();
231  sigcnt2.clear();
232 
233  for (i=0;i<(int)sigvals.size();++i)
234  {
235  if (sigvals[i].first<leftcut) sigcnt[sigvals[i].second]+=1;
236  if (sigvals[i].first>rightcut) sigcnt2[sigvals[i].second]+=1;
237  }
238 
239  t->SetBranchStatus("*",1);
240 
241  double lefteff = sigcnt.size()/Nsigev;
242  double righteff = sigcnt2.size()/Nsigev;
243 
244  cout <<varname<<" "<<flush;
245 
246  bestcut = rightcut;
247  if (lefteff>righteff)
248  {
249  if (dtype==0) cuts[id] = TString::Format("%s<%.3f",varname.Data(),leftcut);
250  else cuts[id] = TString::Format("%s<=%.0f",varname.Data(),leftcut);
251 
252  bestcut = leftcut;
253  return lefteff;
254  }
255 
256  if (dtype==0) cuts[id] = TString::Format("%s>%.3f",varname.Data(),rightcut);
257  else cuts[id] = TString::Format("%s>=%.0f",varname.Data(),rightcut);
258 
259  return righteff;
260 }
261 
262 // ---------------------------------------------------------------
263 
264 double bestSuppressionEvt(TTree *t, TString varname, TEventList &els, TEventList &elb, double &bestcut, double eff, int id)
265 {
266  CountMap bgcnt, bgcnt2, sigcnt;
267  ValueMap sigvals, bgvals;
268 
269  int Nsigval = els.GetN();
270  int dtype = gettype(t,varname);
271 
272  makeMaps(t, varname, els, elb, sigvals, bgvals, id);
273 
274  // sort signals by first element = variable value
275  sort(sigvals.begin(), sigvals.end());
276 
277  int i=0;
278  sigcnt.clear();
279  // find cut for current efficiency requirement left cut
280  while ( (sigcnt.size()/(double)Nsigev)<eff && i<Nsigval ) sigcnt[sigvals[i++].second]+=1;
281  double leftcut = sigvals[i].first;
282 
283  sigcnt.clear();
284  i=sigvals.size()-1;
285  // find cut for current efficiency requirement right cut
286  while ( (sigcnt.size()/(double)Nsigev)<eff && i>=0) sigcnt[sigvals[i--].second]+=1;
287  double rightcut = sigvals[i].first;
288 
289  bgcnt.clear();
290  bgcnt2.clear();
291 
292  for (i=0;i<(int)bgvals.size();++i)
293  {
294  if (bgvals[i].first<=leftcut) bgcnt[bgvals[i].second];
295  if (bgvals[i].first>=rightcut) bgcnt2[bgvals[i].second];
296  }
297 
298  t->SetBranchStatus("*",1);
299 
300 
301  double leftsupr = (Nbgev-bgcnt.size())/Nbgev;
302  double rightsupr = (Nbgev-bgcnt2.size())/Nbgev;
303 
304  cout <<varname<<" "<<flush;
305 
306  if (leftcut!=leftcut || rightcut!=rightcut) return 0;
307 
308  bestcut = rightcut;
309  if (leftsupr>rightsupr)
310  {
311  if (dtype==0) cuts[id] = TString::Format("%s<%.5f",varname.Data(),leftcut);
312  else cuts[id] = TString::Format("%s<=%.0f",varname.Data(),leftcut);
313 
314  bestcut = leftcut;
315  return leftsupr;
316  }
317 
318  if (dtype==0) cuts[id] = TString::Format("%s>%.5f",varname.Data(),rightcut);
319  else cuts[id] = TString::Format("%s>=%.0f",varname.Data(),rightcut);
320 
321  return rightsupr;
322 }
323 // ---------------------------------------------------------------
324 
325 double bestCombiEvt(TTree *t, TString varname, TEventList &els, TEventList &elb, double &bestcut, double eff, int id)
326 {
327  CountMap bgcntl, bgcntr, sigcntl, sigcntr;
328  ValueMap sigvals, bgvals;
329 
330  int dtype = gettype(t,varname);
331  int Ns = els.GetN(), Nb = elb.GetN();
332 
333  makeMaps(t, varname, els, elb, sigvals, bgvals, id);
334 
335  // sort signals by first element = variable value
336  sort(sigvals.begin(), sigvals.end());
337  sort(bgvals.begin(), bgvals.end());
338 
339  double step = (maxh[id]-minh[id])/(double)BINS;
340 
341  int isigl=0, ibgl = 0, isigr = sigvals.size()-1, ibgr = bgvals.size()-1;
342 
343  double leftcut = 0., rightcut = 0., qal = 1000., qar = 1000.;
344  cout <<varname<<" "<<flush;
345 
346  for (int i=0; i<BINS; ++i)
347  {
348  double lcut = (i+1)*step + minh[id];
349  double rcut = maxh[id] - (i+1)*step;
350 
351  while (sigvals[isigl].first<lcut && isigl<Ns) sigcntl[sigvals[isigl++].second]+=1;
352  while (bgvals[ibgl].first<lcut && ibgl<Nb) bgcntl[bgvals[ibgl++].second]+=1;
353 
354  while (sigvals[isigr].first>rcut && isigr>=0) sigcntr[sigvals[isigr--].second]+=1;
355  while (bgvals[ibgr].first>rcut && ibgr>=0) bgcntr[bgvals[ibgr--].second]+=1;
356 
357  double effl = (double) sigcntl.size()/Nsigev;
358  double effr = (double) sigcntr.size()/Nsigev;
359  double supl = 1. - (double) bgcntl.size()/Nbgev;
360  double supr = 1. - (double) bgcntr.size()/Nbgev;
361 
362  double distl = sqrt((1.-effl)*(1.-effl) + (1.-supl)*(1.-supl));
363  double distr = sqrt((1.-effr)*(1.-effr) + (1.-supr)*(1.-supr));
364 
365  printf("%8s (%3d) : left(%5.3f / %5.3f : %5.3f) right(%5.3f / %5.3f : %5.3f)\n", varname.Data(), i, effl, supl, distl, effr, supr, distr);
366  //cout <<i<<" "<<"el:"<<effl<<" eb:"<<effr<<" sl:"<<supl<<" sr:"<<supr<<" dl:"<<distl<<" dr:"<<distr<<endl;
367 
368  grl[id].SetPoint(i, effl, supl);
369  grr[id].SetPoint(i, effr, supr);
370 
371  if (distl<qal) {qal=distl; leftcut=lcut;}
372  if (distr<qar) {qar=distr; rightcut=rcut;}
373  }
374 
375  t->SetBranchStatus("*",1);
376 
377  if (leftcut!=leftcut || rightcut!=rightcut) return 0;
378 
379  bestcut = rightcut;
380  if (qal<qar)
381  {
382  if (dtype==0) cuts[id] = TString::Format("%s<%.5f",varname.Data(),leftcut);
383  else cuts[id] = TString::Format("%s<=%.0f",varname.Data(),leftcut);
384 
385  bestcut = leftcut;
386  return qal;
387  }
388 
389  if (dtype==0) cuts[id] = TString::Format("%s>%.5f",varname.Data(),rightcut);
390  else cuts[id] = TString::Format("%s>=%.0f",varname.Data(),rightcut);
391 
392  return qar;
393 }
394 
395 // ---------------------------------------------------------------
396 
397 int cutfinderx(TString fname, TString precut="", double supr=0.95, int evmult=10000, double norm=1.0, int n0s=-1)
398 {
399  isbt = gROOT->IsBatch();
400 
401  gStyle->SetTitleX(0.2);
402  gStyle->SetTitleY(0.993);
403  gStyle->SetTitleH(0.07);
404  gStyle->SetPadTopMargin(0.08);
405  gStyle->SetTitleBorderSize(0);
406  gStyle->SetOptStat(0);
407 
408  TString tagcut = "tag";
409 
410  TRegexp regntp("n[0-9]+");
411  TRegexp regS("[0-9]+S");
412  TRegexp regB("[0-9]+B");
413 
414  TString s=fname(regS);
415  s = s(0,s.Length()-1);
416  TString b=fname(regB);
417  b = b(0,b.Length()-1);
418 
419  TString ntp = fname(regntp);
420  TString smode = ntp(1,ntp.Length());
421  int mode = smode.Atoi();
422  dstarmode = (mode>=110 && mode<=119) || (mode>=130 && mode<=138) || (mode>=150 && mode<=151);
423 
424  if (n0s<0) n0s = s.Atoi();
425  int n0b = b.Atoi();
426 
427  N0_sig = n0s*evmult;
428  N0_bg = n0b*evmult;
429 
430  cout <<"tree '"<<ntp.Data()<<"' with N_sig = "<<N0_sig<<", N_bg = "<<N0_bg<<endl;
431 
432  int i,j;
433 
434  TFile *f=new TFile(fname,"READ");
435  TTree *t=(TTree*)f->Get(ntp);
436 
437  init(t);
438 
439  TEventList els("els"), elsall("elsall"), elb("elb"), elball("elball");
440 
441  TCanvas *c1=new TCanvas("c1","c1",10,10,1800,550);
442  c1->Divide(7,3);
443 
444  if (precut=="") precut = tagcut;
445  else precut = tagcut+"&&"+precut;
446 
447  TString sigcut = "xmct";
448  TString bgcut = "!xmct";
449 
450  t->Draw(">>elsall",tagcut+"&&"+sigcut);
451  t->Draw(">>elball",tagcut+"&&"+bgcut);
452 
453  if (precut!="")
454  {
455  sigcut+="&&"+precut;
456  bgcut +="&&"+precut;
457  }
458 
459  cout <<sigcut <<" "<<bgcut<<endl;
460  t->Draw(">>els",sigcut);
461  t->Draw(">>elb",bgcut);
462 
463  float nsig = countEvents(t,elsall);
464  float nbg = countEvents(t,elball);
465  Nbgev = countEvents(t, elb);
466  Nsigev = countEvents(t, els);
467 
468  cout <<"SIG EVT: "<<Nsigev<<" ev "<<els.GetN()<<" cn BG: "<<Nbgev<<" ev "<<elb.GetN()<<" cn"<<endl;
469 
470  int cnt=1;
471 
472  for(i=0;i<nbranch;++i) { idx[i]=i; signi[i]=0.; grl[i].Set(BINS); grr[i].Set(BINS);}
473 
474  for(i=0; i<nbranch; ++i)
475  {
476  if (supr<0)
477  signi[i] = bestEffEvt(t, vars[i], els, elb, cut[i], -supr, i);
478  else if (supr>0)
479  signi[i] = bestSuppressionEvt(t, vars[i], els, elb, cut[i], supr, i);
480  else
481  signi[i] = bestCombiEvt(t, vars[i], els, elb, cut[i], supr, i);
482 
483  }
484  cout <<endl;
485 
486  std::vector<int> myidx (idx, idx+nbranch);
487  std::sort(myidx.begin(), myidx.end(), mycompare);
488 
489  if (!isbt)
490  {
491 
492  cout <<"\n\nBEST 20 vars:"<<endl<<endl;
493 
494  TLine l;
495  l.SetLineColor(6);
496  l.SetLineStyle(2);
497  l.SetLineWidth(2);
498  TLatex lt;
499  lt.SetTextSize(0.06);
500  TString target="supr";
501  if (supr<0) target="eff";
502  if (supr==0) target="qa";
503 
504  for (j=0;j<21;++j)
505  {
506  i=myidx[j];
507 
508  if (j<20)
509  {
510  printf("%2d) %-15s : %4s = %5.3f cut = %.4f ( %s )\n",j,vars[i].Data(), target.Data(), signi[i], cut[i], cuts[i].Data());
511 
512  TH1F h1("h1",vars[i],BINS,minh[i],maxh[i]);
513  TH1F h2("h2",vars[i],BINS,minh[i],maxh[i]);
514  h2.SetLineColor(2);
515 
516  t->SetEventList(&els);
517  t->Project("h1",vars[i]);
518  t->SetEventList(&elb);
519  t->Project("h2",vars[i]);
520 
521  h1.Scale(1.0/h1.Integral());
522  h2.Scale(1.0/h2.Integral());
523  h1.SetTitleSize(0.05);
524 
525  c1->cd(j+1);
526  double maxi = h1.GetMaximum();
527  if (h2.GetMaximum()>maxi ) maxi = h2.GetMaximum();
528  maxi*=1.1;
529  h1.SetMaximum(maxi);
530  h2.SetMaximum(maxi);
531 
532  h1.DrawNormalized();
533  h2.DrawNormalized("same");
534 
535  double axmin = h1.GetXaxis()->GetXmin(),axmax = h1.GetXaxis()->GetXmax();
536 
537  l.DrawLine(cut[i],0, cut[i], maxi*0.9);
538 
539  lt.DrawLatex(axmin+(axmax-axmin)*0.6,1.01*maxi,TString::Format("%s = %6.4f",target.Data(),signi[i]));
540 
541  }
542  else
543  {
544  TH1F hr("hr","recmodes",10,0,10);
545  c1->cd(21);
546 
547  for (int k=0;k<10;++k) hr.SetBinContent(k+1,(double) evcntrec[k].size()/evmult*norm);
548 
549  hr.DrawCopy();
550  }
551  }
552 
553  //t->Draw(">>els",sigcut+"&&"+cuts[myidx[0]]);
554 
555 
556  c1->Update();
557  }
558 
559  cout <<"\nCUT : "<<precut.Data()<<endl;
560  cout <<"SIG EVT: "<<Nsigev<<" ev "<<els.GetN()<<" cn BG: "<<Nbgev<<" ev "<<elb.GetN()<<" cn"<<endl;
561 
562  printf("SIG EFF : %6.1f%% BG EFF : %7.3f%%\n",Nsigev/N0_sig*100.,Nbgev/N0_bg*100.);
563  printf("SIG REL : %6.1f%% BG REL : %7.3f%%\n\n",Nsigev/nsig*100.,Nbgev/nbg*100.);
564 
565  cout<<"Recoil : "; for (int k=0;k<10;++k) printf(" %02d ",k);cout <<endl;
566  cout<<"Eff : "; for (int k=0;k<10;++k) printf("%6.1f%%",(double) evcntrec[k].size()/evmult*norm*100.);cout <<endl<<endl;
567 
568  t->SetEventList(0);
569 
570  if (isbt)
571  {
572  cout <<"Vars : ";
573  for (j=0;j<20;++j) cout <<vars[myidx[j]].Data()<<" ";
574  cout <<endl;
575  }
576  return 0;
577 }
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
TEventList * els
Definition: findcuts.C:52
void makeMaps(TTree *t, TString varname, TEventList &els, TEventList &elb, ValueMap &sig, ValueMap &bg, int id)
Definition: autocutx.C:153
Int_t run
Definition: autocutx.C:47
Int_t i
Definition: run_full.C:25
TTree * b
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
#define MAX
Definition: cutfinderx.C:23
TLorentzVector s
Definition: Pnd2DStar.C:50
int ev
bool dstarmode
Definition: autocutx.C:42
TString cuts[MAX]
Definition: autocutx.C:35
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:26
double N0_sig
Definition: autocutx.C:40
double Nsigev
Definition: autocutx.C:41
CountMap evcntrec[10]
Definition: autocutx.C:31
double minh[MAX]
Definition: cutfinderx.C:43
__m128 v
Definition: P4_F32vec4.h:4
int Ns
Definition: findcuts.C:56
bool mycompare(int i, int j)
Definition: autocutx.C:64
static int init
Definition: ranlxd.cxx:374
double bestSuppressionEvt(TTree *t, TString varname, TEventList &els, TEventList &elb, double &bestcut, double eff, int id)
Definition: autocutx.C:259
double cut[MAX]
Definition: autocutx.C:36
Int_t mode
Definition: autocutx.C:47
int idx[MAX]
Definition: autocutx.C:38
int uid(int lev, int lrun, int lmode)
Definition: autocutx.C:122
std::vector< pair< double, int > > ValueMap
Definition: autocutx.C:28
Float_t fbranch[MAX]
Definition: autocutx.C:44
double N0_bg
Definition: autocutx.C:40
int gettype(TTree *t, TString varname)
Definition: autocutx.C:51
float_m ok
double bestCombiEvt(TTree *t, TString varname, TEventList &els, TEventList &elb, double &bestcut, double eff, int id)
Definition: cutfinderx.C:325
int Nb
Definition: findcuts.C:56
TString vars[MAX]
Definition: autocutx.C:34
double maxh[MAX]
Definition: cutfinderx.C:43
TFile * f
Definition: bump_analys.C:12
Bool_t bbranch[MAX]
Definition: autocutx.C:46
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
#define BINS
Definition: autocutx.C:23
c1
Definition: plot_dirc.C:35
TGraph grl[MAX]
Definition: cutfinderx.C:40
std::map< TString, TString > mctvar
Definition: autocutx.C:32
int cutfinderx(TString fname, TString precut="", double supr=0.95, int evmult=10000, double norm=1.0, int n0s=-1)
Definition: cutfinderx.C:397
Int_t nbranch
Definition: autocutx.C:47
int nsig
Definition: toy_core.C:46
std::map< int, int > CountMap
Definition: autocutx.C:29
TEventList * elb
Definition: findcuts.C:52
Int_t cnt
Definition: hist-t7.C:106
int countEvents(TTree *t, TEventList &el)
Definition: autocutx.C:129
TTree * t
Definition: bump_analys.C:13
double Nbgev
Definition: autocutx.C:41
TGraph grr[MAX]
Definition: cutfinderx.C:41
Int_t ibranch[MAX]
Definition: autocutx.C:45
double signi[MAX]
Definition: cutfinderx.C:35
Int_t rec
Definition: autocutx.C:47
CountMap evcnt
Definition: autocutx.C:31
double bestEffEvt(TTree *t, TString varname, TEventList &els, TEventList &elb, double &bestcut, double supr, int id)
Definition: autocutx.C:189
Bool_t isbt
Definition: cutfinderx.C:33