Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
PlotFoams.cxx
Go to the documentation of this file.
1 #include "TMVA/PlotFoams.h"
2 
3 #include "TControlBar.h"
4 #include "TMap.h"
5 #include "TVectorT.h"
6 #include "TLine.h"
7 #include "TPaveText.h"
10 
11 #include <sstream>
12 #include <string>
13 #include <cfloat>
14 
15 #include "TMVA/PDEFoam.h"
16 
17 void TMVA::PlotFoams( TString fileName,
18  bool useTMVAStyle )
19 {
20  cout << "read file: " << fileName << endl;
21  cout << "kValue = " << kValue << endl;
22  TFile *file = TFile::Open(fileName);
23 
24  // set style and remove existing canvas'
25  TMVAGlob::Initialize( useTMVAStyle );
26 
27  // create control bar
28  TControlBar* cbar = new TControlBar( "vertical", "Choose cell value for plot:", 50, 50 );
29  if ((gDirectory->Get("SignalFoam") && gDirectory->Get("BgFoam")) ||
30  gDirectory->Get("MultiTargetRegressionFoam")) {
31  TString macro = Form( "TMVA::Plot(\"%s\",%s, \"Event density\", %s)",
32  fileName.Data(), "TMVA::kValueDensity", (useTMVAStyle ? "kTRUE" : "kFALSE") );
33  cbar->AddButton( "Event density", macro, "Plot event density", "button" );
34  } else if (gDirectory->Get("DiscrFoam") || gDirectory->Get("MultiClassFoam0")){
35  TString macro = Form( "TMVA::Plot(\"%s\", %s, \"Discriminator\", %s)",
36  fileName.Data(), "TMVA::kValue", (useTMVAStyle ? "kTRUE" : "kFALSE") );
37  cbar->AddButton( "Discriminator", macro, "Plot discriminator", "button" );
38  } else if (gDirectory->Get("MonoTargetRegressionFoam")){
39  TString macro = Form( "TMVA::Plot(\"%s\", %s, \"Target\", %s)",
40  fileName.Data(), "TMVA::kValue", (useTMVAStyle ? "kTRUE" : "kFALSE") );
41  cbar->AddButton( "Target", macro, "Plot target", "button" );
42  } else {
43  cout << "Error: no foams found in file: " << fileName << endl;
44  return;
45  }
46 
47  TString macro_rms = Form( "TMVA::Plot(\"%s\", %s, \"Variance\", %s)",
48  fileName.Data(), "TMVA::kRms", (useTMVAStyle ? "kTRUE" : "kFALSE") );
49  cbar->AddButton( "Variance", macro_rms, "Plot variance", "button" );
50  TString macro_rms_ov_mean = Form( "TMVA::Plot(\"%s\", %s, \"Variance/Mean\", %s)",
51  fileName.Data(), "TMVA::kRmsOvMean", (useTMVAStyle ? "kTRUE" : "kFALSE") );
52  cbar->AddButton( "Variance/Mean", macro_rms_ov_mean, "Plot variance over mean", "button" );
53  TString macro_cell_tree = Form( "TMVA::PlotCellTree(\"%s\", \"Cell tree\", %s)",
54  fileName.Data(), (useTMVAStyle ? "kTRUE" : "kFALSE") );
55  cbar->AddButton( "Cell tree", macro_cell_tree, "Plot cell tree", "button" );
56 
57  cbar->Show();
58  file->Close();
59 }
60 
61 // foam plotting macro
62 void TMVA::Plot(TString fileName, TMVA::ECellValue cv, TString cv_long, bool useTMVAStyle )
63 {
64  cout << "read file: " << fileName << endl;
65  TFile *file = TFile::Open(fileName);
66 
67  gStyle->SetNumberContours(999);
68  if (useTMVAStyle) TMVAGlob::SetTMVAStyle();
69 
70  // fileNamed foams and foam type
71  TMVA::PDEFoam* SignalFoam = (TMVA::PDEFoam*) gDirectory->Get("SignalFoam");
72  TMVA::PDEFoam* BgFoam = (TMVA::PDEFoam*) gDirectory->Get("BgFoam");
73  TMVA::PDEFoam* DiscrFoam = (TMVA::PDEFoam*) gDirectory->Get("DiscrFoam");
74  TMVA::PDEFoam* MultiClassFoam0 = (TMVA::PDEFoam*) gDirectory->Get("MultiClassFoam0");
75  TMVA::PDEFoam* MonoTargetRegressionFoam = (TMVA::PDEFoam*) gDirectory->Get("MonoTargetRegressionFoam");
76  TMVA::PDEFoam* MultiTargetRegressionFoam = (TMVA::PDEFoam*) gDirectory->Get("MultiTargetRegressionFoam");
77  TList foam_list; // the foams and their captions
78  if (SignalFoam && BgFoam) {
79  foam_list.Add(new TPair(SignalFoam, new TObjString("Signal Foam")));
80  foam_list.Add(new TPair(BgFoam, new TObjString("Background Foam")));
81  } else if (DiscrFoam) {
82  foam_list.Add(new TPair(DiscrFoam, new TObjString("Discriminator Foam")));
83  } else if (MultiClassFoam0) {
84  UInt_t cls = 0;
85  TMVA::PDEFoam *fm = NULL;
86  while ((fm = (TMVA::PDEFoam*) gDirectory->Get(Form("MultiClassFoam%u", cls)))) {
87  foam_list.Add(new TPair(fm, new TObjString(Form("Discriminator Foam %u",cls))));
88  cls++;
89  }
90  } else if (MonoTargetRegressionFoam) {
91  foam_list.Add(new TPair(MonoTargetRegressionFoam,
92  new TObjString("MonoTargetRegression Foam")));
93  } else if (MultiTargetRegressionFoam) {
94  foam_list.Add(new TPair(MultiTargetRegressionFoam,
95  new TObjString("MultiTargetRegression Foam")));
96  } else {
97  cout << "ERROR: no Foams found in file: " << fileName << endl;
98  return;
99  }
100 
101  // loop over all foams and print out a debug message
102  TListIter foamIter(&foam_list);
103  TPair *fm_pair = NULL;
104  Int_t kDim = 0; // foam dimensions
105  while ((fm_pair = (TPair*) foamIter())) {
106  kDim = ((TMVA::PDEFoam*) fm_pair->Key())->GetTotDim();
107  cout << "Foam loaded: " << ((TObjString*) fm_pair->Value())->String()
108  << " (dimension = " << kDim << ")" << endl;
109  }
110 
111  // kernel to use for the projection
112  TMVA::PDEFoamKernelBase *kernel = new TMVA::PDEFoamKernelTrivial();
113 
114  // plot foams
115  if (kDim == 1) {
116  Plot1DimFoams(foam_list, cv, cv_long, kernel);
117  } else {
118  PlotNDimFoams(foam_list, cv, cv_long, kernel);
119  }
120 
121  file->Close();
122 }
123 
124 
125 void TMVA::Plot1DimFoams(TList& foam_list, TMVA::ECellValue cell_value,
126  const TString& cell_value_description,
127  TMVA::PDEFoamKernelBase* kernel)
128 {
129  // visualize a 1 dimensional PDEFoam via a histogram
130  TCanvas* canvas = NULL;
131  TH1D* projection = NULL;
132 
133  // loop over all foams and draw the histogram
134  TListIter it(&foam_list);
135  TPair* fm_pair = NULL; // the (foam, caption) pair
136  while ((fm_pair = (TPair*) it())) {
137  TMVA::PDEFoam* foam = (TMVA::PDEFoam*) fm_pair->Key();
138  if (!foam) continue;
139  TString foam_caption(((TObjString*) fm_pair->Value())->String());
140  TString variable_name(foam->GetVariableName(0)->String());
141 
142  canvas = new TCanvas(Form("canvas_%p",foam),
143  "1-dimensional PDEFoam", 400, 400);
144 
145  projection = foam->Draw1Dim(cell_value, 100, kernel);
146  projection->SetTitle(cell_value_description + " of " + foam_caption
147  + ";" + variable_name);
148  projection->Draw();
149  projection->SetDirectory(0);
150 
151  canvas->Update();
152  }
153 }
154 
155 
156 void TMVA::PlotNDimFoams(TList& foam_list, TMVA::ECellValue cell_value,
157  const TString& cell_value_description,
158  TMVA::PDEFoamKernelBase* kernel)
159 {
160  // draw 2 dimensional PDEFoam projections
161  TCanvas* canvas = NULL;
162  TH2D* projection = NULL;
163 
164  // loop over all foams and draw the projection
165  TListIter it(&foam_list);
166  TPair* fm_pair = NULL; // the (foam, caption) pair
167  while ((fm_pair = (TPair*) it())) {
168  TMVA::PDEFoam* foam = (TMVA::PDEFoam*) fm_pair->Key();
169  if (!foam) continue;
170  TString foam_caption(((TObjString*) fm_pair->Value())->String());
171  const Int_t kDim = ((TMVA::PDEFoam*) fm_pair->Key())->GetTotDim();
172 
173  // draw all possible projections (kDim*(kDim-1)/2)
174  for (Int_t i = 0; i < kDim; ++i) {
175  for (Int_t k = i + 1; k < kDim; ++k) {
176 
177  canvas = new TCanvas(Form("canvas_%p_%i:%i", foam, i, k),
178  Form("Foam projections %i:%i", i, k),
179  (Int_t)(400/(1.-0.2)), 400);
180  canvas->SetRightMargin(0.2);
181 
182  TString title = Form("%s of %s: Projection %s:%s;%s;%s",
183  cell_value_description.Data(),
184  foam_caption.Data(),
185  foam->GetVariableName(i)->String().Data(),
186  foam->GetVariableName(k)->String().Data(),
187  foam->GetVariableName(i)->String().Data(),
188  foam->GetVariableName(k)->String().Data());
189 
190  projection = foam->Project2(i, k, cell_value, kernel);
191  projection->SetTitle(title);
192  projection->Draw("COLZ");
193  projection->SetDirectory(0);
194 
195  canvas->Update();
196  }
197  }
198  } // loop over foams
199 }
200 
201 
202 void TMVA::PlotCellTree(TString fileName, TString cv_long, bool useTMVAStyle )
203 {
204  // Draw the PDEFoam cell tree
205 
206  cout << "read file: " << fileName << endl;
207  TFile *file = TFile::Open(fileName);
208 
209  if (useTMVAStyle) TMVAGlob::SetTMVAStyle();
210 
211  // find foams
212  TListIter foamIter(gDirectory->GetListOfKeys());
213  TKey *foam_key = NULL; // the foam key
214  TCanvas *canv = NULL; // the canvas
215  while ((foam_key = (TKey*) foamIter())) {
216  TString name(foam_key->GetName());
217  TString class_name(foam_key->GetClassName());
218  if (!class_name.Contains("PDEFoam"))
219  continue;
220  cout << "PDEFoam found: " << class_name << " " << name << endl;
221 
222  // read the foam
223  TMVA::PDEFoam *foam = (TMVA::PDEFoam*) foam_key->ReadObj();
224  canv = new TCanvas(Form("canvas_%s",name.Data()),
225  Form("%s of %s",cv_long.Data(),name.Data()), 640, 480);
226  canv->cd();
227  // get cell tree depth
228  const UInt_t depth = foam->GetRootCell()->GetTreeDepth();
229  const Double_t ystep = 1.0 / depth;
230  DrawCell(foam->GetRootCell(), foam, 0.5, 1.-0.5*ystep, 0.25, ystep);
231  }
232 
233  file->Close();
234 }
235 
236 void TMVA::DrawCell( TMVA::PDEFoamCell *cell, TMVA::PDEFoam *foam,
237  Double_t x, Double_t y,
238  Double_t xscale, Double_t yscale )
239 {
240  // recursively draw cell and it's daughters
241 
242  Float_t xsize = xscale*1.5;
243  Float_t ysize = yscale/3;
244  if (xsize > 0.15) xsize=0.1; //xscale/2;
245  if (cell->GetDau0() != NULL) {
246  TLine *a1 = new TLine(x-xscale/4, y-ysize, x-xscale, y-ysize*2);
247  a1->SetLineWidth(2);
248  a1->Draw();
249  DrawCell(cell->GetDau0(), foam, x-xscale, y-yscale, xscale/2, yscale);
250  }
251  if (cell->GetDau1() != NULL){
252  TLine *a1 = new TLine(x+xscale/4, y-ysize, x+xscale, y-ysize*2);
253  a1->SetLineWidth(2);
254  a1->Draw();
255  DrawCell(cell->GetDau1(), foam, x+xscale, y-yscale, xscale/2, yscale);
256  }
257 
258  TPaveText *t = new TPaveText(x-xsize, y-ysize, x+xsize, y+ysize, "NDC");
259 
260  t->SetBorderSize(1);
261  t->SetFillStyle(1);
262 
263  // draw all cell elements
264  t->AddText( Form("Intg=%.5f", cell->GetIntg()) );
265  t->AddText( Form("Var=%.5f", cell->GetDriv()) );
266  TVectorD *vec = (TVectorD*) cell->GetElement();
267  if (vec) {
268  for (Int_t i = 0; i < vec->GetNrows(); ++i)
269  t->AddText( Form("E[%i]=%.5f", i, (*vec)[i]) );
270  }
271 
272  if (cell->GetStat() != 1) {
273  // cell is inactive --> draw split point
274  t->SetFillColor( TColor::GetColor("#BBBBBB") );
275  t->SetTextColor( TColor::GetColor("#000000") );
276 
277  // cell position and size
278  TMVA::PDEFoamVect cellPosi(foam->GetTotDim()), cellSize(foam->GetTotDim());
279  cell->GetHcub(cellPosi, cellSize);
280  Int_t kBest = cell->GetBest(); // best division variable
281  Double_t xBest = cell->GetXdiv(); // best division point
282  t->AddText( Form("dim=%i", kBest) );
283  t->AddText( Form("cut=%.5g", foam->VarTransformInvers(kBest,cellPosi[kBest] + xBest*cellSize[kBest])) );
284  } else {
285  t->SetFillColor( TColor::GetColor("#DD0033") );
286  t->SetTextColor( TColor::GetColor("#FFFFFF") );
287  }
288 
289  t->Draw();
290 }