Logo ROOT  
Reference Guide
mvasMulticlass.cxx
Go to the documentation of this file.
2#include "TMVA/Types.h"
3#include "TLegend.h"
4#include "TText.h"
5#include "TH2.h"
6
7
8
9// this macro plots the resulting MVA distributions (Signal and
10// Background overlayed) of different MVA methods run in TMVA
11// (e.g. running TMVAnalysis.C).
12
13
14// input: - Input file (result from TMVA)
15// - use of TMVA plotting TStyle
16void TMVA::mvasMulticlass(TString dataset, TString fin , HistType htype , Bool_t useTMVAStyle )
17{
18 // set style and remove existing canvas'
19 TMVAGlob::Initialize( useTMVAStyle );
20
21 // switches
22 const Bool_t Save_Images = kTRUE;
23
24 // checks if file with name "fin" is already open, and if not opens one
26
27 // find directory from input variable transformation and extract class names
28 TIter keys_dir = file->GetDirectory(dataset.Data())->GetListOfKeys();
29 TKey *key;
30 TDirectory *tempdir = 0;
31 while ((key = (TKey *)keys_dir())) {
32 TString name = key->GetName();
33 if (name.BeginsWith("InputVariables_")) {
34 tempdir = (TDirectory *)(file->GetDirectory(dataset.Data())->Get(name));
35 break;
36 }
37 }
38 std::vector<TString> classnames(TMVAGlob::GetClassNames(tempdir));
39
40 // define Canvas layout here!
41 // Int_t xPad = 1; // no of plots in x
42 // Int_t yPad = 1; // no of plots in y
43 // Int_t noPad = xPad * yPad ;
44 const Int_t width = 600; // size of canvas
45
46 // this defines how many canvases we need
47 TCanvas *c = 0;
48
49 // counter variables
50 Int_t countCanvas = 0;
51
52 // search for the right histograms in full list of keys
53 TIter next(file->GetDirectory(dataset.Data())->GetListOfKeys());
54 key = 0;
55 while ((key = (TKey*)next())) {
56
57 if (!TString(key->GetName()).BeginsWith("Method_")) continue;
58 if (!gROOT->GetClass(key->GetClassName())->InheritsFrom("TDirectory")) continue;
59
60 TString methodName;
61 TMVAGlob::GetMethodName(methodName,key);
62
63 TDirectory* mDir = (TDirectory*)key->ReadObj();
64
65 TIter keyIt(mDir->GetListOfKeys());
66 TKey *titkey;
67 while ((titkey = (TKey*)keyIt())) {
68
69 if (!gROOT->GetClass(titkey->GetClassName())->InheritsFrom("TDirectory")) continue;
70
71 TDirectory *titDir = (TDirectory *)titkey->ReadObj();
72 TString methodTitle;
73 TMVAGlob::GetMethodTitle(methodTitle,titDir);
74
75 cout << "--- Found directory for method: " << methodName << "::" << methodTitle << endl;
76 TString hname = "MVA_" + methodTitle;
77 for(UInt_t icls = 0; icls < classnames.size(); ++icls){
78 TObjArray hists;
79
80 std::vector<TString>::iterator classiter = classnames.begin();
81 for(; classiter!=classnames.end(); ++classiter){
82 TString name(hname+"_Test_"+ classnames.at(icls)
83 + "_prob_for_" + *classiter);
84 TH1 *hist = (TH1*)titDir->Get(name);
85 if (hist==0){
86 cout << ":\t mva distribution not available (this is normal for Cut classifier)" << endl;
87 continue;
88 }
89 hists.Add(hist);
90 }
91
92
93 // chop off useless stuff
94 ((TH1*)hists.First())->SetTitle( Form("TMVA response for classifier: %s", methodTitle.Data() ));
95
96 // create new canvas
97 //cout << "Create canvas..." << endl;
98 TString ctitle = ((htype == kMVAType) ?
99 Form("TMVA response for class %s %s", classnames.at(icls).Data(),methodTitle.Data()) :
100 Form("TMVA comparison for class %s %s", classnames.at(icls).Data(),methodTitle.Data())) ;
101
102 c = new TCanvas( Form("canvas%d", countCanvas+1), ctitle,
103 countCanvas*50+200, countCanvas*20, width, (Int_t)width*0.78 );
104
105 // set the histogram style
106 //cout << "Set histogram style..." << endl;
108
109 // normalise all histograms and find maximum
110 Float_t histmax = -1;
111 for(Int_t i=0; i<hists.GetEntriesFast(); ++i){
112 TMVAGlob::NormalizeHist((TH1*)hists[i] );
113 if(((TH1*)hists[i])->GetMaximum() > histmax)
114 histmax = ((TH1*)hists[i])->GetMaximum();
115 }
116
117 // frame limits (between 0 and 1 per definition)
118 Float_t xmin = 0;
119 Float_t xmax = 1;
120 Float_t ymin = 0;
121 Float_t maxMult = (htype == kCompareType) ? 1.3 : 1.2;
122 Float_t ymax = histmax*maxMult;
123 // build a frame
124 Int_t nb = 500;
125 TString hFrameName(TString("frame") + methodTitle);
126 TObject *o = gROOT->FindObject(hFrameName);
127 if(o) delete o;
128 TH2F* frame = new TH2F( hFrameName, ((TH1*)hists.First())->GetTitle(),
129 nb, xmin, xmax, nb, ymin, ymax );
130 frame->GetXaxis()->SetTitle( methodTitle + " response for "+classnames.at(icls));
131 frame->GetYaxis()->SetTitle("(1/N) dN^{ }/^{ }dx");
133
134 // eventually: draw the frame
135 frame->Draw();
136
137 c->GetPad(0)->SetLeftMargin( 0.105 );
138 frame->GetYaxis()->SetTitleOffset( 1.2 );
139
140 // Draw legend
141 TLegend *legend= new TLegend( c->GetLeftMargin(), 1 - c->GetTopMargin() - 0.12,
142 c->GetLeftMargin() + (htype == kCompareType ? 0.40 : 0.3), 1 - c->GetTopMargin() );
143 legend->SetFillStyle( 1 );
144 classiter = classnames.begin();
145
146 for(Int_t i=0; i<hists.GetEntriesFast(); ++i, ++classiter){
147 legend->AddEntry(((TH1*)hists[i]),*classiter,"F");
148 }
149
150 legend->SetBorderSize(1);
151 legend->SetMargin( 0.3 );
152 legend->Draw("same");
153
154
155 for(Int_t i=0; i<hists.GetEntriesFast(); ++i){
156
157 ((TH1*)hists[i])->Draw("histsame");
158 TString ytit = TString("(1/N) ") + ((TH1*)hists[i])->GetYaxis()->GetTitle();
159 ((TH1*)hists[i])->GetYaxis()->SetTitle( ytit ); // histograms are normalised
160
161 }
162
163
164 if (htype == kCompareType) {
165
166 TObjArray othists;
167 // if overtraining check, load additional histograms
168 classiter = classnames.begin();
169 for(; classiter!=classnames.end(); ++classiter){
170 TString name(hname+"_Train_"+ classnames.at(icls)
171 + "_prob_for_" + *classiter);
172 TH1 *hist = (TH1*)titDir->Get(name);
173 if (hist==0){
174 cout << ":\t comparison histogram for overtraining check not available!" << endl;
175 continue;
176 }
177 othists.Add(hist);
178 }
179
180 TLegend *legend2= new TLegend( 1 - c->GetRightMargin() - 0.42, 1 - c->GetTopMargin() - 0.12,
181 1 - c->GetRightMargin(), 1 - c->GetTopMargin() );
182 legend2->SetFillStyle( 1 );
183 legend2->SetBorderSize(1);
184
185 classiter = classnames.begin();
186 for(Int_t i=0; i<othists.GetEntriesFast(); ++i, ++classiter){
187 legend2->AddEntry(((TH1*)othists[i]),*classiter+" (training sample)","P");
188 }
189 legend2->SetMargin( 0.1 );
190 legend2->Draw("same");
191
192 // normalise all histograms and get maximum
193 for(Int_t i=0; i<othists.GetEntriesFast(); ++i){
194 TMVAGlob::NormalizeHist((TH1*)othists[i] );
195 if(((TH1*)othists[i])->GetMaximum() > histmax)
196 histmax = ((TH1*)othists[i])->GetMaximum();
197 }
198
200 for(Int_t i=0; i<othists.GetEntriesFast(); ++i){
201 Int_t col = ((TH1*)hists[i])->GetLineColor();
202 ((TH1*)othists[i])->SetMarkerSize( 0.7 );
203 ((TH1*)othists[i])->SetMarkerStyle( 20 );
204 ((TH1*)othists[i])->SetMarkerColor( col );
205 ((TH1*)othists[i])->SetLineWidth( 1 );
206 ((TH1*)othists[i])->Draw("e1same");
207 }
208
209 ymax = histmax*maxMult;
210 frame->GetYaxis()->SetLimits( 0, ymax );
211
212 // for better visibility, plot thinner lines
214 for(Int_t i=0; i<hists.GetEntriesFast(); ++i){
215 ((TH1*)hists[i])->SetLineWidth( 1 );
216 }
217
218
219 // perform K-S test
220
221 cout << "--- Perform Kolmogorov-Smirnov tests" << endl;
222 cout << "--- Goodness of consistency for class " << classnames.at(icls)<< endl;
223 //TString probatext("Kolmogorov-Smirnov test: ");
224 for(Int_t j=0; j<othists.GetEntriesFast(); ++j){
225 Float_t kol = ((TH1*)hists[j])->KolmogorovTest(((TH1*)othists[j]),"X");
226 cout << classnames.at(j) << ": " << kol << endl;
227 //probatext.Append(classnames.at(j)+Form(" %.3f ",kol));
228 }
229
230
231
232 //TText* tt = new TText( 0.12, 0.74, probatext );
233 //tt->SetNDC(); tt->SetTextSize( 0.032 ); tt->AppendPad();
234
235 }
236
237
238 // redraw axes
239 frame->Draw("sameaxis");
240
241 // text for overflows
242 //Int_t nbin = sig->GetNbinsX();
243 //Double_t dxu = sig->GetBinWidth(0);
244 //Double_t dxo = sig->GetBinWidth(nbin+1);
245 //TString uoflow = Form( "U/O-flow (S,B): (%.1f, %.1f)%% / (%.1f, %.1f)%%",
246 // sig->GetBinContent(0)*dxu*100, bgd->GetBinContent(0)*dxu*100,
247 // sig->GetBinContent(nbin+1)*dxo*100, bgd->GetBinContent(nbin+1)*dxo*100 );
248 //TText* t = new TText( 0.975, 0.115, uoflow );
249 //t->SetNDC();
250 //t->SetTextSize( 0.030 );
251 //t->SetTextAngle( 90 );
252 //t->AppendPad();
253
254 // update canvas
255 c->Update();
256
257 // save canvas to file
258
259 TMVAGlob::plot_logo(1.058);
260 if (Save_Images) {
261 if (htype == kMVAType) TMVAGlob::imgconv( c, Form("%s/plots/mva_%s_%s",dataset.Data(),classnames.at(icls).Data(), methodTitle.Data()) );
262 else if (htype == kCompareType) TMVAGlob::imgconv( c, Form("%s/plots/overtrain_%s_%s",dataset.Data(),classnames.at(icls).Data(), methodTitle.Data()) );
263
264 }
265 countCanvas++;
266 }
267 }
268 cout << "";
269 }
270}
271
#define c(i)
Definition: RSha256.hxx:101
int Int_t
Definition: RtypesCore.h:41
unsigned int UInt_t
Definition: RtypesCore.h:42
bool Bool_t
Definition: RtypesCore.h:59
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
include TDocParser_001 C image html pict1_TDocParser_001 png width
Definition: TDocParser.cxx:121
char name[80]
Definition: TGX11.cxx:109
float xmin
Definition: THbookFile.cxx:93
float ymin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
float ymax
Definition: THbookFile.cxx:93
#define gROOT
Definition: TROOT.h:415
char * Form(const char *fmt,...)
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:294
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition: TAttFill.h:39
virtual void SetLimits(Double_t xmin, Double_t xmax)
Definition: TAxis.h:154
The Canvas class.
Definition: TCanvas.h:31
TIter begin() const
Definition: TCollection.h:283
Describe directory structure in memory.
Definition: TDirectory.h:34
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
Definition: TDirectory.cxx:805
virtual TList * GetListOfKeys() const
Definition: TDirectory.h:160
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:48
The TH1 histogram class.
Definition: TH1.h:56
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:316
TAxis * GetYaxis()
Definition: TH1.h:317
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:251
Book space in a file, create I/O buffers, to fill them, (un)compress them.
Definition: TKey.h:24
virtual const char * GetClassName() const
Definition: TKey.h:72
virtual TObject * ReadObj()
To read a TObject* from the file.
Definition: TKey.cxx:729
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:23
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:330
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:423
void SetMargin(Float_t margin)
Definition: TLegend.h:69
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
An array of TObjects.
Definition: TObjArray.h:37
Int_t GetEntriesFast() const
Definition: TObjArray.h:64
void Add(TObject *obj)
Definition: TObjArray.h:74
TObject * First() const
Return the object in the first slot.
Definition: TObjArray.cxx:495
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void SetBorderSize(Int_t bordersize=4)
Definition: TPave.h:73
Basic string class.
Definition: TString.h:131
const char * Data() const
Definition: TString.h:364
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
Definition: TString.h:610
std::string GetMethodName(TCppMethod_t)
Definition: Cppyy.cxx:757
void Initialize(Bool_t useTMVAStyle=kTRUE)
Definition: tmvaglob.cxx:176
void GetMethodTitle(TString &name, TKey *ikey)
Definition: tmvaglob.cxx:341
void plot_logo(Float_t v_scale=1.0, Float_t skew=1.0)
Definition: tmvaglob.cxx:263
TFile * OpenFile(const TString &fin)
Definition: tmvaglob.cxx:192
void SetFrameStyle(TH1 *frame, Float_t scale=1.0)
Definition: tmvaglob.cxx:77
void SetMultiClassStyle(TObjArray *hists)
Definition: tmvaglob.cxx:47
void NormalizeHist(TH1 *h)
Definition: tmvaglob.cxx:308
std::vector< TString > GetClassNames(TDirectory *dir)
Definition: tmvaglob.cxx:462
void imgconv(TCanvas *c, const TString &fname)
Definition: tmvaglob.cxx:212
void mvasMulticlass(TString dataset, TString fin="TMVAMulticlass.root", HistType htype=kMVAType, Bool_t useTMVAStyle=kTRUE)
Definition: file.py:1