ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 "TMVA/ResultsRegression.h"
29 
30 #include "TMVA/DataSet.h"
31 #include "TMVA/DataSetInfo.h"
32 #include "TMVA/MsgLogger.h"
33 #include "TMVA/Types.h"
34 #include "TMVA/VariableInfo.h"
35 
36 #include "TH1F.h"
37 #include "TH2F.h"
38 #include "TString.h"
39 
40 #include <vector>
41 
42 ////////////////////////////////////////////////////////////////////////////////
43 /// constructor
44 
46  : Results( dsi, resultsName ),
47  fLogger( new MsgLogger(Form("ResultsRegression%s",resultsName.Data()) , kINFO) )
48 {
49 }
50 
51 ////////////////////////////////////////////////////////////////////////////////
52 /// destructor
53 
55 {
56  delete fLogger;
57 }
58 
59 ////////////////////////////////////////////////////////////////////////////////
60 
61 void TMVA::ResultsRegression::SetValue( std::vector<Float_t>& value, Int_t ievt )
62 {
63  if (ievt >= (Int_t)fRegValues.size()) fRegValues.resize( ievt+1 );
64  fRegValues[ievt] = value;
65 }
66 
67 ////////////////////////////////////////////////////////////////////////////////
68 
70 {
71  DataSet* ds = GetDataSet();
72  ds->SetCurrentType( GetTreeType() );
73  const DataSetInfo* dsi = GetDataSetInfo();
74  TString name( Form("tgt_%d",tgtNum) );
75  VariableInfo vinf = dsi->GetTargetInfo(tgtNum);
76  Float_t xmin=0., xmax=0.;
77  if (truncate){
78  xmax = truncvalue;
79  }
80  else{
81  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
82  const Event* ev = ds->GetEvent(ievt);
83  std::vector<Float_t> regVal = fRegValues.at(ievt);
84  Float_t val = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
85  val *= val;
86  xmax = val> xmax? val: xmax;
87  }
88  }
89  xmax *= 1.1;
90  Int_t nbins = 500;
91  TH1F* h = new TH1F( name, name, nbins, xmin, xmax);
92  h->SetDirectory(0);
93  h->GetXaxis()->SetTitle("Quadratic Deviation");
94  h->GetYaxis()->SetTitle("Weighted Entries");
95 
96  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
97  const Event* ev = ds->GetEvent(ievt);
98  std::vector<Float_t> regVal = fRegValues.at(ievt);
99  Float_t val = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
100  val *= val;
101  Float_t weight = ev->GetWeight();
102  if (!truncate || val<=truncvalue ) h->Fill( val, weight);
103  }
104  return h;
105 }
106 
107 ////////////////////////////////////////////////////////////////////////////////
108 
110 {
111  DataSet* ds = GetDataSet();
112  ds->SetCurrentType( GetTreeType() );
113 
114  TString name( Form("tgt_%d_var_%d",tgtNum, varNum) );
115  const DataSetInfo* dsi = GetDataSetInfo();
116  Float_t xmin, xmax;
117  Bool_t takeTargets = kFALSE;
118  if (varNum >= dsi->GetNVariables()) {
119  takeTargets = kTRUE;
120  varNum -= dsi->GetNVariables();
121  }
122  if (!takeTargets) {
123  VariableInfo vinf = dsi->GetVariableInfo(varNum);
124  xmin = vinf.GetMin();
125  xmax = vinf.GetMax();
126 
127  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
128  const Event* ev = ds->GetEvent(ievt);
129  Float_t val = ev->GetValue(varNum);
130 
131  if (val < xmin ) xmin = val;
132  if (val > xmax ) xmax = val;
133  }
134 
135  }
136  else {
137  VariableInfo vinf = dsi->GetTargetInfo(varNum);
138  xmin = vinf.GetMin();
139  xmax = vinf.GetMax();
140 
141  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
142  const Event* ev = ds->GetEvent(ievt);
143  Float_t val = ev->GetTarget(varNum);
144 
145  if (val < xmin ) xmin = val;
146  if (val > xmax ) xmax = val;
147  }
148  }
149 
150  Float_t ymin = FLT_MAX;
151  Float_t ymax = -FLT_MAX;
152 
153  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
154  const Event* ev = ds->GetEvent(ievt);
155  std::vector<Float_t> regVal = fRegValues.at(ievt);
156 
157  Float_t diff = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
158  if (diff < ymin ) ymin = diff;
159  if (diff > ymax ) ymax = diff;
160  }
161 
162  Int_t nxbins = 50;
163  Int_t nybins = 50;
164 
165  Float_t epsilon = TMath::Abs(xmax-xmin)/((Float_t)nxbins-1);
166  xmin -= 1.01*epsilon;
167  xmax += 1.01*epsilon;
168 
169  epsilon = (ymax-ymin)/((Float_t)nybins-1);
170  ymin -= 1.01*epsilon;
171  ymax += 1.01*epsilon;
172 
173 
174  TH2F* h = new TH2F( name, name, nxbins, xmin, xmax, nybins, ymin, ymax );
175  h->SetDirectory(0);
176 
177  h->GetXaxis()->SetTitle( (takeTargets ? dsi->GetTargetInfo(varNum).GetTitle() : dsi->GetVariableInfo(varNum).GetTitle() ) );
178  TString varName( dsi->GetTargetInfo(tgtNum).GetTitle() );
179  TString yName( varName+TString("_{regression} - ") + varName+TString("_{true}") );
180  h->GetYaxis()->SetTitle( yName );
181 
182  for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
183  const Event* ev = ds->GetEvent(ievt);
184  std::vector<Float_t> regVal = fRegValues.at(ievt);
185 
186  Float_t xVal = (takeTargets?ev->GetTarget( varNum ):ev->GetValue( varNum ));
187  Float_t yVal = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
188 
189  h->Fill( xVal, yVal );
190  }
191 
192  return h;
193 }
194 
195 ////////////////////////////////////////////////////////////////////////////////
196 
198 {
199  Log() << kINFO << "Create variable histograms" << Endl;
200  const DataSetInfo* dsi = GetDataSetInfo();
201 
202  for (UInt_t ivar = 0; ivar < dsi->GetNVariables(); ivar++) {
203  for (UInt_t itgt = 0; itgt < dsi->GetNTargets(); itgt++) {
204  TH2F* h = DeviationAsAFunctionOf( ivar, itgt );
205  TString name( Form("%s_reg_var%d_rtgt%d",prefix.Data(),ivar,itgt) );
206  h->SetName( name );
207  h->SetTitle( name );
208  Store( h );
209  }
210  }
211  Log() << kINFO << "Create regression target histograms" << Endl;
212  for (UInt_t ivar = 0; ivar < dsi->GetNTargets(); ivar++) {
213  for (UInt_t itgt = 0; itgt < dsi->GetNTargets(); itgt++) {
214  TH2F* h = DeviationAsAFunctionOf( dsi->GetNVariables()+ivar, itgt );
215  TString name( Form("%s_reg_tgt%d_rtgt%d",prefix.Data(),ivar,itgt) );
216  h->SetName( name );
217  h->SetTitle( name );
218  Store( h );
219  }
220  }
221 
222  Log() << kINFO << "Create regression average deviation" << Endl;
223  for (UInt_t itgt = 0; itgt < dsi->GetNTargets(); itgt++) {
224  TH1F* h = QuadraticDeviation(itgt);
225  TString name( Form("%s_Quadr_Deviation_target_%d_",prefix.Data(),itgt) );
226  h->SetName( name );
227  h->SetTitle( name );
228  Double_t yq[1], xq[]={0.9};
229  h->GetQuantiles(1,yq,xq);
230  Store( h );
231 
232  TH1F* htrunc = QuadraticDeviation(itgt, true, yq[0]);
233  TString name2( Form("%s_Quadr_Dev_best90perc_target_%d_",prefix.Data(),itgt) );
234  htrunc->SetName( name2 );
235  htrunc->SetTitle( name2 );
236  Store( htrunc );
237  }
238  Log() << kINFO << "Results created" << Endl;
239 }
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
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
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
Double_t GetWeight() const
return the event weight - depending on whether the flag IgnoreNegWeightsInTraining is or not...
Definition: Event.cxx:376
TH2F * DeviationAsAFunctionOf(UInt_t varNum, UInt_t tgtNum)
int nbins[3]
void SetValue(std::vector< Float_t > &value, Int_t ievt)
UInt_t GetNVariables() const
Definition: DataSetInfo.h:128
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
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:186
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
Float_t GetValue(UInt_t ivar) const
return value of i'th variable
Definition: Event.cxx:231
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:287
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