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