Logo ROOT   6.18/05
Reference Guide
RuleFitParams.h
Go to the documentation of this file.
1// @(#)root/tmva $Id$
2// Author: Andreas Hoecker, Joerg Stelzer, Fredrik Tegenfeldt, Helge Voss
3
4/**********************************************************************************
5 * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6 * Package: TMVA *
7 * Class : RuleFitParams *
8 * Web : http://tmva.sourceforge.net *
9 * *
10 * Description: *
11 * A class doing the actual fitting of a linear model using rules as *
12 * base functions. *
13 * Reference paper: 1.Gradient Directed Regularization *
14 * Friedman, Popescu, 2004 *
15 * 2.Predictive Learning with Rule Ensembles *
16 * Friedman, Popescu, 2005 *
17 * *
18 * *
19 * Authors (alphabetical): *
20 * Fredrik Tegenfeldt <Fredrik.Tegenfeldt@cern.ch> - Iowa State U., USA *
21 * Helge Voss <Helge.Voss@cern.ch> - MPI-KP Heidelberg, Ger. *
22 * *
23 * Copyright (c) 2005: *
24 * CERN, Switzerland *
25 * Iowa State U. *
26 * MPI-K Heidelberg, Germany *
27 * *
28 * Redistribution and use in source and binary forms, with or without *
29 * modification, are permitted according to the terms listed in LICENSE *
30 * (http://tmva.sourceforge.net/LICENSE) *
31 **********************************************************************************/
32
33#ifndef ROOT_TMVA_RuleFitParams
34#define ROOT_TMVA_RuleFitParams
35
36// #if ROOT_VERSION_CODE >= 364802
37#include "TMathBase.h"
38// #else
39// #ifndef ROOT_TMath
40// #include "TMath.h"
41// #endif
42// #endif
43
44#include "TMVA/Event.h"
45
46class TTree;
47
48namespace TMVA {
49
50 class RuleEnsemble;
51 class MsgLogger;
52 class RuleFit;
54
55 public:
56
58 virtual ~RuleFitParams();
59
60 void Init();
61
62 // set message type
63 void SetMsgType( EMsgType t );
64
65 // set RuleFit ptr
66 void SetRuleFit( RuleFit *rf ) { fRuleFit = rf; }
67 //
68 // GD path: set N(path steps)
69 void SetGDNPathSteps( Int_t np ) { fGDNPathSteps = np; }
70
71 // GD path: set path step size
73
74 // GD path: set tau search range
76 {
77 fGDTauMin = (t0>1.0 ? 1.0:(t0<0.0 ? 0.0:t0));
78 fGDTauMax = (t1>1.0 ? 1.0:(t1<0.0 ? 0.0:t1));
80 }
81
82 // GD path: set number of steps in tau search range
84
85 // GD path: set tau
86 void SetGDTau( Double_t t ) { fGDTau = t; }
87
88
91
92 // return type such that +1 = signal and -1 = background
93 Int_t Type( const Event * e ) const; // return (fRuleFit->GetMethodRuleFit()->DataInfo().IsSignal(e) ? 1:-1); }
94 //
95 UInt_t GetPathIdx1() const { return fPathIdx1; }
96 UInt_t GetPathIdx2() const { return fPathIdx2; }
97 UInt_t GetPerfIdx1() const { return fPerfIdx1; }
98 UInt_t GetPerfIdx2() const { return fPerfIdx2; }
99
100 // Loss function; Huber loss eq 33
101 Double_t LossFunction( const Event& e ) const;
102
103 // same but using evt idx (faster)
104 Double_t LossFunction( UInt_t evtidx ) const;
105 Double_t LossFunction( UInt_t evtidx, UInt_t itau ) const;
106
107 // Empirical risk
108 Double_t Risk(UInt_t ind1, UInt_t ind2, Double_t neff) const;
109 Double_t Risk(UInt_t ind1, UInt_t ind2, Double_t neff, UInt_t itau) const;
110
111 // Risk evaluation for fPathIdx and fPerfInd
115
116 // Risk evaluation for all tau
118
119 // Penalty function; Lasso function (eq 8)
120 Double_t Penalty() const;
121
122 // initialize GD path
123 void InitGD();
124
125 // find best tau and return the number of scan steps used
127
128 // make path for binary classification (squared-error ramp, sect 6 in ref 1)
129 void MakeGDPath();
130
131 protected:
132
133 // typedef of an Event const iterator
134 typedef std::vector<const TMVA::Event *>::const_iterator EventItr;
135
136 // init ntuple
137 void InitNtuple();
138
139 // calculate N(tau) in scan - limit to 100000.
140 void CalcGDNTau() { fGDNTau = static_cast<UInt_t>(1.0/fGDTauPrec)+1; if (fGDNTau>100000) fGDNTau=100000; }
141
142 // fill ntuple with coefficient info
143 void FillCoefficients();
144
145 // estimate the optimum scoring function
146 void CalcFStar();
147
148 // estimate of binary error rate
150
151 // estimate of scale average error rate
153
154 // estimate 1-area under ROC
155 Double_t ErrorRateRocRaw( std::vector<Double_t> & sFsig, std::vector<Double_t> & sFbkg );
157 void ErrorRateRocTst();
158
159 // estimate optimism
161
162 // make gradient vector (eq 44 in ref 1)
163 void MakeGradientVector();
164
165 // Calculate the direction in parameter space (eq 25, ref 1) and update coeffs (eq 22, ref 1)
166 void UpdateCoefficients();
167
168 // calculate average of responses of F
171
172 // calculate average of true response (initial estimate of a0)
174
175 // calculate the average of each variable over the range
176 void EvaluateAverage(UInt_t ind1, UInt_t ind2,
177 std::vector<Double_t> &avsel,
178 std::vector<Double_t> &avrul);
179
180 // evaluate using fPathIdx1,2
182
183 // evaluate using fPerfIdx1,2
185
186 // the same as above but for the various tau
190
191
192 RuleFit * fRuleFit; // rule fit
193 RuleEnsemble * fRuleEnsemble; // rule ensemble
194 //
195 UInt_t fNRules; // number of rules
196 UInt_t fNLinear; // number of linear terms
197 //
198 // Event indices for path/validation - TODO: should let the user decide
199 // Now it is just a simple one-fold cross validation.
200 //
201 UInt_t fPathIdx1; // first event index for path search
202 UInt_t fPathIdx2; // last event index for path search
203 UInt_t fPerfIdx1; // first event index for performance evaluation
204 UInt_t fPerfIdx2; // last event index for performance evaluation
205 Double_t fNEveEffPath; // sum of weights for Path events
206 Double_t fNEveEffPerf; // idem for Perf events
207
208 std::vector<Double_t> fAverageSelectorPath; // average of each variable over the range fPathIdx1,2
209 std::vector<Double_t> fAverageRulePath; // average of each rule, same range
210 std::vector<Double_t> fAverageSelectorPerf; // average of each variable over the range fPerfIdx1,2
211 std::vector<Double_t> fAverageRulePerf; // average of each rule, same range
212
213 std::vector<Double_t> fGradVec; // gradient vector - dimension = number of rules in ensemble
214 std::vector<Double_t> fGradVecLin; // gradient vector - dimension = number of variables
215
216 std::vector< std::vector<Double_t> > fGradVecTst; // gradient vector - one per tau
217 std::vector< std::vector<Double_t> > fGradVecLinTst; // gradient vector, linear terms - one per tau
218 //
219 std::vector<Double_t> fGDErrTst; // error rates per tau
220 std::vector<Char_t> fGDErrTstOK; // error rate is sufficiently low <--- stores boolean
221 std::vector< std::vector<Double_t> > fGDCoefTst; // rule coeffs - one per tau
222 std::vector< std::vector<Double_t> > fGDCoefLinTst; // linear coeffs - one per tau
223 std::vector<Double_t> fGDOfsTst; // offset per tau
224 std::vector< Double_t > fGDTauVec; // the tau's
225 UInt_t fGDNTauTstOK; // number of tau in the test-phase that are ok
226 UInt_t fGDNTau; // number of tau-paths - calculated in SetGDTauPrec
227 Double_t fGDTauPrec; // precision in tau
228 UInt_t fGDTauScan; // number scan for tau-paths
229 Double_t fGDTauMin; // min threshold parameter (tau in eq 26, ref 1)
230 Double_t fGDTauMax; // max threshold parameter (tau in eq 26, ref 1)
231 Double_t fGDTau; // selected threshold parameter (tau in eq 26, ref 1)
232 Double_t fGDPathStep; // step size along path (delta nu in eq 22, ref 1)
233 Int_t fGDNPathSteps; // number of path steps
234 Double_t fGDErrScale; // stop scan at error = scale*errmin
235 //
236 Double_t fAverageTruth; // average truth, ie sum(y)/N, y=+-1
237 //
238 std::vector<Double_t> fFstar; // vector of F*() - filled in CalcFStar()
239 Double_t fFstarMedian; // median value of F*() using
240 //
241 TTree *fGDNtuple; // Gradient path ntuple, contains params for each step along the path
242 Double_t fNTRisk; // GD path: risk
243 Double_t fNTErrorRate; // GD path: error rate (or performance)
244 Double_t fNTNuval; // GD path: value of nu
245 Double_t fNTCoefRad; // GD path: 'radius' of all rulecoeffs
246 Double_t fNTOffset; // GD path: model offset
247 Double_t *fNTCoeff; // GD path: rule coefficients
248 Double_t *fNTLinCoeff; // GD path: linear coefficients
249
250 Double_t fsigave; // Sigma of current signal score function F(sig)
251 Double_t fsigrms; // Rms of F(sig)
252 Double_t fbkgave; // Average of F(bkg)
253 Double_t fbkgrms; // Rms of F(bkg)
254
255 private:
256
257 mutable MsgLogger* fLogger; //! message logger
258 MsgLogger& Log() const { return *fLogger; }
259
260 };
261
262 // --------------------------------------------------------
263
264 class AbsValue {
265
266 public:
267
269 };
270}
271
272
273#endif
#define e(i)
Definition: RSha256.hxx:103
int Int_t
Definition: RtypesCore.h:41
unsigned int UInt_t
Definition: RtypesCore.h:42
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
Bool_t operator()(Double_t first, Double_t second) const
ostringstream derivative to redirect and format output
Definition: MsgLogger.h:59
A class doing the actual fitting of a linear model using rules as base functions.
Definition: RuleFitParams.h:53
void CalcTstAverageResponse()
calc average response for all test paths - TODO: see comment under CalcAverageResponse() note that 0 ...
void MakeGDPath()
The following finds the gradient directed path in parameter space.
std::vector< std::vector< Double_t > > fGDCoefLinTst
void EvaluateAverage(UInt_t ind1, UInt_t ind2, std::vector< Double_t > &avsel, std::vector< Double_t > &avrul)
evaluate the average of each variable and f(x) in the given range
Double_t * fNTLinCoeff
UInt_t RiskPerfTst()
Estimates the error rate with the current set of parameters.
Double_t Risk(UInt_t ind1, UInt_t ind2, Double_t neff) const
risk assessment
Double_t RiskPerf(UInt_t itau) const
Double_t Optimism()
implementation of eq.
Int_t FindGDTau()
This finds the cutoff parameter tau by scanning several different paths.
void SetGDPathStep(Double_t s)
Definition: RuleFitParams.h:72
std::vector< Double_t > fAverageSelectorPath
std::vector< Double_t > fGradVecLin
virtual ~RuleFitParams()
destructor
RuleFitParams()
constructor
void SetRuleFit(RuleFit *rf)
Definition: RuleFitParams.h:66
std::vector< constTMVA::Event * >::const_iterator EventItr
void Init()
Initializes all parameters using the RuleEnsemble and the training tree.
Double_t CalcAverageResponse()
calculate the average response - TODO : rewrite bad dependancy on EvaluateAverage() !
std::vector< Double_t > fAverageRulePerf
void SetMsgType(EMsgType t)
Double_t Penalty() const
This is the "lasso" penalty To be used for regression.
std::vector< std::vector< Double_t > > fGradVecLinTst
void SetGDTauRange(Double_t t0, Double_t t1)
Definition: RuleFitParams.h:75
std::vector< std::vector< Double_t > > fGDCoefTst
Double_t RiskPath() const
void FillCoefficients()
helper function to store the rule coefficients in local arrays
std::vector< Double_t > fGDTauVec
Double_t LossFunction(const Event &e) const
Implementation of squared-error ramp loss function (eq 39,40 in ref 1) This is used for binary Classi...
void SetGDErrScale(Double_t s)
Definition: RuleFitParams.h:89
void InitGD()
Initialize GD path search.
Int_t Type(const Event *e) const
std::vector< Double_t > fGDErrTst
void ErrorRateRocTst()
Estimates the error rate with the current set of parameters.
void SetGDTauPrec(Double_t p)
Definition: RuleFitParams.h:90
std::vector< Double_t > fGradVec
std::vector< Double_t > fAverageSelectorPerf
std::vector< std::vector< Double_t > > fGradVecTst
std::vector< Char_t > fGDErrTstOK
RuleEnsemble * fRuleEnsemble
void CalcFStar()
Estimates F* (optimum scoring function) for all events for the given sets.
MsgLogger & Log() const
message logger
UInt_t GetPathIdx2() const
Definition: RuleFitParams.h:96
void MakeTstGradientVector()
make test gradient vector for all tau same algorithm as MakeGradientVector()
UInt_t GetPerfIdx2() const
Definition: RuleFitParams.h:98
Double_t CalcAverageTruth()
calculate the average truth
void SetGDTauScan(UInt_t n)
Definition: RuleFitParams.h:83
std::vector< Double_t > fFstar
void UpdateTstCoefficients()
Establish maximum gradient for rules, linear terms and the offset for all taus TODO: do not need inde...
Double_t RiskPerf() const
void MakeGradientVector()
make gradient vector
void UpdateCoefficients()
Establish maximum gradient for rules, linear terms and the offset.
UInt_t GetPerfIdx1() const
Definition: RuleFitParams.h:97
void InitNtuple()
initializes the ntuple
Double_t CalcAverageResponseOLD()
Double_t ErrorRateBin()
Estimates the error rate with the current set of parameters It uses a binary estimate of (y-F*(x)) (y...
void SetGDTau(Double_t t)
Definition: RuleFitParams.h:86
Double_t ErrorRateReg()
Estimates the error rate with the current set of parameters This code is pretty messy at the moment.
std::vector< Double_t > fGDOfsTst
UInt_t GetPathIdx1() const
Definition: RuleFitParams.h:95
Double_t ErrorRateRoc()
Estimates the error rate with the current set of parameters.
Double_t ErrorRateRocRaw(std::vector< Double_t > &sFsig, std::vector< Double_t > &sFbkg)
Estimates the error rate with the current set of parameters.
std::vector< Double_t > fAverageRulePath
void SetGDNPathSteps(Int_t np)
Definition: RuleFitParams.h:69
A class implementing various fits of rule ensembles.
Definition: RuleFit.h:45
A TTree represents a columnar dataset.
Definition: TTree.h:71
const Int_t n
Definition: legend1.C:16
static constexpr double s
static constexpr double second
create variable transformations
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
Definition: first.py:1
auto * t1
Definition: textangle.C:20