Logo ROOT  
Reference Guide
ROCCurve.h
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Omar Zapata, Lorenzo Moneta, Sergei Gleyzer, Kim Albertsson
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : ROCCurve *
8  * *
9  * Description: *
10  * This is class to compute ROC Integral (AUC) *
11  * *
12  * Authors : *
13  * Omar Zapata <Omar.Zapata@cern.ch> - UdeA/ITM Colombia *
14  * Lorenzo Moneta <Lorenzo.Moneta@cern.ch> - CERN, Switzerland *
15  * Sergei Gleyzer <Sergei.Gleyzer@cern.ch> - U of Florida & CERN *
16  * Kim Albertsson <kim.albertsson@cern.ch> - LTU & CERN *
17  * *
18  * Copyright (c) 2015: *
19  * CERN, Switzerland *
20  * UdeA/ITM, Colombia *
21  * U. of Florida, USA *
22  **********************************************************************************/
23 #ifndef ROOT_TMVA_ROCCurve
24 #define ROOT_TMVA_ROCCurve
25 
26 #include "RtypesCore.h"
27 
28 #include <vector>
29 #include <utility>
30 
31 class TList;
32 class TTree;
33 class TString;
34 class TH1;
35 class TH2;
36 class TH2F;
37 class TSpline;
38 class TSpline1;
39 class TGraph;
40 
41 namespace TMVA {
42 
43 class MsgLogger;
44 
45 class ROCCurve {
46 
47 public:
48  ROCCurve(const std::vector<std::tuple<Float_t, Float_t, Bool_t>> &mvas);
49 
50  ROCCurve(const std::vector<Float_t> &mvaValues, const std::vector<Bool_t> &mvaTargets,
51  const std::vector<Float_t> &mvaWeights);
52 
53  ROCCurve(const std::vector<Float_t> &mvaValues, const std::vector<Bool_t> &mvaTargets);
54 
55  ROCCurve(const std::vector<Float_t> &mvaSignal, const std::vector<Float_t> &mvaBackground,
56  const std::vector<Float_t> &mvaSignalWeights, const std::vector<Float_t> &mvaBackgroundWeights);
57 
58  ROCCurve(const std::vector<Float_t> &mvaSignal, const std::vector<Float_t> &mvaBackground);
59 
60  ~ROCCurve();
61 
62  Double_t GetEffSForEffB(Double_t effB, const UInt_t num_points = 41);
63 
65  TGraph *GetROCCurve(const UInt_t points = 100); // n divisions = #points -1
66 
67  const std::vector<std::tuple<Float_t, Float_t, Bool_t>> GetMvas() const { return fMva; }
68 private:
69  mutable MsgLogger *fLogger; //! message logger
70  MsgLogger &Log() const;
71 
73 
74  std::vector<std::tuple<Float_t, Float_t, Bool_t>> fMva;
75 
76  std::vector<Double_t> ComputeSensitivity(const UInt_t num_points);
77  std::vector<Double_t> ComputeSpecificity(const UInt_t num_points);
78 };
79 }
80 #endif
TH2F
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:251
TMVA::ROCCurve::fLogger
MsgLogger * fLogger
Definition: ROCCurve.h:69
TTree
A TTree represents a columnar dataset.
Definition: TTree.h:79
TMVA::ROCCurve::GetROCIntegral
Double_t GetROCIntegral(const UInt_t points=41)
Calculates the ROC integral (AUC)
Definition: ROCCurve.cxx:250
TSpline
Base class for spline implementation containing the Draw/Paint methods.
Definition: TSpline.h:31
TMVA::ROCCurve::fGraph
TGraph * fGraph
Definition: ROCCurve.h:72
TString
Basic string class.
Definition: TString.h:136
TMVA::ROCCurve::~ROCCurve
~ROCCurve()
destructor
Definition: ROCCurve.cxx:125
TMVA::ROCCurve::ComputeSpecificity
std::vector< Double_t > ComputeSpecificity(const UInt_t num_points)
Definition: ROCCurve.cxx:140
TMVA::ROCCurve::ROCCurve
ROCCurve(const std::vector< std::tuple< Float_t, Float_t, Bool_t >> &mvas)
Definition: ROCCurve.cxx:47
TH2
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
TMVA::ROCCurve::ComputeSensitivity
std::vector< Double_t > ComputeSensitivity(const UInt_t num_points)
Definition: ROCCurve.cxx:176
unsigned int
TMVA::ROCCurve::Log
MsgLogger & Log() const
message logger
Definition: ROCCurve.cxx:130
Double_t
double Double_t
Definition: RtypesCore.h:59
TGraph
A TGraph is an object made of two arrays X and Y with npoints each.
Definition: TGraph.h:41
RtypesCore.h
TMVA::MsgLogger
ostringstream derivative to redirect and format output
Definition: MsgLogger.h:59
TMVA::ROCCurve::GetROCCurve
TGraph * GetROCCurve(const UInt_t points=100)
Returns a new TGraph containing the ROC curve.
Definition: ROCCurve.cxx:276
TMVA::ROCCurve::GetMvas
const std::vector< std::tuple< Float_t, Float_t, Bool_t > > GetMvas() const
Definition: ROCCurve.h:67
points
point * points
Definition: X3DBuffer.c:22
TH1
TH1 is the base class of all histogram classes in ROOT.
Definition: TH1.h:58
TMVA::mvas
void mvas(TString dataset, TString fin="TMVA.root", HistType htype=kMVAType, Bool_t useTMVAStyle=kTRUE)
TMVA::ROCCurve
Definition: ROCCurve.h:45
TMVA::ROCCurve::fMva
std::vector< std::tuple< Float_t, Float_t, Bool_t > > fMva
Definition: ROCCurve.h:74
TList
A doubly linked list.
Definition: TList.h:44
TMVA
create variable transformations
Definition: GeneticMinimizer.h:22
TMVA::ROCCurve::GetEffSForEffB
Double_t GetEffSForEffB(Double_t effB, const UInt_t num_points=41)
Calculate the signal efficiency (sensitivity) for a given background efficiency (sensitivity).
Definition: ROCCurve.cxx:219