Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 <tuple>
29#include <utility>
30#include <vector>
31
32class TList;
33class TTree;
34class TString;
35class TH1;
36class TH2;
37class TH2F;
38class TSpline;
39class TSpline1;
40class TGraph;
41
42namespace TMVA {
43
44class MsgLogger;
45
46class ROCCurve {
47
48public:
49 ROCCurve(const std::vector<std::tuple<Float_t, Float_t, Bool_t>> &mvas);
50
51 ROCCurve(const std::vector<Float_t> &mvaValues, const std::vector<Bool_t> &mvaTargets,
52 const std::vector<Float_t> &mvaWeights);
53
54 ROCCurve(const std::vector<Float_t> &mvaValues, const std::vector<Bool_t> &mvaTargets);
55
56 ROCCurve(const std::vector<Float_t> &mvaSignal, const std::vector<Float_t> &mvaBackground,
57 const std::vector<Float_t> &mvaSignalWeights, const std::vector<Float_t> &mvaBackgroundWeights);
58
59 ROCCurve(const std::vector<Float_t> &mvaSignal, const std::vector<Float_t> &mvaBackground);
60
61 ~ROCCurve();
62
63 Double_t GetEffSForEffB(Double_t effB, const UInt_t num_points = 41);
64
66 TGraph *GetROCCurve(const UInt_t points = 100); // n divisions = #points -1
67
68 const std::vector<std::tuple<Float_t, Float_t, Bool_t>> GetMvas() const { return fMva; }
69private:
70 mutable MsgLogger *fLogger; ///<! message logger
71 MsgLogger &Log() const;
72
74
75 std::vector<std::tuple<Float_t, Float_t, Bool_t>> fMva;
76
77 std::vector<Double_t> ComputeSensitivity(const UInt_t num_points);
78 std::vector<Double_t> ComputeSpecificity(const UInt_t num_points);
79};
80}
81#endif
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:307
Service class for 2-D histogram classes.
Definition TH2.h:30
A doubly linked list.
Definition TList.h:38
ostringstream derivative to redirect and format output
Definition MsgLogger.h:57
std::vector< Double_t > ComputeSpecificity(const UInt_t num_points)
Definition ROCCurve.cxx:138
~ROCCurve()
destructor
Definition ROCCurve.cxx:123
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:217
TGraph * fGraph
Definition ROCCurve.h:73
std::vector< Double_t > ComputeSensitivity(const UInt_t num_points)
Definition ROCCurve.cxx:174
Double_t GetROCIntegral(const UInt_t points=41)
Calculates the ROC integral (AUC)
Definition ROCCurve.cxx:248
const std::vector< std::tuple< Float_t, Float_t, Bool_t > > GetMvas() const
Definition ROCCurve.h:68
MsgLogger & Log() const
Definition ROCCurve.cxx:128
std::vector< std::tuple< Float_t, Float_t, Bool_t > > fMva
Definition ROCCurve.h:75
MsgLogger * fLogger
! message logger
Definition ROCCurve.h:70
TGraph * GetROCCurve(const UInt_t points=100)
Returns a new TGraph containing the ROC curve.
Definition ROCCurve.cxx:274
Base class for spline implementation containing the Draw/Paint methods.
Definition TSpline.h:31
Basic string class.
Definition TString.h:139
A TTree represents a columnar dataset.
Definition TTree.h:79
create variable transformations
void mvas(TString dataset, TString fin="TMVA.root", HistType htype=kMVAType, Bool_t useTMVAStyle=kTRUE)