ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ROCCalc.cxx
Go to the documentation of this file.
1 #include <algorithm>
2 #include <cstdlib>
3 #include <errno.h>
4 
5 #include "TObjString.h"
6 #include "TMath.h"
7 #include "TString.h"
8 #include "TTree.h"
9 #include "TLeaf.h"
10 #include "TH1.h"
11 #include "TH2.h"
12 #include "TList.h"
13 #include "TSpline.h"
14 #include "TVector.h"
15 #include "TMatrixD.h"
16 #include "TMatrixDSymEigen.h"
17 #include "TVectorD.h"
18 #include "TTreeFormula.h"
19 #include "TXMLEngine.h"
20 #include "TROOT.h"
21 #include "TMatrixDSymEigen.h"
22 #include "TColor.h"
23 #include "TMVA/Config.h"
24 
25 
26 #ifndef ROOT_TMVA_Tools
27 #include "TMVA/Tools.h"
28 #endif
29 #ifndef ROOT_TMVA_ROCCalc
30 #include "TMVA/ROCCalc.h"
31 #endif
32 #ifndef ROOT_TMVA_Config
33 #include "TMVA/Config.h"
34 #endif
35 #ifndef ROOT_TMVA_Event
36 #include "TMVA/Event.h"
37 #endif
38 #ifndef ROOT_TMVA_Version
39 #include "TMVA/Version.h"
40 #endif
41 #ifndef ROOT_TMVA_PDF
42 #include "TMVA/PDF.h"
43 #endif
44 #ifndef ROOT_TMVA_MsgLogger
45 #include "TMVA/MsgLogger.h"
46 #endif
47 
48 #include "TMVA/PDF.h"
49 #include "TMVA/TSpline1.h"
50 #include "TMVA/TSpline2.h"
51 #include "TMVA/Types.h"
52 
53 using namespace std;
54 
55 ////////////////////////////////////////////////////////////////////////////////
56 
58  fMaxIter(100),
59  fAbsTol(0.0),
60  fmvaS(0),
61  fmvaB(0),
62  fmvaSpdf(0),
63  fmvaBpdf(0),
64  fSplS(0),
65  fSplB(0),
66  fSplmvaCumS(0),
67  fSplmvaCumB(0),
68  fSpleffBvsS(0),
69  fnStot(0),
70  fnBtot(0),
71  fSignificance(0),
72  fPurity(0),
73  fLogger ( new TMVA::MsgLogger("ROCCalc") )
74 {
76  fNbins = 100;
77  // fmvaS = (TH1*) mvaS->Clone("MVA Signal"); fmvaS->SetTitle("MVA Signal");
78  // fmvaB = (TH1*) mvaB->Clone("MVA Backgr"); fmvaB->SetTitle("MVA Backgr");
79  fmvaS = mvaS; fmvaS->SetTitle("MVA Signal");
80  fmvaB = mvaB; fmvaB->SetTitle("MVA Backgr");
81  fXmax = fmvaS->GetXaxis()->GetXmax();
82  fXmin = fmvaS->GetXaxis()->GetXmin();
83 
84  if (TMath::Abs(fXmax-fmvaB->GetXaxis()->GetXmax()) > 0.000001 ||
85  TMath::Abs(fXmin-fmvaB->GetXaxis()->GetXmin()) > 0.000001 ||
86  fmvaB->GetNbinsX() != fmvaS->GetNbinsX()) {
87  Log() << kFATAL << " Cannot cal ROC curve etc, as in put mvaS and mvaB have differen #nbins or range "<<Endl;
88  }
89  if (!strcmp(fmvaS->GetXaxis()->GetTitle(),"")) fmvaS->SetXTitle("MVA-value");
90  if (!strcmp(fmvaB->GetXaxis()->GetTitle(),"")) fmvaB->SetXTitle("MVA-value");
91  if (!strcmp(fmvaS->GetYaxis()->GetTitle(),"")) fmvaS->SetYTitle("#entries");
92  if (!strcmp(fmvaB->GetYaxis()->GetTitle(),"")) fmvaB->SetYTitle("#entries");
94  fmvaSpdf = mvaS->RebinX(mvaS->GetNbinsX()/100,"MVA Signal PDF");
95  fmvaBpdf = mvaB->RebinX(mvaB->GetNbinsX()/100,"MVA Backgr PDF");
96  fmvaSpdf->SetTitle("MVA Signal PDF");
97  fmvaBpdf->SetTitle("MVA Backgr PDF");
103 
104  fCutOrientation = (fmvaS->GetMean() > fmvaB->GetMean()) ? +1 : -1;
105 
106  fNevtS = 0;
107 
108 }
109 
110 ////////////////////////////////////////////////////////////////////////////////
111 /// Int_t c_Canvas = TColor::GetColor( "#f0f0f0" );
112 /// Int_t c_FrameFill = TColor::GetColor( "#fffffd" );
113 /// Int_t c_TitleBox = TColor::GetColor( "#5D6B7D" );
114 /// Int_t c_TitleBorder = TColor::GetColor( "#7D8B9D" );
115 /// Int_t c_TitleText = TColor::GetColor( "#FFFFFF" );
116 
118  Int_t c_SignalLine = TColor::GetColor( "#0000ee" );
119  Int_t c_SignalFill = TColor::GetColor( "#7d99d1" );
120  Int_t c_BackgroundLine = TColor::GetColor( "#ff0000" );
121  Int_t c_BackgroundFill = TColor::GetColor( "#ff0000" );
122  // Int_t c_NovelBlue = TColor::GetColor( "#2244a5" );
123 
124  //signal
125  // const Int_t FillColor__S = 38 + 150; // change of Color Scheme in ROOT-5.16.
126  // convince yourself with gROOT->GetListOfColors()->Print()
127  Int_t FillColor__S = c_SignalFill;
128  Int_t FillStyle__S = 1001;
129  Int_t LineColor__S = c_SignalLine;
130  Int_t LineWidth__S = 2;
131 
132  // background
133  //Int_t icolor = gConfig().fVariablePlotting.fUsePaperStyle ? 2 + 100 : 2;
134  Int_t FillColor__B = c_BackgroundFill;
135  Int_t FillStyle__B = 3554;
136  Int_t LineColor__B = c_BackgroundLine;
137  Int_t LineWidth__B = 2;
138 
139  if (sig != NULL) {
140  sig->SetLineColor( LineColor__S );
141  sig->SetLineWidth( LineWidth__S );
142  sig->SetFillStyle( FillStyle__S );
143  sig->SetFillColor( FillColor__S );
144  }
145 
146  if (bkg != NULL) {
147  bkg->SetLineColor( LineColor__B );
148  bkg->SetLineWidth( LineWidth__B );
149  bkg->SetFillStyle( FillStyle__B );
150  bkg->SetFillColor( FillColor__B );
151  }
152 
153  if (any != NULL) {
154  any->SetLineColor( LineColor__S );
155  any->SetLineWidth( LineWidth__S );
156  any->SetFillStyle( FillStyle__S );
157  any->SetFillColor( FillColor__S );
158  }
159 }
160 
161 
162 ////////////////////////////////////////////////////////////////////////////////
163 /// destructor
164 
166  // delete Splines and all histograms that were created only for internal use
167  if (fSplS) { delete fSplS; fSplS = 0; }
168  if (fSplB) { delete fSplB; fSplB = 0; }
169  if (fSpleffBvsS) { delete fSpleffBvsS; fSpleffBvsS = 0; }
170  if (fSplmvaCumS) { delete fSplmvaCumS; fSplmvaCumS = 0; }
171  if (fSplmvaCumB) { delete fSplmvaCumB; fSplmvaCumB = 0; }
172  if (fmvaScumul) { delete fmvaScumul; }
173  if (fmvaBcumul) { delete fmvaBcumul; }
174 
175  delete fLogger;
176 }
177 
178 ////////////////////////////////////////////////////////////////////////////////
179 /// get the ROC curve
180 
182  // first get the cumulative distributions of the mva distribution
183  // --> efficiencies vs cut value
184  fNevtS = fmvaS->GetSumOfWeights(); // needed to get the error on the eff.. will only be correct if the histogram is not scaled to "integral == 1" Yet;
185  if (fNevtS < 2) {
186  Log() << kWARNING << "I guess the mva distributions fed into ROCCalc were already normalized, therefore the calculated error on the efficiency will be incorrect !! " << Endl;
187  fNevtS = 0; // reset to zero --> no error will be calculated on the efficiencies
188  }
189  fmvaScumul = gTools().GetCumulativeDist(fmvaS);
190  fmvaBcumul = gTools().GetCumulativeDist(fmvaB);
191  fmvaScumul->Scale( 1.0/TMath::Max(std::numeric_limits<double>::epsilon(),fmvaScumul->GetMaximum()) );
192  fmvaBcumul->Scale( 1.0/TMath::Max(std::numeric_limits<double>::epsilon(),fmvaBcumul->GetMaximum()) );
193  fmvaScumul->SetMinimum(0);
194  fmvaBcumul->SetMinimum(0);
195  // fmvaScumul->Draw("hist");
196  // fmvaBcumul->Draw("histsame");
197 
198  // background efficiency versus signal efficiency
199  TH1D* effBvsS = new TH1D("effBvsS", "ROC-Curve", fNbins, 0, 1 );
200  effBvsS->SetXTitle( "Signal eff" );
201  effBvsS->SetYTitle( "Backgr eff" );
202 
203  // background rejection (=1-eff.) versus signal efficiency
204  TH1D* rejBvsS = new TH1D( "rejBvsS", "ROC-Curve", fNbins, 0, 1 );
205  rejBvsS->SetXTitle( "Signal eff" );
206  rejBvsS->SetYTitle( "Backgr rejection (1-eff)" );
207 
208  // inverse background eff (1/eff.) versus signal efficiency
209  TH1D* inveffBvsS = new TH1D("invBeffvsSeff", "ROC-Curve" , fNbins, 0, 1 );
210  inveffBvsS->SetXTitle( "Signal eff" );
211  inveffBvsS->SetYTitle( "Inverse backgr. eff (1/eff)" );
212 
213  // use root finder
214  // spline background efficiency plot
215  // note that there is a bin shift when going from a TH1D object to a TGraph :-(
216  if (fUseSplines) {
217  fSplmvaCumS = new TSpline1( "spline2_signal", new TGraph( fmvaScumul ) );
218  fSplmvaCumB = new TSpline1( "spline2_background", new TGraph( fmvaBcumul ) );
219  // verify spline sanity
220  gTools().CheckSplines( fmvaScumul, fSplmvaCumS );
221  gTools().CheckSplines( fmvaBcumul, fSplmvaCumB );
222  }
223 
224  Double_t effB = 0;
225  for (UInt_t bini=1; bini<=fNbins; bini++) {
226 
227  // find cut value corresponding to a given signal efficiency
228  Double_t effS = effBvsS->GetBinCenter( bini );
229  Double_t cut = Root( effS );
230 
231  // retrieve background efficiency for given cut
232  if (fUseSplines) effB = fSplmvaCumB->Eval( cut );
233  else effB = fmvaBcumul->GetBinContent( fmvaBcumul->FindBin( cut ) );
234 
235  // and fill histograms
236  effBvsS->SetBinContent( bini, effB );
237  rejBvsS->SetBinContent( bini, 1.0-effB );
239  inveffBvsS->SetBinContent( bini, 1.0/effB );
240  }
241 
242  // create splines for histogram
243  fSpleffBvsS = new TSpline1( "effBvsS", new TGraph( effBvsS ) );
244 
245  // search for overlap point where, when cutting on it,
246  // one would obtain: eff_S = rej_B = 1 - eff_B
247 
248  Double_t effS = 0., rejB = 0., effS_ = 0., rejB_ = 0.;
249  Int_t nbins = 5000;
250  for (Int_t bini=1; bini<=nbins; bini++) {
251 
252  // get corresponding signal and background efficiencies
253  effS = (bini - 0.5)/Float_t(nbins);
254  rejB = 1.0 - fSpleffBvsS->Eval( effS );
255 
256  // find signal efficiency that corresponds to required background efficiency
257  if ((effS - rejB)*(effS_ - rejB_) < 0) break;
258  effS_ = effS;
259  rejB_ = rejB;
260  }
261  // find cut that corresponds to signal efficiency and update signal-like criterion
262  fSignalCut = Root( 0.5*(effS + effS_) );
263 
264  return rejBvsS;
265 }
266 
267 ////////////////////////////////////////////////////////////////////////////////
268 /// code to compute the area under the ROC ( rej-vs-eff ) curve
269 
271  Double_t effS = 0, effB = 0;
272  Int_t nbins = 1000;
273  if (fSpleffBvsS == 0) this->GetROC(); // that will make the ROC calculation if not done yet
274 
275  // compute area of rej-vs-eff plot
276  Double_t integral = 0;
277  for (Int_t bini=1; bini<=nbins; bini++) {
278 
279  // get corresponding signal and background efficiencies
280  effS = (bini - 0.5)/Float_t(nbins);
281  effB = fSpleffBvsS->Eval( effS );
282  integral += (1.0 - effB);
283  }
284  integral /= nbins;
285 
286  return integral;
287 }
288 
289 ////////////////////////////////////////////////////////////////////////////////
290 /// get the signal efficiency for a particular backgroud efficiency
291 /// that will be the value of the efficiency retured (does not affect
292 /// the efficiency-vs-bkg plot which is done anyway.
293 
295  // find precise efficiency value
296  Double_t effS=0., effB, effSOld=1., effBOld=0.;
297  Int_t nbins = 1000;
298  if (fSpleffBvsS == 0) this->GetROC(); // that will make the ROC calculation if not done yet
299 
300  Float_t step=1./nbins; // stepsize in efficiency binning
301  for (Int_t bini=1; bini<=nbins; bini++) {
302  // get corresponding signal and background efficiencies
303  effS = (bini - 0.5)*step; // efficiency goes from 0-to-1 in nbins steps of 1/nbins (take middle of the bin)
304  effB = fSpleffBvsS->Eval( effS );
305 
306  // find signal efficiency that corresponds to required background efficiency
307  if ((effB - effBref)*(effBOld - effBref) <= 0) break;
308  effSOld = effS;
309  effBOld = effB;
310  }
311 
312  // take mean between bin above and bin below
313  effS = 0.5*(effS + effSOld);
314 
315 
316  if (fNevtS > 0) effSerr = TMath::Sqrt( effS*(1.0 - effS)/fNevtS );
317  else effSerr = 0;
318 
319  return effS;
320 }
321 
322 ////////////////////////////////////////////////////////////////////////////////
323 /// returns efficiency as function of cut
324 
326 {
327  Double_t retVal=0;
328 
329  // retrieve the class object
330  if (fUseSplines) retVal = fSplmvaCumS->Eval( theCut );
331  else retVal = fmvaScumul->GetBinContent( fmvaScumul->FindBin( theCut ) );
332 
333  // caution: here we take some "forbidden" action to hide a problem:
334  // in some cases, in particular for likelihood, the binned efficiency distributions
335  // do not equal 1, at xmin, and 0 at xmax; of course, in principle we have the
336  // unbinned information available in the trees, but the unbinned minimization is
337  // too slow, and we don't need to do a precision measurement here. Hence, we force
338  // this property.
339  Double_t eps = 1.0e-5;
340  if (theCut-fXmin < eps) retVal = (fCutOrientation > 0) ? 1.0 : 0.0;
341  else if (fXmax-theCut < eps) retVal = (fCutOrientation > 0) ? 0.0 : 1.0;
342 
343 
344  return retVal;
345 }
346 
347 ////////////////////////////////////////////////////////////////////////////////
348 /// Root finding using Brents algorithm; taken from CERNLIB function RZERO
349 
351 {
352  Double_t a = fXmin, b = fXmax;
353  Double_t fa = GetEffForRoot( a ) - refValue;
354  Double_t fb = GetEffForRoot( b ) - refValue;
355  if (fb*fa > 0) {
356  Log() << kWARNING << "<ROCCalc::Root> initial interval w/o root: "
357  << "(a=" << a << ", b=" << b << "),"
358  << " (Eff_a=" << GetEffForRoot( a )
359  << ", Eff_b=" << GetEffForRoot( b ) << "), "
360  << "(fa=" << fa << ", fb=" << fb << "), "
361  << "refValue = " << refValue << Endl;
362  return 1;
363  }
364 
365  Bool_t ac_equal(kFALSE);
366  Double_t fc = fb;
367  Double_t c = 0, d = 0, e = 0;
368  for (Int_t iter= 0; iter <= fMaxIter; iter++) {
369  if ((fb < 0 && fc < 0) || (fb > 0 && fc > 0)) {
370 
371  // Rename a,b,c and adjust bounding interval d
372  ac_equal = kTRUE;
373  c = a; fc = fa;
374  d = b - a; e = b - a;
375  }
376 
377  if (TMath::Abs(fc) < TMath::Abs(fb)) {
378  ac_equal = kTRUE;
379  a = b; b = c; c = a;
380  fa = fb; fb = fc; fc = fa;
381  }
382 
383  Double_t tol = 0.5 * 2.2204460492503131e-16 * TMath::Abs(b);
384  Double_t m = 0.5 * (c - b);
385  if (fb == 0 || TMath::Abs(m) <= tol || TMath::Abs(fb) < fAbsTol) return b;
386 
387  // Bounds decreasing too slowly: use bisection
388  if (TMath::Abs (e) < tol || TMath::Abs (fa) <= TMath::Abs (fb)) { d = m; e = m; }
389  else {
390  // Attempt inverse cubic interpolation
391  Double_t p, q, r;
392  Double_t s = fb / fa;
393 
394  if (ac_equal) { p = 2 * m * s; q = 1 - s; }
395  else {
396  q = fa / fc; r = fb / fc;
397  p = s * (2 * m * q * (q - r) - (b - a) * (r - 1));
398  q = (q - 1) * (r - 1) * (s - 1);
399  }
400  // Check whether we are in bounds
401  if (p > 0) q = -q;
402  else p = -p;
403 
404  Double_t min1 = 3 * m * q - TMath::Abs (tol * q);
405  Double_t min2 = TMath::Abs (e * q);
406  if (2 * p < (min1 < min2 ? min1 : min2)) {
407  // Accept the interpolation
408  e = d; d = p / q;
409  }
410  else { d = m; e = m; } // Interpolation failed: use bisection.
411  }
412  // Move last best guess to a
413  a = b; fa = fb;
414  // Evaluate new trial root
415  if (TMath::Abs(d) > tol) b += d;
416  else b += (m > 0 ? +tol : -tol);
417 
418  fb = GetEffForRoot( b ) - refValue;
419 
420  }
421 
422  // Return our best guess if we run out of iterations
423  Log() << kWARNING << "<ROCCalc::Root> maximum iterations (" << fMaxIter
424  << ") reached before convergence" << Endl;
425 
426  return b;
427 }
428 
429 ////////////////////////////////////////////////////////////////////////////////
430 
432 {
433  if (fnStot!=nStot || fnBtot!=nBtot || !fSignificance) {
434  GetSignificance(nStot, nBtot);
435  fnStot=nStot;
436  fnBtot=nBtot;
437  }
438  return fPurity;
439 }
440 ////////////////////////////////////////////////////////////////////////////////
441 
443 {
444  if (fnStot==nStot && fnBtot==nBtot && !fSignificance) return fSignificance;
445  fnStot=nStot; fnBtot=nBtot;
446 
447  fSignificance = (TH1*) fmvaScumul->Clone("Significance"); fSignificance->SetTitle("Significance");
448  fSignificance->Reset(); fSignificance->SetFillStyle(0);
449  fSignificance->SetXTitle("mva cut value");
450  fSignificance->SetYTitle("Stat. significance S/Sqrt(S+B)");
451  fSignificance->SetLineColor(2);
452  fSignificance->SetLineWidth(5);
453 
454  fPurity = (TH1*) fmvaScumul->Clone("Purity"); fPurity->SetTitle("Purity");
455  fPurity->Reset(); fPurity->SetFillStyle(0);
456  fPurity->SetXTitle("mva cut value");
457  fPurity->SetYTitle("Purity: S/(S+B)");
458  fPurity->SetLineColor(3);
459  fPurity->SetLineWidth(5);
460 
461  Double_t maxSig=0;
462  for (Int_t i=1; i<=fSignificance->GetNbinsX(); i++) {
463  Double_t S = fmvaScumul->GetBinContent( i ) * nStot;
464  Double_t B = fmvaBcumul->GetBinContent( i ) * nBtot;
465  Double_t purity;
466  Double_t sig;
467  if (S+B > 0){
468  purity = S/(S+B);
469  sig = S/TMath::Sqrt(S+B);
470  if (sig > maxSig) {
471  maxSig = sig;
472  }
473  } else {
474  purity=0;
475  sig=0;
476  }
477  cout << "S="<<S<<" B="<<B<< " purity="<<purity<< endl;
478  fPurity->SetBinContent( i, purity );
479  fSignificance->SetBinContent( i, sig );
480  }
481 
482  /*
483  TLatex* line1;
484  TLatex* line2;
485  TLatex tl;
486  tl.SetNDC();
487  tl.SetTextSize( 0.033 );
488  Int_t maxbin = fSignificance->GetMaximumBin();
489  line1 = tl.DrawLatex( 0.15, 0.23, Form("For %1.0f signal and %1.0f background", nStot, nBtot));
490  tl.DrawLatex( 0.15, 0.19, "events the maximum S/Sqrt(S+B) is");
491 
492  line2 = tl.DrawLatex( 0.15, 0.15, Form("%4.2f when cutting at %5.2f",
493  maxSig,
494  fSignificance->GetXaxis()->GetBinCenter(maxbin)) );
495  */
496  return fSignificance;
497 
498 }
499 
500 
501 
502 
503 
504 
505 
506 
virtual void SetLineWidth(Width_t lwidth)
Definition: TAttLine.h:57
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6174
static double B[]
MsgLogger & Log() const
message logger
Definition: ROCCalc.h:75
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
virtual void SetMaximum(Double_t maximum=-1111)
Definition: TH1.h:394
float Float_t
Definition: RtypesCore.h:53
return c
TH1 * fmvaSpdf
Definition: ROCCalc.h:60
ROCCalc(TH1 *mvaS, TH1 *mvaB)
Definition: ROCCalc.cxx:57
static Int_t c_SignalFill
Definition: tmvaglob.h:50
TH1 * fmvaBpdf
Definition: ROCCalc.h:60
UInt_t fNbins
Definition: ROCCalc.h:56
int Int_t
Definition: RtypesCore.h:41
virtual void SetYTitle(const char *title)
Definition: TH1.h:409
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual void SetFillStyle(Style_t fstyle)
Definition: TAttFill.h:52
int nbins[3]
virtual Int_t GetNbinsX() const
Definition: TH1.h:296
TH1 * fmvaS
Definition: ROCCalc.h:59
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
Definition: TH1.cxx:2565
TH1D * GetROC()
get the ROC curve
Definition: ROCCalc.cxx:181
static struct mg_connection * fc(struct mg_context *ctx)
Definition: civetweb.c:839
Tools & gTools()
Definition: Tools.cxx:79
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TObject.cxx:203
TH1 * fmvaB
Definition: ROCCalc.h:59
Double_t GetEffSForEffBof(Double_t effBref, Double_t &effSerr)
get the signal efficiency for a particular backgroud efficiency that will be the value of the efficie...
Definition: ROCCalc.cxx:294
Double_t Root(Double_t)
Root finding using Brents algorithm; taken from CERNLIB function RZERO.
Definition: ROCCalc.cxx:350
Bool_t fUseSplines
Definition: ROCCalc.h:57
int d
Definition: tornado.py:11
Double_t GetEffForRoot(Double_t theCut)
returns efficiency as function of cut
Definition: ROCCalc.cxx:325
Bool_t CheckSplines(const TH1 *, const TSpline *)
Definition: Tools.cxx:487
std::map< std::string, std::string >::const_iterator iter
Definition: TAlienJob.cxx:54
Float_t fXmin
Definition: ROCCalc.h:61
Double_t GetXmin() const
Definition: TAxis.h:137
TH1 * GetPurity(Int_t nStot, Int_t nBtot)
Definition: ROCCalc.cxx:431
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
virtual Double_t GetBinCenter(Int_t bin) const
return bin center for 1D historam Better to use h1.GetXaxis().GetBinCenter(bin)
Definition: TH1.cxx:8470
const double tol
virtual Double_t GetMean(Int_t axis=1) const
For axis = 1,2 or 3 returns the mean value of the histogram along X,Y or Z axis.
Definition: TH1.cxx:7014
ROOT::R::TRInterface & r
Definition: Object.C:4
Float_t fXmax
Definition: ROCCalc.h:61
virtual void SetFillColor(Color_t fcolor)
Definition: TAttFill.h:50
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition: TH1.cxx:8543
TH1 * GetSignificance(Int_t nStot, Int_t nBtot)
Definition: ROCCalc.cxx:442
unsigned int UInt_t
Definition: RtypesCore.h:42
TMarker * m
Definition: textangle.C:8
virtual Double_t GetSumOfWeights() const
Return the sum of weights excluding under/overflows.
Definition: TH1.cxx:7354
static Int_t GetColor(const char *hexcolor)
Static method returning color number for color specified by hex color string of form: "#rrggbb"...
Definition: TColor.cxx:1023
TAxis * GetYaxis()
Definition: TH1.h:320
const char * GetTitle() const
Returns title of object.
Definition: TAxis.h:133
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
TH1 * GetCumulativeDist(TH1 *h)
Definition: Tools.cxx:1759
REAL epsilon
Definition: triangle.c:617
const char * Root
Definition: TXMLSetup.cxx:47
static Int_t c_SignalLine
Definition: tmvaglob.h:49
virtual TH1 * RebinX(Int_t ngroup=2, const char *newname="")
Definition: TH1.h:347
double Double_t
Definition: RtypesCore.h:55
static const float S
Definition: mandel.cpp:113
Double_t GetXmax() const
Definition: TAxis.h:138
Int_t fCutOrientation
Definition: ROCCalc.h:63
The TH1 histogram class.
Definition: TH1.h:80
static Int_t c_BackgroundFill
Definition: tmvaglob.h:52
void ApplySignalAndBackgroundStyle(TH1 *sig, TH1 *bkg, TH1 *any=0)
Int_t c_Canvas = TColor::GetColor( "#f0f0f0" ); Int_t c_FrameFill = TColor::GetColor( "#fffffd" ); In...
Definition: ROCCalc.cxx:117
virtual void SetXTitle(const char *title)
Definition: TH1.h:408
Double_t fNevtS
Definition: ROCCalc.h:62
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
#define NULL
Definition: Rtypes.h:82
virtual void SetTitle(const char *title)
Change (i.e.
Definition: TH1.cxx:6268
Double_t GetROCIntegral()
code to compute the area under the ROC ( rej-vs-eff ) curve
Definition: ROCCalc.cxx:270
static Int_t c_BackgroundLine
Definition: tmvaglob.h:51
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
virtual Double_t GetMaximum(Double_t maxval=FLT_MAX) const
Return maximum value smaller than maxval of bins in the range, unless the value has been overridden b...
Definition: TH1.cxx:7921
const Bool_t kTRUE
Definition: Rtypes.h:91
float * q
Definition: THbookFile.cxx:87
~ROCCalc()
destructor
Definition: ROCCalc.cxx:165
Definition: math.cpp:60
TAxis * GetXaxis()
Definition: TH1.h:319