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 "Rtypes.h"
27
28#include <iomanip>
29#include <iostream>
30#include <sstream>
31#include <vector>
32
33class TList;
34class TTree;
35class TString;
36class TH1;
37class TH2;
38class TH2F;
39class TSpline;
40class TSpline1;
41class TGraph;
42
43namespace TMVA {
44
45class MsgLogger;
46
47class ROCCurve {
48
49public:
50 ROCCurve(const std::vector<std::tuple<Float_t, Float_t, Bool_t>> &mvas);
51
52 ROCCurve(const std::vector<Float_t> &mvaValues, const std::vector<Bool_t> &mvaTargets,
53 const std::vector<Float_t> &mvaWeights);
54
55 ROCCurve(const std::vector<Float_t> &mvaValues, const std::vector<Bool_t> &mvaTargets);
56
57 ROCCurve(const std::vector<Float_t> &mvaSignal, const std::vector<Float_t> &mvaBackground,
58 const std::vector<Float_t> &mvaSignalWeights, const std::vector<Float_t> &mvaBackgroundWeights);
59
60 ROCCurve(const std::vector<Float_t> &mvaSignal, const std::vector<Float_t> &mvaBackground);
61
62 ~ROCCurve();
63
64 Double_t GetEffSForEffB(Double_t effB, const UInt_t num_points = 41);
65
67 TGraph *GetROCCurve(const UInt_t points = 100); // n divisions = #points -1
68
69 const std::vector<std::tuple<Float_t, Float_t, Bool_t>> GetMvas() const { return fMva; }
70private:
71 mutable MsgLogger *fLogger; //! message logger
72 MsgLogger &Log() const;
73
75
76 std::vector<std::tuple<Float_t, Float_t, Bool_t>> fMva;
77
78 std::vector<Double_t> ComputeSensitivity(const UInt_t num_points);
79 std::vector<Double_t> ComputeSpecificity(const UInt_t num_points);
80};
81}
82#endif
unsigned int UInt_t
Definition: RtypesCore.h:42
double Double_t
Definition: RtypesCore.h:55
point * points
Definition: X3DBuffer.c:22
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:41
The TH1 histogram class.
Definition: TH1.h:56
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:251
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
A doubly linked list.
Definition: TList.h:44
ostringstream derivative to redirect and format output
Definition: MsgLogger.h:59
std::vector< Double_t > ComputeSpecificity(const UInt_t num_points)
Definition: ROCCurve.cxx:141
ROCCurve(const std::vector< std::tuple< Float_t, Float_t, Bool_t > > &mvas)
Definition: ROCCurve.cxx:48
~ROCCurve()
destructor
Definition: ROCCurve.cxx:126
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:220
TGraph * fGraph
Definition: ROCCurve.h:74
std::vector< Double_t > ComputeSensitivity(const UInt_t num_points)
Definition: ROCCurve.cxx:177
Double_t GetROCIntegral(const UInt_t points=41)
Calculates the ROC integral (AUC)
Definition: ROCCurve.cxx:251
const std::vector< std::tuple< Float_t, Float_t, Bool_t > > GetMvas() const
Definition: ROCCurve.h:69
MsgLogger & Log() const
message logger
Definition: ROCCurve.cxx:131
std::vector< std::tuple< Float_t, Float_t, Bool_t > > fMva
Definition: ROCCurve.h:76
MsgLogger * fLogger
Definition: ROCCurve.h:71
TGraph * GetROCCurve(const UInt_t points=100)
Returns a new TGraph containing the ROC curve.
Definition: ROCCurve.cxx:277
Base class for spline implementation containing the Draw/Paint methods.
Definition: TSpline.h:22
Basic string class.
Definition: TString.h:131
A TTree represents a columnar dataset.
Definition: TTree.h:72
create variable transformations
void mvas(TString dataset, TString fin="TMVA.root", HistType htype=kMVAType, Bool_t useTMVAStyle=kTRUE)