Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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/*! \class TMVA::ResultsRegression
29\ingroup TMVA
30Class that is the base-class for a vector of result
31*/
33
34#include "TMVA/DataSet.h"
35#include "TMVA/DataSetInfo.h"
36#include "TMVA/MsgLogger.h"
37#include "TMVA/Results.h"
38#include "TMVA/Types.h"
39#include "TMVA/VariableInfo.h"
40
41#include "TH1F.h"
42#include "TH2F.h"
43#include "TString.h"
44
45#include <vector>
46
47
48////////////////////////////////////////////////////////////////////////////////
49/// constructor
50
52 : Results( dsi, resultsName ),
53 fLogger( new MsgLogger(TString::Format("ResultsRegression%s",resultsName.Data()).Data() , kINFO) )
54{
55}
56
57////////////////////////////////////////////////////////////////////////////////
58/// destructor
59
61{
62 delete fLogger;
63}
64
65////////////////////////////////////////////////////////////////////////////////
66
67void TMVA::ResultsRegression::SetValue( std::vector<Float_t>& value, Int_t ievt )
68{
69 if (ievt >= (Int_t)fRegValues.size()) fRegValues.resize( ievt+1 );
70 fRegValues[ievt] = value;
71}
72
73////////////////////////////////////////////////////////////////////////////////
74
76{
77 DataSet* ds = GetDataSet();
78 ds->SetCurrentType( GetTreeType() );
79 const DataSetInfo* dsi = GetDataSetInfo();
80 TString name = TString::Format("tgt_%d",tgtNum);
81 VariableInfo vinf = dsi->GetTargetInfo(tgtNum);
82 Float_t xmin=0., xmax=0.;
83 if (truncate){
84 xmax = truncvalue;
85 }
86 else{
87 for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
88 const Event* ev = ds->GetEvent(ievt);
89 std::vector<Float_t> regVal = fRegValues.at(ievt);
90 Float_t val = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
91 val *= val;
92 xmax = val> xmax? val: xmax;
93 }
94 }
95 xmax *= 1.1;
96 Int_t nbins = 500;
97 TH1F* h = new TH1F( name, name, nbins, xmin, xmax);
98 h->SetDirectory(nullptr);
99 h->GetXaxis()->SetTitle("Quadratic Deviation");
100 h->GetYaxis()->SetTitle("Weighted Entries");
101
102 for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
103 const Event* ev = ds->GetEvent(ievt);
104 std::vector<Float_t> regVal = fRegValues.at(ievt);
105 Float_t val = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
106 val *= val;
107 Float_t weight = ev->GetWeight();
108 if (!truncate || val<=truncvalue ) h->Fill( val, weight);
109 }
110 return h;
111}
112
113////////////////////////////////////////////////////////////////////////////////
114
116{
117 DataSet* ds = GetDataSet();
118 ds->SetCurrentType( GetTreeType() );
119
120 TString name = TString::Format("tgt_%d_var_%d",tgtNum, varNum);
121 const DataSetInfo* dsi = GetDataSetInfo();
123 Bool_t takeTargets = kFALSE;
124 if (varNum >= dsi->GetNVariables()) {
125 takeTargets = kTRUE;
126 varNum -= dsi->GetNVariables();
127 }
128 if (!takeTargets) {
129 VariableInfo vinf = dsi->GetVariableInfo(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->GetValue(varNum);
136
137 if (val < xmin ) xmin = val;
138 if (val > xmax ) xmax = val;
139 }
140
141 }
142 else {
143 VariableInfo vinf = dsi->GetTargetInfo(varNum);
144 xmin = vinf.GetMin();
145 xmax = vinf.GetMax();
146
147 for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
148 const Event* ev = ds->GetEvent(ievt);
149 Float_t val = ev->GetTarget(varNum);
150
151 if (val < xmin ) xmin = val;
152 if (val > xmax ) xmax = val;
153 }
154 }
155
156 Float_t ymin = FLT_MAX;
157 Float_t ymax = -FLT_MAX;
158
159 for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
160 const Event* ev = ds->GetEvent(ievt);
161 std::vector<Float_t> regVal = fRegValues.at(ievt);
162
163 Float_t diff = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
164 if (diff < ymin ) ymin = diff;
165 if (diff > ymax ) ymax = diff;
166 }
167
168 Int_t nxbins = 50;
169 Int_t nybins = 50;
170
171 Float_t epsilon = TMath::Abs(xmax-xmin)/((Float_t)nxbins-1);
172 xmin -= 1.01*epsilon;
173 xmax += 1.01*epsilon;
174
175 epsilon = (ymax-ymin)/((Float_t)nybins-1);
176 ymin -= 1.01*epsilon;
177 ymax += 1.01*epsilon;
178
179
180 TH2F* h = new TH2F( name, name, nxbins, xmin, xmax, nybins, ymin, ymax );
181 h->SetDirectory(nullptr);
182
183 h->GetXaxis()->SetTitle( (takeTargets ? dsi->GetTargetInfo(varNum).GetTitle() : dsi->GetVariableInfo(varNum).GetTitle() ) );
184 TString varName( dsi->GetTargetInfo(tgtNum).GetTitle() );
185 TString yName( varName+TString("_{regression} - ") + varName+TString("_{true}") );
186 h->GetYaxis()->SetTitle( yName );
187
188 for (Int_t ievt=0; ievt<ds->GetNEvents(); ievt++) {
189 const Event* ev = ds->GetEvent(ievt);
190 std::vector<Float_t> regVal = fRegValues.at(ievt);
191
192 Float_t xVal = (takeTargets?ev->GetTarget( varNum ):ev->GetValue( varNum ));
193 Float_t yVal = regVal.at( tgtNum ) - ev->GetTarget( tgtNum );
194
195 h->Fill( xVal, yVal );
196 }
197
198 return h;
199}
200
201////////////////////////////////////////////////////////////////////////////////
202
204{
205 Log() << kINFO << "Create variable histograms" << Endl;
206 const DataSetInfo* dsi = GetDataSetInfo();
207
208 for (UInt_t ivar = 0; ivar < dsi->GetNVariables(); ivar++) {
209 for (UInt_t itgt = 0; itgt < dsi->GetNTargets(); itgt++) {
210 TH2F* h = DeviationAsAFunctionOf( ivar, itgt );
211 TString name = TString::Format("%s_reg_var%d_rtgt%d",prefix.Data(),ivar,itgt);
212 h->SetName( name );
213 h->SetTitle( name );
214 Store( h );
215 }
216 }
217 Log() << kINFO << "Create regression target histograms" << Endl;
218 for (UInt_t ivar = 0; ivar < dsi->GetNTargets(); ivar++) {
219 for (UInt_t itgt = 0; itgt < dsi->GetNTargets(); itgt++) {
220 TH2F* h = DeviationAsAFunctionOf( dsi->GetNVariables()+ivar, itgt );
221 TString name = TString::Format("%s_reg_tgt%d_rtgt%d",prefix.Data(),ivar,itgt);
222 h->SetName( name );
223 h->SetTitle( name );
224 Store( h );
225 }
226 }
227
228 Log() << kINFO << "Create regression average deviation" << Endl;
229 for (UInt_t itgt = 0; itgt < dsi->GetNTargets(); itgt++) {
230 TH1F* h = QuadraticDeviation(itgt);
231 TString name = TString::Format("%s_Quadr_Deviation_target_%d_",prefix.Data(),itgt);
232 h->SetName( name );
233 h->SetTitle( name );
234 Double_t yq[1], xq[]={0.9};
235 h->GetQuantiles(1,yq,xq);
236 Store( h );
237
238 TH1F* htrunc = QuadraticDeviation(itgt, true, yq[0]);
239 TString name2 = TString::Format("%s_Quadr_Dev_best90perc_target_%d_",prefix.Data(),itgt);
240 htrunc->SetName( name2 );
241 htrunc->SetTitle( name2 );
242 Store( htrunc );
243 }
244 Log() << kINFO << "Results created" << Endl;
245}
#define h(i)
Definition RSha256.hxx:106
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Definition TGX11.cxx:110
float xmin
float ymin
float xmax
float ymax
static char * Format(const char *format, va_list ap)
Format a string in a circular formatting buffer (using a printf style format descriptor).
Definition TString.cxx:2420
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:577
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6707
void SetName(const char *name) override
Change the name of this histogram.
Definition TH1.cxx:8877
2-D histogram with a float per channel (see TH1 documentation)}
Definition TH2.h:258
Class that contains all the data information.
Definition DataSetInfo.h:62
UInt_t GetNVariables() const
UInt_t GetNTargets() const
VariableInfo & GetVariableInfo(Int_t i)
VariableInfo & GetTargetInfo(Int_t i)
Class that contains all the data information.
Definition DataSet.h:58
const Event * GetEvent() const
returns event without transformations
Definition DataSet.cxx:202
Long64_t GetNEvents(Types::ETreeType type=Types::kMaxTreeType) const
Definition DataSet.h:206
void SetCurrentType(Types::ETreeType type) const
Definition DataSet.h:89
Float_t GetValue(UInt_t ivar) const
return value of i'th variable
Definition Event.cxx:236
Double_t GetWeight() const
return the event weight - depending on whether the flag IgnoreNegWeightsInTraining is or not.
Definition Event.cxx:389
Float_t GetTarget(UInt_t itgt) const
Definition Event.h:102
ostringstream derivative to redirect and format output
Definition MsgLogger.h:57
void SetValue(std::vector< Float_t > &value, Int_t ievt)
TH2F * DeviationAsAFunctionOf(UInt_t varNum, UInt_t tgtNum)
TH1F * QuadraticDeviation(UInt_t tgtNum, Bool_t truncate=false, Double_t truncvalue=0.)
void CreateDeviationHistograms(TString prefix)
ResultsRegression(const DataSetInfo *dsi, TString resultsName)
constructor
Class that is the base-class for a vector of result.
Definition Results.h:57
Class for type info of MVA input variable.
Double_t GetMin() const
Double_t GetMax() const
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:380
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:2356
MsgLogger & Endl(MsgLogger &ml)
Definition MsgLogger.h:148
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
double epsilon
Definition triangle.c:618