Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
BoostControlPlots.cxx
Go to the documentation of this file.
2 #include <vector>
3 #include <string>
4 #include "TLegend.h"
5 #include "TText.h"
6 
7 
8 
9 // input: - Input file (result from TMVA),
10 // - use of TMVA plotting TStyle
11 // this macro is based on BDTControlPlots.C
12 void TMVA::BoostControlPlots(TString dataset, TString fin , Bool_t useTMVAStyle )
13 {
14  // set style and remove existing canvas'
15  TMVAGlob::Initialize( useTMVAStyle );
16 
17  // checks if file with name "fin" is already open, and if not opens one
18  TFile* file = TMVAGlob::OpenFile( fin );
19 
20  // get all titles of the method Boost
21  TList titles;
22  TString dirname="Method_Boost";
23  UInt_t ninst = TMVA::TMVAGlob::GetListOfTitles(dirname,titles,file->GetDirectory(dataset.Data()));
24  if (ninst==0) {
25  cout << "Could not locate directory 'Method_Boost' in file " << fin << endl;
26  return;
27  }
28  // loop over all titles
29  TIter keyIter(&titles);
30  TDirectory *boostdir;
31  TKey *key;
32  while ((key = TMVAGlob::NextKey(keyIter,"TDirectory"))) {
33  boostdir = (TDirectory *)key->ReadObj();
34  boostcontrolplots(dataset, boostdir );
35  }
36 }
37 
38 void TMVA::boostcontrolplots(TString dataset, TDirectory *boostdir ) {
39 
40  const Int_t nPlots = 6;
41 
42  Int_t width = 900;
43  Int_t height = 900;
44  char cn[100];
45  const TString titName = boostdir->GetName();
46  sprintf( cn, "cv_%s", titName.Data() );
47  TCanvas *c = new TCanvas( cn, Form( "%s Control Plots", titName.Data() ),
48  width, height );
49  c->Divide(2,4);
50 
51 
52  //TString hname[nPlots]={"Booster_BoostWeight","Booster_MethodWeight","Booster_ErrFraction","Booster_OrigErrFraction"};
53 
54  TString hname[nPlots]={"BoostWeight","MethodWeight","ErrFraction","SoverBtotal","SeparationGain", "SeparationGain"};
55 
56  // Note: the ROCIntegral plots are only filled for option "Boost_DetailedMonitoring=ture" currently not filled...
57  // TString hname[nPlots]={"BoostWeight","MethodWeight","ErrFraction","ROCIntegral_test"};
58 
59  for (Int_t i=0; i<nPlots; i++){
60  Int_t color = 4;
61  TH1 *h = (TH1*) boostdir->Get(hname[i]);
62  TString plotname = h->GetName();
63  h->SetMaximum(h->GetMaximum()*1.3);
64  h->SetMinimum( 0 );
65  h->SetMarkerColor(color);
66  h->SetMarkerSize( 0.7 );
67  h->SetMarkerStyle( 24 );
68  h->SetLineWidth(2);
69  h->SetLineColor(color);
70  h->Draw();
71  c->Update();
72  }
73 
74  // draw combined ROC plots
75 
76  TString hname_roctest[2] ={"ROCIntegral_test", "ROCIntegralBoosted_test"};
77  TString hname_roctrain[2]={"ROCIntegral_train", "ROCIntegralBoosted_train"};
78  TString htitle[2] = {"ROC integral of single classifier", "ROC integral of boosted method"};
79 
80  for (Int_t i=0; i<2; i++){
81  Int_t color = 4;
82  TPad * cPad = (TPad*)c->cd(nPlots+i+1);
83  TH1 *htest = (TH1*) boostdir->Get(hname_roctest[i]);
84  TH1 *htrain = (TH1*) boostdir->Get(hname_roctrain[i]);
85 
86  // check if filled
87  // Bool_t histFilled = (htest->GetMaximum() > 0 || htrain->GetMaximum() > 0);
88  Bool_t histFilled = (htest && htrain);
89 
90  if (!htest) htest = new TH1F("htest","",2,0,1);
91  if (!htrain) htrain = new TH1F("htrain","",2,0,1);
92 
93  htest->SetTitle(htitle[i]);
94  htest->SetMaximum(1.0);
95  htest->SetMinimum(0.0);
96  htest->SetMarkerColor(color);
97  htest->SetMarkerSize( 0.7 );
98  htest->SetMarkerStyle( 24 );
99  htest->SetLineWidth(2);
100  htest->SetLineColor(color);
101  htest->Draw();
102  htrain->SetMaximum(1.0);
103  htrain->SetMinimum(0.0);
104  htrain->SetMarkerColor(color-2);
105  htrain->SetMarkerSize( 0.7 );
106  htrain->SetMarkerStyle( 24 );
107  htrain->SetLineWidth(2);
108  htrain->SetLineColor(color-2);
109  htrain->Draw("same");
110 
111  if (histFilled) {
112  TLegend *legend= new TLegend( cPad->GetLeftMargin(),
113  0.2 + cPad->GetBottomMargin(),
114  cPad->GetLeftMargin() + 0.6,
115  cPad->GetBottomMargin() );
116  legend->AddEntry(htest, TString("testing sample"), "L");
117  legend->AddEntry(htrain, TString("training sample (orig. weights)"), "L");
118  legend->SetFillStyle( 1 );
119  legend->SetBorderSize(1);
120  legend->SetMargin( 0.3 );
121  legend->Draw("same");
122  }
123  else {
124  TText* t = new TText();
125  t->SetTextSize( 0.056 );
126  t->SetTextColor( 2 );
127  t->DrawTextNDC( .2, 0.6, "Use MethodBoost option: \"Boost_DetailedMonitoring\" " );
128  t->DrawTextNDC( .2, 0.51, "to fill this histograms" );
129  }
130 
131  c->Update();
132  }
133 
134  // write to file
135  TString fname = dataset+Form( "/plots/%s_ControlPlots", titName.Data() );
136  TMVAGlob::imgconv( c, fname );
137 
138 }
139 
140