Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
regression_averagedevs.cxx
Go to the documentation of this file.
1#include <limits>
2
4
5#include "TLatex.h"
6#include "TGraphErrors.h"
7#include "TFrame.h"
8
9/*
10 this macro plots the quadratic deviation of the estimated from the target value, averaged over the first nevt events in test sample (all if Nevt=-1)
11 a) normal average
12 b) truncated average, using best 90%
13 created January 2009, Eckhard von Toerne, University of Bonn, Germany
14*/
15
16void TMVA::regression_averagedevs(TString dataset,TString fin, Int_t Nevt, Bool_t useTMVAStyle )
17{
18 bool debug=false;
19 if (Nevt <0) Nevt=1000000;
20 TMVAGlob::Initialize( useTMVAStyle );
21 // checks if file with name "fin" is already open, and if not opens one
22 TFile* file = TMVAGlob::OpenFile( fin );
23 TList jobDirList;
24 TMVAGlob::GetListOfJobs((TFile*)file->GetDirectory(dataset.Data()),jobDirList);
25 if (jobDirList.GetSize()==0) {
26 cout << "error could not find jobs" << endl;
27 return;
28 }
29
30 Bool_t __PLOT_LOGO__ = kTRUE;
31 Bool_t __SAVE_IMAGE__ = kTRUE;
32
33 TDirectory* dir0 = (TDirectory*) (jobDirList.At(0));
34 //TDirectory* dir0 = (TDirectory*) (file->Get("InputVariables_Id"));
35 Int_t nTargets = TMVAGlob::GetNumberOfTargets( dir0);
36
37 if (debug) cout << "found targets " << nTargets<<endl;
38 TCanvas* c=0;
39 for (Int_t itrgt = 0 ; itrgt < nTargets; itrgt++){
40 if (debug) cout << "loop targets " << itrgt<<endl;
41 TString xtit = "Method";
42 TString ytit = "Average Quadratic Deviation";
43 TString ftit = ytit + " versus " + xtit + TString::Format(" for target %d",itrgt);
44 c = new TCanvas( TString::Format("c%d",itrgt), ftit , 50+20*itrgt, 10*itrgt, 750, 650 );
45
46 // global style settings
47 c->SetGrid();
48 c->SetTickx(1);
49 c->SetTicky(0);
50 c->SetTopMargin(0.28);
51 c->SetBottomMargin(0.1);
52
53 TString hNameRef = TString::Format("regression_average_devs_target%d",itrgt);
54
55 const Int_t maxMethods = 100;
56 // const Int_t maxTargets = 100;
57 Float_t m[4][maxMethods]; // h0 train-all, h1 train-90%, h2 test-all, h3 test-90%
58 Float_t em[4][maxMethods];
59 Float_t x[4][maxMethods];
60 Float_t ex[4][maxMethods];
61
62 TIter next(&jobDirList);
63 Float_t mymax=0., mymin=std::numeric_limits<float>::max();
64 TString mvaNames[maxMethods];
65 TDirectory *jobDir;
66 Int_t nMethods = 0;
67 // loop over all methods
68 while ( (jobDir = (TDirectory*)next()) ) {
69 TString methodTitle;
70 TMVAGlob::GetMethodTitle(methodTitle,jobDir);
71 mvaNames[nMethods]=methodTitle;
72 if (debug) cout << "--- Found directory for method: " << methodTitle << endl;
73 TIter keyIt(jobDir->GetListOfKeys());
74 TKey *histKey;
75 while ( (histKey = (TKey*)keyIt()) ) {
76 if (histKey->ReadObj()->InheritsFrom("TH1F") ){
77 TString s(histKey->ReadObj()->GetName());
78 if( !s.Contains("Quadr_Dev") ) continue;
79 if( !s.Contains(TString::Format("target_%d_",itrgt))) continue;
80 Int_t ihist = 0 ;
81 if( !s.Contains("best90perc") && s.Contains("train")) ihist=0;
82 if( s.Contains("best90perc") && s.Contains("train")) ihist=1;
83 if( !s.Contains("best90perc") && s.Contains("test")) ihist=2;
84 if( s.Contains("best90perc") && s.Contains("test")) ihist=3;
85 if (debug) cout <<"using histogram" << s << ", ihist="<<ihist<<endl;
86 TH1F* h = (TH1F*) (histKey->ReadObj());
87 m[ihist][nMethods] = sqrt(h->GetMean());
88 em[ihist][nMethods] = h->GetRMS()/(sqrt(h->GetEntries())*2.*h->GetMean());
89 x[ihist][nMethods] = nMethods+0.44+0.12*ihist;
90 ex[ihist][nMethods] = 0.001;
91 mymax= m[ihist][nMethods] > mymax ? m[ihist][nMethods] : mymax;
92 mymin= m[ihist][nMethods] < mymin ? m[ihist][nMethods] : mymin;
93 if (debug) cout << "m"<< ihist << "="<<m[ihist][nMethods]<<endl;
94 }
95 }
96 nMethods++;
97 }
98 TH1F* haveragedevs= new TH1F(TString::Format("haveragedevs%d",itrgt),ftit,nMethods,0.,nMethods);
99 for (int i=0;i<nMethods;i++) haveragedevs->GetXaxis()->SetBinLabel(i+1, mvaNames[i]);
100 haveragedevs->SetStats(0);
101 TGraphErrors* graphTrainAv= new TGraphErrors(nMethods,x[0],m[0],ex[0],em[0]);
102 TGraphErrors* graphTruncTrainAv= new TGraphErrors(nMethods,x[1],m[1],ex[1],em[1]);
103 TGraphErrors* graphTestAv= new TGraphErrors(nMethods,x[2],m[2],ex[2],em[2]);
104 TGraphErrors* graphTruncTestAv= new TGraphErrors(nMethods,x[3],m[3],ex[3],em[3]);
105
106 Double_t xmax = 1.2 * mymax;
107 Double_t xmin = 0.8 * mymin - (mymax - mymin)*0.05;
108 Double_t xheader = 0.2;
109 Double_t yheader = xmax*0.92;
110 xmin = xmin > 0.? xmin : 0.;
111 if (mymin > 1.e-20 && log10(mymax/mymin)>1.5){
112 c->SetLogy();
113 cout << "--- result differ significantly using log scale for display of regression results"<< endl;
114 xmax = 1.5 * xmax;
115 xmin = 0.75 * mymin;
116 yheader = xmax*0.78;
117 }
118 Float_t x0L = 0.03, y0H = 0.91;
119 Float_t dxL = 0.457-x0L, dyH = 0.14;
120 // TLegend *legend = new TLegend( x0L, y0H-dyH, x0L+dxL, y0H , "Average Deviation = (#sum_{evts} (f_{MVA} - f_{target})^{2} )^{1/2}");
121 TLegend *legend = new TLegend( x0L, y0H-dyH, x0L+dxL, y0H );
122 legend->SetTextSize( 0.035 );
123 legend->SetTextAlign(12);
124 legend->SetMargin( 0.1 );
125
126 TH1F *hr = c->DrawFrame(-1.,0.,nMethods+1, xmax);
127 cout << endl;
128 cout << "Training: Average Deviation between target " << itrgt <<" and estimate" << endl;
129 cout << TString::Format("%-15s%-15s%-15s", "Method","Average Dev.","trunc. Aver.(90%)") <<endl;
130 for (int i=0;i<nMethods;i++){
131 cout << TString::Format("%-15s:%#10.3g%#10.3g",
132 (const char*)mvaNames[i], m[0][i],m[1][i])<<endl;
133 // cout << mvaNames[i] << " " << m[0][i]<< " "<< m[1][i]<<endl;
134 hr->GetXaxis()->SetBinLabel(i+1," ");
135 }
136 cout << endl;
137 cout << "Testing: Average Deviation between target " << itrgt <<" and estimate" << endl;
138 cout << TString::Format("%-15s%-15s%-15s", "Method","Average Dev.","trunc. Aver.(90%)") <<endl;
139 for (int i=0;i<nMethods;i++){
140 cout << TString::Format("%-15s:%#10.3g%#10.3g",
141 (const char*)mvaNames[i], m[2][i],m[3][i])<<endl;
142 //cout << mvaNames[i] << " " << m[2][i]<< " "<< m[3][i]<<endl;
143 }
144
145 haveragedevs->SetMinimum(xmin);
146 haveragedevs->SetMaximum(xmax);
147 haveragedevs->SetXTitle("Method");
148 haveragedevs->SetYTitle("Deviation from target");
149 haveragedevs->Draw();
150 c->GetFrame()->SetFillColor(21);
151 c->GetFrame()->SetBorderSize(12);
152 graphTrainAv->SetMarkerSize(1.);
153 graphTrainAv->SetMarkerColor(kBlue);
154 graphTrainAv->SetMarkerStyle(25);
155 graphTrainAv->Draw("P");
156
157 graphTruncTrainAv->SetMarkerSize(1.);
158 graphTruncTrainAv->SetMarkerColor(kBlack);
159 graphTruncTrainAv->SetMarkerStyle(25);
160 graphTruncTrainAv->Draw("P");
161
162 graphTestAv->SetMarkerSize(1.);
163 graphTestAv->SetMarkerColor(kBlue);
164 graphTestAv->SetMarkerStyle(21);
165 graphTestAv->Draw("P");
166
167 graphTruncTestAv->SetMarkerSize(1.);
168 graphTruncTestAv->SetMarkerColor(kBlack);
169 graphTruncTestAv->SetMarkerStyle(21);
170 graphTruncTestAv->Draw("P");
171 legend->AddEntry(graphTrainAv,TString("Training Sample, Average Deviation"),"p");
172 legend->AddEntry(graphTruncTrainAv,TString("Training Sample, truncated Average Dev. (best 90%)"),"p");
173 legend->AddEntry(graphTestAv,TString("Test Sample, Average Deviation"),"p");
174 legend->AddEntry(graphTruncTestAv,TString("Test Sample, truncated Average Dev. (best 90%)"),"p");
175
176 legend->Draw();
177 TLatex legHeader;
178 legHeader.SetTextSize(0.035);
179 legHeader.SetTextAlign(12);
180 //legHeader.DrawLatex(x0L, y0H+0.01, "Average Deviation = (#sum (_{ } f_{MVA} - f_{target})^{2} )^{1/2}");
181 legHeader.DrawLatex(xheader, yheader, "Average Deviation = (#sum (_{ } f_{MVA} - f_{target})^{2} )^{1/2}");
182 // ============================================================
183
184 if (__PLOT_LOGO__) TMVAGlob::plot_logo();
185 // ============================================================
186
187 c->Update();
188 TString fname = dataset+"/plots/" + hNameRef;
189 if (__SAVE_IMAGE__) TMVAGlob::imgconv( c, fname );
190 } // end loop itrgt
191 return;
192}
193
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
bool Bool_t
Definition RtypesCore.h:63
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
@ kBlack
Definition Rtypes.h:65
@ kBlue
Definition Rtypes.h:66
float xmin
float xmax
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:45
virtual void SetTextAlign(Short_t align=11)
Set the text alignment.
Definition TAttText.h:42
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:47
virtual void SetBinLabel(Int_t bin, const char *label)
Set label for bin.
Definition TAxis.cxx:886
The Canvas class.
Definition TCanvas.h:23
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
TDirectory * GetDirectory(const char *apath, Bool_t printError=false, const char *funcname="GetDirectory") override
Find a directory named "apath".
Describe directory structure in memory.
Definition TDirectory.h:45
virtual TList * GetListOfKeys() const
Definition TDirectory.h:223
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:53
A TGraphErrors is a TGraph with error bars.
void Draw(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition TGraph.cxx:814
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:621
virtual void SetXTitle(const char *title)
Definition TH1.h:418
TAxis * GetXaxis()
Definition TH1.h:324
virtual void SetMaximum(Double_t maximum=-1111)
Definition TH1.h:403
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3066
virtual void SetMinimum(Double_t minimum=-1111)
Definition TH1.h:404
virtual void SetYTitle(const char *title)
Definition TH1.h:419
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TH1.cxx:8958
Book space in a file, create I/O buffers, to fill them, (un)compress them.
Definition TKey.h:28
virtual TObject * ReadObj()
To read a TObject* from the file.
Definition TKey.cxx:762
To draw Mathematical Formula.
Definition TLatex.h:18
TLatex * DrawLatex(Double_t x, Double_t y, const char *text)
Make a copy of this object with the new parameters And copy object attributes.
Definition TLatex.cxx:1943
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:317
void Draw(Option_t *option="") override
Draw this legend with its current attributes.
Definition TLegend.cxx:422
void SetMargin(Float_t margin)
Definition TLegend.h:69
A doubly linked list.
Definition TList.h:38
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:355
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:439
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:525
Basic string class.
Definition TString.h:139
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2378
RVec< PromoteType< T > > log10(const RVec< T > &v)
Definition RVec.hxx:1842
Double_t x[n]
Definition legend1.C:17
Double_t ex[n]
Definition legend1.C:17
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
UInt_t GetListOfJobs(TFile *file, TList &jobdirs)
Definition tmvaglob.cxx:615
void Initialize(Bool_t useTMVAStyle=kTRUE)
Definition tmvaglob.cxx:176
void GetMethodTitle(TString &name, TKey *ikey)
Definition tmvaglob.cxx:348
void plot_logo(Float_t v_scale=1.0, Float_t skew=1.0)
Definition tmvaglob.cxx:270
TFile * OpenFile(const TString &fin)
Definition tmvaglob.cxx:192
Int_t GetNumberOfTargets(TDirectory *dir)
Definition tmvaglob.cxx:403
void imgconv(TCanvas *c, const TString &fname)
Definition tmvaglob.cxx:212
void regression_averagedevs(TString dataset, TString fin, Int_t Nevt=-1, Bool_t useTMVAStyle=kTRUE)
TMarker m
Definition textangle.C:8