ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
likelihoodrefs.cxx
Go to the documentation of this file.
1 #include "TMVA/likelihoodrefs.h"
2 #include <vector>
3 #include <string>
4 
5 
6 
7 // this macro plots the reference distribuions for the Likelihood
8 // methods for the various input variables used in TMVA (e.g. running
9 // TMVAnalysis.C). Signal and Background are plotted separately
10 
11 // input: - Input file (result from TMVA),
12 // - use of TMVA plotting TStyle
13 
14 
16  Bool_t newCanvas = kTRUE;
17 
18  const UInt_t maxCanvas = 200;
19  TCanvas** c = new TCanvas*[maxCanvas];
20  Int_t width = 670;
21  Int_t height = 380;
22 
23  // avoid duplicated printing
24  std::vector<std::string> hasBeenUsed;
25  const TString titName = lhdir->GetName();
26  UInt_t ic = -1;
27 
28  TIter next(lhdir->GetListOfKeys());
29  TKey *key;
30  while ((key = TMVAGlob::NextKey(next,"TH1"))) { // loop over all TH1
31  TH1 *h = (TH1*)key->ReadObj();
32  TH1F *b( 0 );
33  TString hname( h->GetName() );
34 
35  // avoid duplicated plotting
36  Bool_t found = kFALSE;
37  for (UInt_t j = 0; j < hasBeenUsed.size(); j++) {
38  if (hasBeenUsed[j] == hname.Data()) found = kTRUE;
39  }
40  if (!found) {
41 
42  // draw original plots
43  if (hname.EndsWith("_sig_nice")) {
44 
45  if (newCanvas) {
46  char cn[20];
47  sprintf( cn, "cv%d_%s", ic+1, titName.Data() );
48  ++ic;
49  TString n = hname;
50  c[ic] = new TCanvas( cn, Form( "%s reference for variable: %s",
51  titName.Data(),(n.ReplaceAll("_sig","")).Data() ),
52  ic*50+50, ic*20, width, height );
53  c[ic]->Divide(2,1);
54  newCanvas = kFALSE;
55  }
56 
57  // signal
58  Int_t color = 4;
59  TPad * cPad = (TPad*)c[ic]->cd(1);
60  TString plotname = hname;
61 
62  h->SetMaximum(h->GetMaximum()*1.3);
63  h->SetMinimum( 0 );
64  h->SetMarkerColor(color);
65  h->SetMarkerSize( 0.7 );
66  h->SetMarkerStyle( 24 );
67  h->SetLineWidth(1);
68  h->SetLineColor(color);
69  color++;
70  h->Draw("e1");
71  Double_t hSscale = 1.0/(h->GetSumOfWeights()*h->GetBinWidth(1));
72 
73  TLegend *legS= new TLegend( cPad->GetLeftMargin(),
74  1-cPad->GetTopMargin()-.14,
75  cPad->GetLeftMargin()+.77,
76  1-cPad->GetTopMargin() );
77  legS->SetBorderSize(1);
78  legS->AddEntry(h,"Input data (signal)","p");
79 
80  // background
81  TString bname( hname );
82  b = (TH1F*)lhdir->Get( bname.ReplaceAll("_sig","_bgd") );
83  cPad = (TPad*)c[ic]->cd(2);
84  color = 2;
85  b->SetMaximum(b->GetMaximum()*1.3);
86  b->SetMinimum( 0 );
87  b->SetLineWidth(1);
88  b->SetLineColor(color);
89  b->SetMarkerColor(color);
90  b->SetMarkerSize( 0.7 );
91  b->SetMarkerStyle( 24 );
92  b->Draw("e1");
93  Double_t hBscale = 1.0/(b->GetSumOfWeights()*b->GetBinWidth(1));
94  TLegend *legB= new TLegend( cPad->GetLeftMargin(),
95  1-cPad->GetTopMargin()-.14,
96  cPad->GetLeftMargin()+.77,
97  1-cPad->GetTopMargin() );
98  legB->SetBorderSize(1);
99  legB->AddEntry(b,"Input data (backgr.)","p");
100 
101  // register
102  hasBeenUsed.push_back( bname.Data() );
103 
104  // the PDFs --------------
105 
106  // check for splines
107  h = 0;
108  b = 0;
109  TString pname = hname; pname.ReplaceAll("_nice","");
110  for (int i=0; i<= 5; i++) {
111  TString hspline = pname + Form( "_smoothed_hist_from_spline%i", i );
112  h = (TH1F*)lhdir->Get( hspline );
113  if (h) {
114  b = (TH1F*)lhdir->Get( hspline.ReplaceAll("_sig","_bgd") );
115  break;
116  }
117  }
118 
119  // check for KDE
120  if (h == 0 && b == 0) {
121  TString hspline = pname +"_smoothed_hist_from_KDE";
122  h = (TH1F*)lhdir->Get( hspline );
123  if (h) {
124  b = (TH1F*)lhdir->Get( hspline.ReplaceAll("_sig","_bgd") );
125  }
126  }
127 
128  // found something ?
129  if (h == 0 || b == 0) {
130  cout << "--- likelihoodrefs.C: did not find spline for histogram: " << pname.Data() << endl;
131  }
132  else {
133 
134  Double_t pSscale = 1.0/(h->GetSumOfWeights()*h->GetBinWidth(1));
135  h->Scale( pSscale/hSscale );
136  color = 4;
137  c[ic]->cd(1);
138  h->SetLineWidth(2);
139  h->SetLineColor(color);
140  legS->AddEntry(h,"Estimated PDF (norm. signal)","l");
141  h->Draw("histsame");
142  legS->Draw();
143 
144  Double_t pBscale = 1.0/(b->GetSumOfWeights()*b->GetBinWidth(1));
145  b->Scale( pBscale/hBscale );
146  color = 2;
147  c[ic]->cd(2);
148  b->SetLineColor(color);
149  b->SetLineWidth(2);
150  legB->AddEntry(b,"Estimated PDF (norm. backgr.)","l");
151  b->Draw("histsame");
152 
153  // draw the legends
154  legB->Draw();
155 
156  hasBeenUsed.push_back( pname.Data() );
157  }
158 
159  c[ic]->Update();
160 
161  // write to file
162  TString fname = Form( "plots/%s_refs_c%i", titName.Data(), ic+1 );
163  TMVAGlob::imgconv( c[ic], fname );
164  // c[ic]->Update();
165 
166  newCanvas = kTRUE;
167  hasBeenUsed.push_back( hname.Data() );
168  }
169  }
170  }
171 }
172 
173 void TMVA::likelihoodrefs( TString fin , Bool_t useTMVAStyle )
174 {
175  // set style and remove existing canvas'
176  TMVAGlob::Initialize( useTMVAStyle );
177 
178  // checks if file with name "fin" is already open, and if not opens one
179  TMVAGlob::OpenFile( fin );
180 
181  // get all titles of the method likelihood
182  TList titles;
183  TString metlike="Method_Likelihood";
184  UInt_t ninst = TMVAGlob::GetListOfTitles(metlike,titles);
185  if (ninst==0) {
186  cout << "Could not locate directory 'Method_Likelihood' in file " << fin << endl;
187  return;
188  }
189  // loop over all titles
190  TIter keyIter(&titles);
191  TDirectory *lhdir;
192  TKey *key;
193  while ((key = TMVAGlob::NextKey(keyIter,"TDirectory"))) {
194  lhdir = (TDirectory *)key->ReadObj();
195  likelihoodrefs( lhdir );
196  }
197 }
198 
virtual void SetLineWidth(Width_t lwidth)
Definition: TAttLine.h:57
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6174
void imgconv(TCanvas *c, const TString &fname)
Definition: tmvaglob.cxx:212
virtual void SetMaximum(Double_t maximum=-1111)
Definition: TH1.h:394
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:35
virtual TList * GetListOfKeys() const
Definition: TDirectory.h:158
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
Definition: TDirectory.cxx:727
ClassImp(TSeqCollection) Int_t TSeqCollection TIter next(this)
Return index of object in collection.
return c
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
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:373
TH1 * h
Definition: legend2.C:5
tuple pname
Definition: tree.py:131
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
virtual void SetMinimum(Double_t minimum=-1111)
Definition: TH1.h:395
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
TKey * NextKey(TIter &keyIter, TString className)
Definition: tmvaglob.cxx:357
virtual Double_t GetBinWidth(Int_t bin) const
return bin width for 1D historam Better to use h1.GetXaxis().GetBinWidth(bin)
Definition: TH1.cxx:8492
const char * Data() const
Definition: TString.h:349
Float_t GetTopMargin() const
Definition: TAttPad.h:56
virtual void SetMarkerColor(Color_t mcolor=1)
Definition: TAttMarker.h:51
std::vector< std::vector< double > > Data
Book space in a file, create I/O buffers, to fill them, (un)compress them.
Definition: TKey.h:30
A doubly linked list.
Definition: TList.h:47
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
void Initialize(Bool_t useTMVAStyle=kTRUE)
Definition: tmvaglob.cxx:176
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
unsigned int UInt_t
Definition: RtypesCore.h:42
The most important graphics class in the ROOT system.
Definition: TPad.h:46
char * Form(const char *fmt,...)
virtual Double_t GetSumOfWeights() const
Return the sum of weights excluding under/overflows.
Definition: TH1.cxx:7354
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
virtual void SetMarkerStyle(Style_t mstyle=1)
Definition: TAttMarker.h:53
void likelihoodrefs(TDirectory *lhdir)
virtual void SetMarkerSize(Size_t msize=1)
Definition: TAttMarker.h:54
Float_t GetLeftMargin() const
Definition: TAttPad.h:54
The Canvas class.
Definition: TCanvas.h:48
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:280
Describe directory structure in memory.
Definition: TDirectory.h:44
The TH1 histogram class.
Definition: TH1.h:80
virtual TObject * ReadObj()
To read a TObject* from the file.
Definition: TKey.cxx:727
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:1073
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2179
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:7921
const Bool_t kTRUE
Definition: Rtypes.h:91
const char * cd(char *path=0)
Definition: rootalias.C:45
const Int_t n
Definition: legend1.C:16
virtual void SetBorderSize(Int_t bordersize=4)
Definition: TPave.h:82
UInt_t GetListOfTitles(TDirectory *rfdir, TList &titles)
Definition: tmvaglob.cxx:634