ROOT  6.06/09
Reference Guide
ResultsRegression.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Andreas Hoecker, Peter Speckmayer, Joerg Stelzer, Helge Voss
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : ResultsRegression *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Implementation (see header for description) *
12  * *
13  * Authors (alphabetical): *
14  * Andreas Hoecker <Andreas.Hocker@cern.ch> - CERN, Switzerland *
15  * Peter Speckmayer <Peter.Speckmayer@cern.ch> - CERN, Switzerland *
16  * Joerg Stelzer <Joerg.Stelzer@cern.ch> - CERN, Switzerland *
17  * Helge Voss <Helge.Voss@cern.ch> - MPI-K Heidelberg, Germany *
18  * *
19  * Copyright (c) 2006: *
20  * CERN, Switzerland *
21  * MPI-K Heidelberg, Germany *
22  * *
23  * Redistribution and use in source and binary forms, with or without *
24  * modification, are permitted according to the terms listed in LICENSE *
25  * (http://tmva.sourceforge.net/LICENSE) *
26  **********************************************************************************/
27 
28 #include <vector>
29 
30 #include "TMVA/ResultsRegression.h"
31 #include "TMVA/MsgLogger.h"
32 #include "TMVA/DataSet.h"
33 
34 ////////////////////////////////////////////////////////////////////////////////
35 /// constructor
36 
38  : Results( dsi, resultsName ),
39  fLogger( new MsgLogger(Form("ResultsRegression%s",resultsName.Data()) , kINFO) )
40 {
41 }
42 
43 ////////////////////////////////////////////////////////////////////////////////
44 /// destructor
45 
47 {
48  delete fLogger;
49 }
50 
51 ////////////////////////////////////////////////////////////////////////////////
52 
53 void TMVA::ResultsRegression::SetValue( std::vector<Float_t>& value, Int_t ievt )
54 {
55  if (ievt >= (Int_t)fRegValues.size()) fRegValues.resize( ievt+1 );
56  fRegValues[ievt] = value;
57 }
58 
59 ////////////////////////////////////////////////////////////////////////////////
60 
62 {
63  DataSet* ds = GetDataSet();
64  ds->SetCurrentType( GetTreeType() );
65  const DataSetInfo* dsi = GetDataSetInfo();
66  TString name( Form("tgt_%d",tgtNum) );
67  VariableInfo vinf = dsi->GetTargetInfo(tgtNum);
68  Float_t xmin=0., xmax=0.;
69  if (truncate){
70  xmax = truncvalue;
71  }
72  else{
73  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
74  const Event* ev = ds->GetEvent(ievt);
75  std::vector<Float_t> regVal = fRegValues.at(ievt);
76  Float_t val = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
77  val *= val;
78  xmax = val> xmax? val: xmax;
79  }
80  }
81  xmax *= 1.1;
82  Int_t nbins = 500;
83  TH1F* h = new TH1F( name, name, nbins, xmin, xmax);
84  h->SetDirectory(0);
85  h->GetXaxis()->SetTitle("Quadratic Deviation");
86  h->GetYaxis()->SetTitle("Weighted Entries");
87 
88  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
89  const Event* ev = ds->GetEvent(ievt);
90  std::vector<Float_t> regVal = fRegValues.at(ievt);
91  Float_t val = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
92  val *= val;
93  Float_t weight = ev->GetWeight();
94  if (!truncate || val<=truncvalue ) h->Fill( val, weight);
95  }
96  return h;
97 }
98 
99 ////////////////////////////////////////////////////////////////////////////////
100 
102 {
103  DataSet* ds = GetDataSet();
104  ds->SetCurrentType( GetTreeType() );
105 
106  TString name( Form("tgt_%d_var_%d",tgtNum, varNum) );
107  const DataSetInfo* dsi = GetDataSetInfo();
108  Float_t xmin, xmax;
109  Bool_t takeTargets = kFALSE;
110  if (varNum >= dsi->GetNVariables()) {
111  takeTargets = kTRUE;
112  varNum -= dsi->GetNVariables();
113  }
114  if (!takeTargets) {
115  VariableInfo vinf = dsi->GetVariableInfo(varNum);
116  xmin = vinf.GetMin();
117  xmax = vinf.GetMax();
118 
119  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
120  const Event* ev = ds->GetEvent(ievt);
121  Float_t val = ev->GetValue(varNum);
122 
123  if (val < xmin ) xmin = val;
124  if (val > xmax ) xmax = val;
125  }
126 
127  }
128  else {
129  VariableInfo vinf = dsi->GetTargetInfo(varNum);
130  xmin = vinf.GetMin();
131  xmax = vinf.GetMax();
132 
133  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
134  const Event* ev = ds->GetEvent(ievt);
135  Float_t val = ev->GetTarget(varNum);
136 
137  if (val < xmin ) xmin = val;
138  if (val > xmax ) xmax = val;
139  }
140  }
141 
142  Float_t ymin = FLT_MAX;
143  Float_t ymax = -FLT_MAX;
144 
145  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
146  const Event* ev = ds->GetEvent(ievt);
147  std::vector<Float_t> regVal = fRegValues.at(ievt);
148 
149  Float_t diff = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
150  if (diff < ymin ) ymin = diff;
151  if (diff > ymax ) ymax = diff;
152  }
153 
154  Int_t nxbins = 50;
155  Int_t nybins = 50;
156 
157  Float_t epsilon = TMath::Abs(xmax-xmin)/((Float_t)nxbins-1);
158  xmin -= 1.01*epsilon;
159  xmax += 1.01*epsilon;
160 
161  epsilon = (ymax-ymin)/((Float_t)nybins-1);
162  ymin -= 1.01*epsilon;
163  ymax += 1.01*epsilon;
164 
165 
166  TH2F* h = new TH2F( name, name, nxbins, xmin, xmax, nybins, ymin, ymax );
167  h->SetDirectory(0);
168 
169  h->GetXaxis()->SetTitle( (takeTargets ? dsi->GetTargetInfo(varNum).GetTitle() : dsi->GetVariableInfo(varNum).GetTitle() ) );
170  TString varName( dsi->GetTargetInfo(tgtNum).GetTitle() );
171  TString yName( varName+TString("_{regression} - ") + varName+TString("_{true}") );
172  h->GetYaxis()->SetTitle( yName );
173 
174  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
175  const Event* ev = ds->GetEvent(ievt);
176  std::vector<Float_t> regVal = fRegValues.at(ievt);
177 
178  Float_t xVal = (takeTargets?ev->GetTarget( varNum ):ev->GetValue( varNum ));
179  Float_t yVal = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
180 
181  h->Fill( xVal, yVal );
182  }
183 
184  return h;
185 }
186 
187 ////////////////////////////////////////////////////////////////////////////////
188 
190 {
191  Log() << kINFO << "Create variable histograms" << Endl;
192  const DataSetInfo* dsi = GetDataSetInfo();
193 
194  for (UInt_t ivar = 0; ivar < dsi->GetNVariables(); ivar++) {
195  for (UInt_t itgt = 0; itgt < dsi->GetNTargets(); itgt++) {
196  TH2F* h = DeviationAsAFunctionOf( ivar, itgt );
197  TString name( Form("%s_reg_var%d_rtgt%d",prefix.Data(),ivar,itgt) );
198  h->SetName( name );
199  h->SetTitle( name );
200  Store( h );
201  }
202  }
203  Log() << kINFO << "Create regression target histograms" << Endl;
204  for (UInt_t ivar = 0; ivar < dsi->GetNTargets(); ivar++) {
205  for (UInt_t itgt = 0; itgt < dsi->GetNTargets(); itgt++) {
206  TH2F* h = DeviationAsAFunctionOf( dsi->GetNVariables()+ivar, itgt );
207  TString name( Form("%s_reg_tgt%d_rtgt%d",prefix.Data(),ivar,itgt) );
208  h->SetName( name );
209  h->SetTitle( name );
210  Store( h );
211  }
212  }
213 
214  Log() << kINFO << "Create regression average deviation" << Endl;
215  for (UInt_t itgt = 0; itgt < dsi->GetNTargets(); itgt++) {
216  TH1F* h = QuadraticDeviation(itgt);
217  TString name( Form("%s_Quadr_Deviation_target_%d_",prefix.Data(),itgt) );
218  h->SetName( name );
219  h->SetTitle( name );
220  Double_t yq[1], xq[]={0.9};
221  h->GetQuantiles(1,yq,xq);
222  Store( h );
223 
224  TH1F* htrunc = QuadraticDeviation(itgt, true, yq[0]);
225  TString name2( Form("%s_Quadr_Dev_best90perc_target_%d_",prefix.Data(),itgt) );
226  htrunc->SetName( name2 );
227  htrunc->SetTitle( name2 );
228  Store( htrunc );
229  }
230  Log() << kINFO << "Results created" << Endl;
231 }
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
float xmin
Definition: THbookFile.cxx:93
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
float Float_t
Definition: RtypesCore.h:53
virtual void SetDirectory(TDirectory *dir)
By default when an histogram is created, it is added to the list of histogram objects in the current ...
Definition: TH1.cxx:8266
THist< 2, float > TH2F
Definition: THist.h:321
float ymin
Definition: THbookFile.cxx:93
Double_t GetMin() const
Definition: VariableInfo.h:69
TH1 * h
Definition: legend2.C:5
virtual Int_t GetQuantiles(Int_t nprobSum, Double_t *q, const Double_t *probSum=0)
Compute Quantiles for this histogram Quantile x_q of a probability distribution Function F is defined...
Definition: TH1.cxx:4182
UInt_t GetNTargets() const
Definition: DataSetInfo.h:129
THist< 1, float > TH1F
Definition: THist.h:315
Double_t GetMax() const
Definition: VariableInfo.h:70
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
TH2F * DeviationAsAFunctionOf(UInt_t varNum, UInt_t tgtNum)
int nbins[3]
void SetValue(std::vector< Float_t > &value, Int_t ievt)
Double_t GetWeight() const
return the event weight - depending on whether the flag IgnoreNegWeightsInTraining is or not...
Definition: Event.cxx:376
UInt_t GetNVariables() const
Definition: DataSetInfo.h:128
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Float_t GetValue(UInt_t ivar) const
return value of i'th variable
Definition: Event.cxx:231
const char * Data() const
Definition: TString.h:349
std::vector< std::vector< double > > Data
TFileCollection * GetDataSet(const char *ds, const char *server="")
GetDataSet wrapper.
Definition: pq2wrappers.cxx:87
float ymax
Definition: THbookFile.cxx:93
const TString & GetTitle() const
Definition: VariableInfo.h:65
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:256
VariableInfo & GetTargetInfo(Int_t i)
Definition: DataSetInfo.h:119
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
TAxis * GetYaxis()
Definition: TH1.h:320
float xmax
Definition: THbookFile.cxx:93
ResultsRegression(const DataSetInfo *dsi, TString resultsName)
constructor
REAL epsilon
Definition: triangle.c:617
virtual void SetName(const char *name)
Change the name of this histogram.
Definition: TH1.cxx:8288
const Event * GetEvent() const
Definition: DataSet.cxx:180
void SetCurrentType(Types::ETreeType type) const
Definition: DataSet.h:111
double Double_t
Definition: RtypesCore.h:55
Long64_t GetNEvents(Types::ETreeType type=Types::kMaxTreeType) const
Definition: DataSet.h:225
VariableInfo & GetVariableInfo(Int_t i)
Definition: DataSetInfo.h:114
#define name(a, b)
Definition: linkTestLib0.cpp:5
Float_t GetTarget(UInt_t itgt) const
Definition: Event.h:101
void CreateDeviationHistograms(TString prefix)
virtual void SetTitle(const char *title)
Change (i.e.
Definition: TH1.cxx:6268
const Bool_t kTRUE
Definition: Rtypes.h:91
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:285
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
Definition: TNamed.cxx:152
float value
Definition: math.cpp:443
Definition: math.cpp:60
TH1F * QuadraticDeviation(UInt_t tgtNum, Bool_t truncate=false, Double_t truncvalue=0.)
TAxis * GetXaxis()
Definition: TH1.h:319