Logo ROOT   6.14/05
Reference Guide
rulevisHists.cxx
Go to the documentation of this file.
1 #include "TMVA/rulevisHists.h"
2 #include "TH2.h"
3 
4 
5 // This macro plots the distributions of the different input variables overlaid on
6 // the sum of importance per bin.
7 // The scale goes from violett (no importance) to red (high importance).
8 // Areas where many important rules are active, will thus be very red.
9 //
10 // input: - Input file (result from TMVA),
11 // - normal/decorrelated/PCA
12 // - use of TMVA plotting TStyle
13 void TMVA::rulevisHists(TString fin, TMVAGlob::TypeOfPlot type, bool useTMVAStyle)
14 {
15  // set style and remove existing canvas'
16  TMVAGlob::Initialize( useTMVAStyle );
17 
18  // checks if file with name "fin" is already open, and if not opens one
19  //TFile *file =
20  TMVAGlob::OpenFile( fin );
21 
22  // get all titles of the method rulefit
23  TList titles;
24  TString dirname = "Method_RuleFit";
25  UInt_t ninst = TMVAGlob::GetListOfTitles(dirname,titles);
26  if (ninst==0) return;
27 
28  // get top dir containing all hists of the variables
30  if (vardir==0) return;
31 
32  TDirectory* corrdir = TMVAGlob::GetCorrelationPlotsDir( type, vardir );
33  if (corrdir==0) return;
34 
35  // loop over all titles
36  TIter keyIter(&titles);
37  TDirectory *rfdir;
38  TKey *rfkey;
39  while ((rfkey = TMVAGlob::NextKey(keyIter,"TDirectory"))) {
40  rfdir = (TDirectory *)rfkey->ReadObj();
41  rulevisHists( rfdir, vardir, corrdir, type );
42  }
43 }
44 
45 void TMVA::rulevisHists( TDirectory *rfdir, TDirectory *vardir, TDirectory *corrdir, TMVAGlob::TypeOfPlot type) {
46  //
47  if (rfdir==0) return;
48  if (vardir==0) return;
49  if (corrdir==0) return;
50  //
51  const TString rfName = rfdir->GetName();
52  const TString maintitle = rfName + " : Rule Importance";
53  const TString rfNameOpt = "_RF2D_";
54  const TString outfname[TMVAGlob::kNumOfMethods] = { "rulevisHists",
55  "rulevisHists_decorr",
56  "rulevisCorr_pca",
57  "rulevisCorr_gaussdecorr" };
58 
59  const TString outputName = outfname[type]+"_"+rfdir->GetName();
60  //
61  TIter rfnext(rfdir->GetListOfKeys());
62  TKey *rfkey;
63  Double_t rfmax=0;//keep compiler quiet .. it IS in any case
64  Double_t rfmin=0;// later initialized
65  // Bool_t allEmpty=kTRUE;
67  while ((rfkey = (TKey*)rfnext())) {
68  // make sure, that we only look at histograms
69  TClass *cl = gROOT->GetClass(rfkey->GetClassName());
70  if (!cl->InheritsFrom("TH2F")) continue;
71  TH2F *hrf = (TH2F*)rfkey->ReadObj();
72  TString hname= hrf->GetName();
73  if (hname.Contains("__RF_")){ // found a new RF plot
74  Double_t valmin = hrf->GetMinimum();
75  Double_t valmax = hrf->GetMaximum();
76  if (first) {
77  rfmin=valmin;
78  rfmax=valmax;
79  first = kFALSE;
80  } else {
81  if (valmax>rfmax) rfmax=valmax;
82  if (valmin<rfmin) rfmin=valmin;
83  }
84  // if (hrf->GetEntries()>0) allEmpty=kFALSE;
85  }
86  }
87  if (first) {
88  cout << "ERROR: no RF plots found..." << endl;
89  return;
90  }
91 
92  const Int_t nContours = 100;
93  Double_t contourLevels[nContours];
94  Double_t dcl = (rfmax-rfmin)/Double_t(nContours-1);
95  //
96  for (Int_t i=0; i<nContours; i++) {
97  contourLevels[i] = rfmin+dcl*Double_t(i);
98  }
99 
100  ///////////////////////////
101  vardir->cd();
102 
103  // how many plots are in the directory?
104  Int_t noPlots = ((vardir->GetListOfKeys())->GetEntries()) / 2;
105 
106  // define Canvas layout here!
107  // default setting
108  Int_t xPad; // no of plots in x
109  Int_t yPad; // no of plots in y
110  Int_t width; // size of canvas
111  Int_t height;
112  switch (noPlots) {
113  case 1:
114  xPad = 1; yPad = 1; width = 500; height = 0.7*width; break;
115  case 2:
116  xPad = 2; yPad = 1; width = 600; height = 0.7*width; break;
117  case 3:
118  xPad = 3; yPad = 1; width = 900; height = 0.4*width; break;
119  case 4:
120  xPad = 2; yPad = 2; width = 600; height = width; break;
121  default:
122  xPad = 3; yPad = 2; width = 800; height = 0.7*width; break;
123  }
124  Int_t noPad = xPad * yPad ;
125 
126  // this defines how many canvases we need
127  const Int_t noCanvas = 1 + (Int_t)((noPlots - 0.001)/noPad);
128  TCanvas **c = new TCanvas*[noCanvas];
129  for (Int_t ic=0; ic<noCanvas; ic++) c[ic] = 0;
130 
131  // counter variables
132  Int_t countCanvas = 0;
133  Int_t countPad = 1;
134 
135  // loop over all objects in directory
136  TIter next(vardir->GetListOfKeys());
137  TKey *key;
138  //
139  first = kTRUE;
140 
141  while ((key = (TKey*)next())) {
142 
143  // make sure, that we only look at histograms
144  TClass *cl = gROOT->GetClass(key->GetClassName());
145  if (!cl->InheritsFrom("TH1")) continue;
146  TH1F* sig = (TH1F*)key->ReadObj();
147  TString hname= sig->GetName();
148 
149  // check for all signal histograms
150  if (hname.Contains("__S")){ // found a new signal plot
151  // create new canvas
152  if ((c[countCanvas]==NULL) || (countPad>noPad)) {
153  char cn[20];
154  sprintf( cn, "rulehist%d_", countCanvas+1 );
155  TString cname(cn);
156  cname += rfdir->GetName();
157  c[countCanvas] = new TCanvas( cname, maintitle,
158  countCanvas*50+200, countCanvas*20, width, height );
159  // style
160  c[countCanvas]->Divide(xPad,yPad);
161  countPad = 1;
162  }
163 
164  // save canvas to file
165  TPad *cPad = (TPad *)(c[countCanvas]->GetPad(countPad));
166  c[countCanvas]->cd(countPad);
167  countPad++;
168 
169  // find the corredponding background histo
170  TString bgname = hname;
171  bgname.ReplaceAll("__S","__B");
172  TKey *hkey = vardir->GetKey(bgname);
173  TH1F* bgd = (TH1F*)hkey->ReadObj();
174  if (bgd == NULL) {
175  cout << "ERROR!!! couldn't find backgroung histo for" << hname << endl;
176  //exit(1);
177  delete [] c;
178  return;
179  }
180 
181  TString rfname = hname;
182  rfname.ReplaceAll("__S","__RF");
183  TKey *hrfkey = rfdir->GetKey(rfname);
184  TH2F *hrf = (TH2F*)hrfkey->ReadObj();
185  // Double_t wv = hrf->GetMaximum();
186  // if (rfmax>0.0)
187  // hrf->Scale(1.0/rfmax);
188  hrf->SetMinimum(rfmin); // make sure it's zero -> for palette axis
189  hrf->SetMaximum(rfmax); // make sure max is 1.0 -> idem
190  hrf->SetContour(nContours,&contourLevels[0]);
191 
192  // this is set but not stored during plot creation in MVA_Factory
193  // TMVAGlob::SetSignalAndBackgroundStyle( sigK, bgd );
194  sig->SetFillStyle(3002);
195  sig->SetFillColor(1);
196  sig->SetLineColor(1);
197  sig->SetLineWidth(2);
198 
199  bgd->SetFillStyle(3554);
200  bgd->SetFillColor(1);
201  bgd->SetLineColor(1);
202  bgd->SetLineWidth(2);
203 
204  // chop off "signal"
205  TString title(hrf->GetTitle());
206  title.ReplaceAll("signal","");
207 
208  // finally plot and overlay
209  Float_t sc = 1.1;
210  if (countPad==2) sc = 1.3;
211  sig->SetMaximum( TMath::Max( sig->GetMaximum(), bgd->GetMaximum() )*sc );
212  Double_t smax = sig->GetMaximum();
213 
214  if (first) {
215  hrf->SetTitle( maintitle );
216  first = kFALSE;
217  } else {
218  hrf->SetTitle( "" );
219  }
220  hrf->Draw("colz ah");
221  TMVAGlob::SetFrameStyle( hrf, 1.2 );
222 
223  sig->Draw("same ah");
224  bgd->Draw("same ah");
225  // draw axis using range [0,smax]
226  hrf->GetXaxis()->SetTitle( title );
227  hrf->GetYaxis()->SetTitleOffset( 1.30 );
228  hrf->GetYaxis()->SetTitle("Events");
229  hrf->GetYaxis()->SetLimits(0,smax);
230  hrf->Draw("same axis");
231 
232  cPad->SetRightMargin(0.13);
233  cPad->Update();
234 
235  // Draw legend
236  if (countPad==2){
237  TLegend *legend= new TLegend( cPad->GetLeftMargin(),
238  1-cPad->GetTopMargin()-.18,
239  cPad->GetLeftMargin()+.4,
240  1-cPad->GetTopMargin() );
241  legend->AddEntry(sig,"Signal","F");
242  legend->AddEntry(bgd,"Background","F");
243  legend->Draw("same");
244  legend->SetBorderSize(1);
245  legend->SetMargin( 0.3 );
246  legend->SetFillColor(19);
247  legend->SetFillStyle(1);
248  }
249 
250  // save canvas to file
251  if (countPad > noPad) {
252  c[countCanvas]->Update();
253  TString fname = Form( "plots/%s_c%i", outputName.Data(), countCanvas+1 );
254  TMVAGlob::imgconv( c[countCanvas], fname );
255  // TMVAGlob::plot_logo(); // don't understand why this doesn't work ... :-(
256  countCanvas++;
257  }
258  }
259  }
260 
261  if (countPad <= noPad) {
262  c[countCanvas]->Update();
263  TString fname = Form( "plots/%s_c%i", outputName.Data(), countCanvas+1 );
264  TMVAGlob::imgconv( c[countCanvas], fname );
265  }
266  delete[] c;
267 }
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:294
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:43
virtual Double_t GetMaximum(Double_t maxval=FLT_MAX) const
Return maximum value smaller than maxval of bins in the range, unless the value has been overridden b...
Definition: TH1.cxx:7844
void imgconv(TCanvas *c, const TString &fname)
Definition: tmvaglob.cxx:212
Float_t GetLeftMargin() const
Definition: TAttPad.h:44
virtual TList * GetListOfKeys() const
Definition: TDirectory.h:150
virtual void SetMaximum(Double_t maximum=-1111)
Definition: TH1.h:390
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:23
virtual void SetLimits(Double_t xmin, Double_t xmax)
Definition: TAxis.h:154
float Float_t
Definition: RtypesCore.h:53
image html pict1_TGaxis_012 png width
Define new text attributes for the label number "labNum".
Definition: TGaxis.cxx:2551
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:687
virtual void SetContour(Int_t nlevels, const Double_t *levels=0)
Set the number and values of contour levels.
Definition: TH1.cxx:7785
TFile * OpenFile(const TString &fin)
Definition: tmvaglob.cxx:192
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:423
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:688
virtual void SetMinimum(Double_t minimum=-1111)
Definition: TH1.h:391
#define gROOT
Definition: TROOT.h:410
void SetFrameStyle(TH1 *frame, Float_t scale=1.0)
Definition: tmvaglob.cxx:77
Basic string class.
Definition: TString.h:131
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:567
virtual void Update()
Update pad.
Definition: TPad.cxx:2815
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
void SetMargin(Float_t margin)
Definition: TLegend.h:69
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition: TAttFill.h:39
TKey * NextKey(TIter &keyIter, TString className)
Definition: tmvaglob.cxx:357
TDirectory * GetCorrelationPlotsDir(TMVAGlob::TypeOfPlot type, TDirectory *dir=0)
Definition: tmvaglob.cxx:718
Book space in a file, create I/O buffers, to fill them, (un)compress them.
Definition: TKey.h:24
A doubly linked list.
Definition: TList.h:44
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
TDirectory * GetInputVariablesDir(TMVAGlob::TypeOfPlot type, TDirectory *dir=0)
Definition: tmvaglob.cxx:700
void Initialize(Bool_t useTMVAStyle=kTRUE)
Definition: tmvaglob.cxx:176
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2974
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:37
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:250
unsigned int UInt_t
Definition: RtypesCore.h:42
The most important graphics class in the ROOT system.
Definition: TPad.h:29
char * Form(const char *fmt,...)
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:75
TAxis * GetYaxis()
Definition: TH1.h:316
Bool_t InheritsFrom(const char *cl) const
Return kTRUE if this class inherits from a class with name "classname".
Definition: TClass.cxx:4688
virtual Double_t GetMinimum(Double_t minval=-FLT_MAX) const
Return minimum value larger than minval of bins in the range, unless the value has been overridden by...
Definition: TH1.cxx:7929
const Bool_t kFALSE
Definition: RtypesCore.h:88
The Canvas class.
Definition: TCanvas.h:31
virtual TKey * GetKey(const char *, Short_t=9999) const
Definition: TDirectory.h:148
double Double_t
Definition: RtypesCore.h:55
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:330
Describe directory structure in memory.
Definition: TDirectory.h:34
int type
Definition: TGX11.cxx:120
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:619
virtual void SetRightMargin(Float_t rightmargin)
Set Pad right margin in fraction of the pad width.
Definition: TAttPad.cxx:120
virtual TObject * ReadObj()
To read a TObject* from the file.
Definition: TKey.cxx:722
Float_t GetTopMargin() const
Definition: TAttPad.h:46
virtual Bool_t cd(const char *path=0)
Change current directory to "this" directory.
Definition: TDirectory.cxx:497
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1162
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:200
#define c(i)
Definition: RSha256.hxx:101
void rulevisHists(TString fin="TMVA.root", TMVAGlob::TypeOfPlot type=TMVAGlob::kNorm, bool useTMVAStyle=kTRUE)
virtual void SetTitle(const char *title)
See GetStatOverflows for more information.
Definition: TH1.cxx:6192
Definition: first.py:1
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2248
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
const Bool_t kTRUE
Definition: RtypesCore.h:87
virtual void SetBorderSize(Int_t bordersize=4)
Definition: TPave.h:70
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:315
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
const char * Data() const
Definition: TString.h:364
UInt_t GetListOfTitles(TDirectory *rfdir, TList &titles)
Definition: tmvaglob.cxx:636