Logo ROOT   6.08/07
Reference Guide
MethodBDT.h
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Andreas Hoecker, Joerg Stelzer, Helge Voss, Kai Voss, Jan Therhaag
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : MethodBDT (Boosted Decision Trees) *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Analysis of Boosted Decision Trees *
12  * *
13  * Authors (alphabetical): *
14  * Andreas Hoecker <Andreas.Hocker@cern.ch> - CERN, Switzerland *
15  * Helge Voss <Helge.Voss@cern.ch> - MPI-K Heidelberg, Germany *
16  * Kai Voss <Kai.Voss@cern.ch> - U. of Victoria, Canada *
17  * Doug Schouten <dschoute@sfu.ca> - Simon Fraser U., Canada *
18  * Jan Therhaag <jan.therhaag@cern.ch> - U. of Bonn, Germany *
19  * *
20  * Copyright (c) 2005-2011: *
21  * CERN, Switzerland *
22  * U. of Victoria, Canada *
23  * MPI-K Heidelberg, Germany *
24  * U. of Bonn, Germany *
25  * *
26  * Redistribution and use in source and binary forms, with or without *
27  * modification, are permitted according to the terms listed in LICENSE *
28  * (http://tmva.sourceforge.net/LICENSE) *
29  **********************************************************************************/
30 
31 #ifndef ROOT_TMVA_MethodBDT
32 #define ROOT_TMVA_MethodBDT
33 
34 //////////////////////////////////////////////////////////////////////////
35 // //
36 // MethodBDT //
37 // //
38 // Analysis of Boosted Decision Trees //
39 // //
40 //////////////////////////////////////////////////////////////////////////
41 
42 #include <vector>
43 #ifndef ROOT_TH2
44 #include "TH2.h"
45 #endif
46 #ifndef ROOT_TTree
47 #include "TTree.h"
48 #endif
49 #ifndef ROOT_TMVA_MethodBase
50 #include "TMVA/MethodBase.h"
51 #endif
52 #ifndef ROOT_TMVA_DecisionTree
53 #include "TMVA/DecisionTree.h"
54 #endif
55 #ifndef ROOT_TMVA_Event
56 #include "TMVA/Event.h"
57 #endif
58 #include "TMVA/LossFunction.h"
59 
60 namespace TMVA {
61 
62  class SeparationBase;
63 
64  class MethodBDT : public MethodBase {
65 
66  public:
67  // constructor for training and reading
68  MethodBDT( const TString& jobName,
69  const TString& methodTitle,
70  DataSetInfo& theData,
71  const TString& theOption = "");
72 
73  // constructor for calculating BDT-MVA using previously generatad decision trees
74  MethodBDT( DataSetInfo& theData,
75  const TString& theWeightFile);
76 
77  virtual ~MethodBDT( void );
78 
79  virtual Bool_t HasAnalysisType( Types::EAnalysisType type, UInt_t numberClasses, UInt_t numberTargets );
80 
81 
82  // write all Events from the Tree into a vector of Events, that are
83  // more easily manipulated
84  void InitEventSample();
85 
86  // optimize tuning parameters
87  virtual std::map<TString,Double_t> OptimizeTuningParameters(TString fomType="ROCIntegral", TString fitType="FitGA");
88  virtual void SetTuneParameters(std::map<TString,Double_t> tuneParameters);
89 
90  // training method
91  void Train( void );
92 
93  // revoke training
94  void Reset( void );
95 
97 
98  // write weights to file
99  void AddWeightsXMLTo( void* parent ) const;
100 
101  // read weights from file
102  void ReadWeightsFromStream( std::istream& istr );
103  void ReadWeightsFromXML(void* parent);
104 
105  // write method specific histos to target file
106  void WriteMonitoringHistosToFile( void ) const;
107 
108  // calculate the MVA value
109  Double_t GetMvaValue( Double_t* err = 0, Double_t* errUpper = 0);
110 
111  // get the actual forest size (might be less than fNTrees, the requested one, if boosting is stopped early
112  UInt_t GetNTrees() const {return fForest.size();}
113  private:
114  Double_t GetMvaValue( Double_t* err, Double_t* errUpper, UInt_t useNTrees );
115  Double_t PrivateGetMvaValue( const TMVA::Event *ev, Double_t* err=0, Double_t* errUpper=0, UInt_t useNTrees=0 );
116  void BoostMonitor(Int_t iTree);
117 
118  public:
119  const std::vector<Float_t>& GetMulticlassValues();
120 
121  // regression response
122  const std::vector<Float_t>& GetRegressionValues();
123 
124  // apply the boost algorithm to a tree in the collection
125  Double_t Boost( std::vector<const TMVA::Event*>&, DecisionTree *dt, UInt_t cls = 0);
126 
127  // ranking of input variables
128  const Ranking* CreateRanking();
129 
130  // the option handling methods
131  void DeclareOptions();
132  void ProcessOptions();
134  void SetMinNodeSize(Double_t sizeInPercent);
135  void SetMinNodeSize(TString sizeInPercent);
136 
137  void SetNTrees(Int_t d){fNTrees = d;}
143 
144 
145  // get the forest
146  inline const std::vector<TMVA::DecisionTree*> & GetForest() const;
147 
148  // get the forest
149  inline const std::vector<const TMVA::Event*> & GetTrainingEvents() const;
150 
151  inline const std::vector<double> & GetBoostWeights() const;
152 
153  //return the individual relative variable importance
154  std::vector<Double_t> GetVariableImportance();
156 
158 
159  // make ROOT-independent C++ class for classifier response (classifier-specific implementation)
160  void MakeClassSpecific( std::ostream&, const TString& ) const;
161 
162  // header and auxiliary classes
163  void MakeClassSpecificHeader( std::ostream&, const TString& ) const;
164 
165  void MakeClassInstantiateNode( DecisionTreeNode *n, std::ostream& fout,
166  const TString& className ) const;
167 
168  void GetHelpMessage() const;
169 
170  protected:
172 
173  private:
174  // Init used in the various constructors
175  void Init( void );
176 
178 
179  // boosting algorithm (adaptive boosting)
180  Double_t AdaBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
181 
182  // boosting algorithm (adaptive boosting with cost matrix)
183  Double_t AdaCost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
184 
185  // boosting as a random re-weighting
186  Double_t Bagging( );
187 
188  // boosting special for regression
189  Double_t RegBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
190 
191  // adaboost adapted to regression
192  Double_t AdaBoostR2( std::vector<const TMVA::Event*>&, DecisionTree *dt );
193 
194  // binomial likelihood gradient boost for classification
195  // (see Friedman: "Greedy Function Approximation: a Gradient Boosting Machine"
196  // Technical report, Dept. of Statistics, Stanford University)
197  Double_t GradBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt, UInt_t cls = 0);
198  Double_t GradBoostRegression(std::vector<const TMVA::Event*>&, DecisionTree *dt );
199  void InitGradBoost( std::vector<const TMVA::Event*>&);
200  void UpdateTargets( std::vector<const TMVA::Event*>&, UInt_t cls = 0);
201  void UpdateTargetsRegression( std::vector<const TMVA::Event*>&,Bool_t first=kFALSE);
202  Double_t GetGradBoostMVA(const TMVA::Event *e, UInt_t nTrees);
203  void GetBaggedSubSample(std::vector<const TMVA::Event*>&);
204 
205  std::vector<const TMVA::Event*> fEventSample; // the training events
206  std::vector<const TMVA::Event*> fValidationSample;// the Validation events
207  std::vector<const TMVA::Event*> fSubSample; // subsample for bagged grad boost
208  std::vector<const TMVA::Event*> *fTrainSample; // pointer to sample actually used in training (fEventSample or fSubSample) for example
209 
210  Int_t fNTrees; // number of decision trees requested
211  std::vector<DecisionTree*> fForest; // the collection of decision trees
212  std::vector<double> fBoostWeights; // the weights applied in the individual boosts
213  Double_t fSigToBkgFraction;// Signal to Background fraction assumed during training
214  TString fBoostType; // string specifying the boost type
215  Double_t fAdaBoostBeta; // beta parameter for AdaBoost algorithm
216  TString fAdaBoostR2Loss; // loss type used in AdaBoostR2 (Linear,Quadratic or Exponential)
217  //Double_t fTransitionPoint; // break-down point for gradient regression
218  Double_t fShrinkage; // learning rate for gradient boost;
219  Bool_t fBaggedBoost; // turn bagging in combination with boost on/off
220  Bool_t fBaggedGradBoost; // turn bagging in combination with grad boost on/off
221  //Double_t fSumOfWeights; // sum of all event weights
222  //std::map< const TMVA::Event*, std::pair<Double_t, Double_t> > fWeightedResiduals; // weighted regression residuals
223  std::map< const TMVA::Event*, LossFunctionEventInfo> fLossFunctionEventInfo; // map event to true value, predicted value, and weight
224  // used by different loss functions for BDT regression
225  std::map< const TMVA::Event*,std::vector<double> > fResiduals; // individual event residuals for gradient boost
226 
227  //options for the decision Tree
228  SeparationBase *fSepType; // the separation used in node splitting
229  TString fSepTypeS; // the separation (option string) used in node splitting
230  Int_t fMinNodeEvents; // min number of events in node
231  Float_t fMinNodeSize; // min percentage of training events in node
232  TString fMinNodeSizeS; // string containing min percentage of training events in node
233 
234  Int_t fNCuts; // grid used in cut applied in node splitting
235  Bool_t fUseFisherCuts; // use multivariate splits using the Fisher criterium
236  Double_t fMinLinCorrForFisher; // the minimum linear correlation between two variables demanded for use in fisher criterium in node splitting
237  Bool_t fUseExclusiveVars; // individual variables already used in fisher criterium are not anymore analysed individually for node splitting
238  Bool_t fUseYesNoLeaf; // use sig or bkg classification in leave nodes or sig/bkg
239  Double_t fNodePurityLimit; // purity limit for sig/bkg nodes
240  UInt_t fNNodesMax; // max # of nodes
241  UInt_t fMaxDepth; // max depth
242 
243  DecisionTree::EPruneMethod fPruneMethod; // method used for prunig
244  TString fPruneMethodS; // prune method option String
245  Double_t fPruneStrength; // a parameter to set the "amount" of pruning..needs to be adjusted
246  Double_t fFValidationEvents; // fraction of events to use for pruning
247  Bool_t fAutomatic; // use user given prune strength or automatically determined one using a validation sample
248  Bool_t fRandomisedTrees; // choose a random subset of possible cut variables at each node during training
249  UInt_t fUseNvars; // the number of variables used in the randomised tree splitting
250  Bool_t fUsePoissonNvars; // use "fUseNvars" not as fixed number but as mean of a possion distr. in each split
251  UInt_t fUseNTrainEvents; // number of randomly picked training events used in randomised (and bagged) trees
252 
253  Double_t fBaggedSampleFraction; // relative size of bagged event sample to original sample size
254  TString fNegWeightTreatment; // variable that holds the option of how to treat negative event weights in training
255  Bool_t fNoNegWeightsInTraining; // ignore negative event weights in the training
256  Bool_t fInverseBoostNegWeights; // boost ev. with neg. weights with 1/boostweight rathre than boostweight
257  Bool_t fPairNegWeightsGlobal; // pair ev. with neg. and pos. weights in traning sample and "annihilate" them
258  Bool_t fTrainWithNegWeights; // yes there are negative event weights and we don't ignore them
259  Bool_t fDoBoostMonitor; //create control plot with ROC integral vs tree number
260 
261 
262  //some histograms for monitoring
263  TTree* fMonitorNtuple; // monitoring ntuple
264  Int_t fITree; // ntuple var: ith tree
265  Double_t fBoostWeight; // ntuple var: boost weight
266  Double_t fErrorFraction; // ntuple var: misclassification error fraction
267 
268  Double_t fCss; // Cost factor
269  Double_t fCts_sb; // Cost factor
270  Double_t fCtb_ss; // Cost factor
271  Double_t fCbb; // Cost factor
272 
273  Bool_t fDoPreselection; // do or do not perform automatic pre-selection of 100% eff. cuts
274 
275  Bool_t fSkipNormalization; // true for skipping normalization at initialization of trees
276 
277  std::vector<Double_t> fVariableImportance; // the relative importance of the different variables
278 
279 
280  void DeterminePreselectionCuts(const std::vector<const TMVA::Event*>& eventSample);
282 
283  std::vector<Double_t> fLowSigCut;
284  std::vector<Double_t> fLowBkgCut;
285  std::vector<Double_t> fHighSigCut;
286  std::vector<Double_t> fHighBkgCut;
287 
288  std::vector<Bool_t> fIsLowSigCut;
289  std::vector<Bool_t> fIsLowBkgCut;
290  std::vector<Bool_t> fIsHighSigCut;
291  std::vector<Bool_t> fIsHighBkgCut;
292 
293  Bool_t fHistoricBool; //historic variable, only needed for "CompatibilityOptions"
294 
295  TString fRegressionLossFunctionBDTGS; // the option string determining the loss function for BDT regression
296  Double_t fHuberQuantile; // the option string determining the quantile for the Huber Loss Function
297  // in BDT regression.
299 
300  // debugging flags
301  static const Int_t fgDebugLevel; // debug level determining some printout/control plots etc.
302 
303  // for backward compatibility
304 
305  ClassDef(MethodBDT,0); // Analysis of Boosted Decision Trees
306  };
307 
308 } // namespace TMVA
309 
310 const std::vector<TMVA::DecisionTree*>& TMVA::MethodBDT::GetForest() const { return fForest; }
311 const std::vector<const TMVA::Event*> & TMVA::MethodBDT::GetTrainingEvents() const { return fEventSample; }
312 const std::vector<double>& TMVA::MethodBDT::GetBoostWeights() const { return fBoostWeights; }
313 
314 #endif
Bool_t fUseYesNoLeaf
Definition: MethodBDT.h:238
void Train(void)
BDT training.
Definition: MethodBDT.cxx:1134
void PreProcessNegativeEventWeights()
o.k.
Definition: MethodBDT.cxx:921
void GetBaggedSubSample(std::vector< const TMVA::Event *> &)
fills fEventSample with fBaggedSampleFraction*NEvents random training events
Definition: MethodBDT.cxx:2024
std::vector< Bool_t > fIsLowSigCut
Definition: MethodBDT.h:288
Double_t RegBoost(std::vector< const TMVA::Event *> &, DecisionTree *dt)
a special boosting only for Regression ...
Definition: MethodBDT.cxx:2059
void DeclareCompatibilityOptions()
options that are used ONLY for the READER to ensure backward compatibility
Definition: MethodBDT.cxx:444
std::map< const TMVA::Event *, LossFunctionEventInfo > fLossFunctionEventInfo
Definition: MethodBDT.h:223
Bool_t fPairNegWeightsGlobal
Definition: MethodBDT.h:257
void SetUseNvars(Int_t n)
Definition: MethodBDT.h:141
Bool_t fRandomisedTrees
Definition: MethodBDT.h:248
TString fMinNodeSizeS
Definition: MethodBDT.h:232
const Ranking * CreateRanking()
Compute ranking of input variables.
Definition: MethodBDT.cxx:2544
float Float_t
Definition: RtypesCore.h:53
TString fPruneMethodS
Definition: MethodBDT.h:244
TString fSepTypeS
Definition: MethodBDT.h:229
TTree * fMonitorNtuple
Definition: MethodBDT.h:263
Double_t fAdaBoostBeta
Definition: MethodBDT.h:215
std::vector< Bool_t > fIsHighSigCut
Definition: MethodBDT.h:290
void DeclareOptions()
define the options (their key words) that can be set in the option string know options: nTrees number...
Definition: MethodBDT.cxx:323
std::vector< Double_t > fVariableImportance
Definition: MethodBDT.h:277
void DeterminePreselectionCuts(const std::vector< const TMVA::Event *> &eventSample)
find useful preselection cuts that will be applied before and Decision Tree training.
Definition: MethodBDT.cxx:2854
EAnalysisType
Definition: Types.h:129
void MakeClassInstantiateNode(DecisionTreeNode *n, std::ostream &fout, const TString &className) const
recursively descends a tree and writes the node instance to the output streem
Definition: MethodBDT.cxx:2808
Double_t fMinLinCorrForFisher
Definition: MethodBDT.h:236
std::vector< const TMVA::Event * > fEventSample
Definition: MethodBDT.h:205
Double_t Bagging()
call it boot-strapping, re-sampling or whatever you like, in the end it is nothing else but applying ...
Definition: MethodBDT.cxx:2013
Bool_t fBaggedGradBoost
Definition: MethodBDT.h:220
Bool_t fDoBoostMonitor
Definition: MethodBDT.h:259
Basic string class.
Definition: TString.h:137
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
Double_t AdaBoost(std::vector< const TMVA::Event *> &, DecisionTree *dt)
the AdaBoost implementation.
Definition: MethodBDT.cxx:1714
void ProcessOptions()
the option string is decoded, for available options see "DeclareOptions"
Definition: MethodBDT.cxx:463
Bool_t fBaggedBoost
Definition: MethodBDT.h:219
void GetHelpMessage() const
Get help message text.
Definition: MethodBDT.cxx:2564
std::vector< Bool_t > fIsHighBkgCut
Definition: MethodBDT.h:291
void SetShrinkage(Double_t s)
Definition: MethodBDT.h:140
Double_t AdaCost(std::vector< const TMVA::Event *> &, DecisionTree *dt)
the AdaCost boosting algorithm takes a simple cost Matrix (currently fixed for all events...
Definition: MethodBDT.cxx:1895
Bool_t fAutomatic
Definition: MethodBDT.h:247
void MakeClassSpecific(std::ostream &, const TString &) const
make ROOT-independent C++ class for classifier response (classifier-specific implementation) ...
Definition: MethodBDT.cxx:2621
const std::vector< const TMVA::Event * > & GetTrainingEvents() const
Definition: MethodBDT.h:311
Double_t fCts_sb
Definition: MethodBDT.h:269
TString fRegressionLossFunctionBDTGS
Definition: MethodBDT.h:295
Double_t fBoostWeight
Definition: MethodBDT.h:265
Double_t GetMvaValue(Double_t *err=0, Double_t *errUpper=0)
Definition: MethodBDT.cxx:2317
Double_t fPruneStrength
Definition: MethodBDT.h:245
std::vector< Double_t > fHighBkgCut
Definition: MethodBDT.h:286
Double_t GetGradBoostMVA(const TMVA::Event *e, UInt_t nTrees)
returns MVA value: -1 for background, 1 for signal
Definition: MethodBDT.cxx:1410
Double_t fBaggedSampleFraction
Definition: MethodBDT.h:253
Bool_t fInverseBoostNegWeights
Definition: MethodBDT.h:256
Double_t GradBoostRegression(std::vector< const TMVA::Event *> &, DecisionTree *dt)
Implementation of M_TreeBoost using any loss function as desribed by Friedman 1999.
Definition: MethodBDT.cxx:1502
virtual void SetTuneParameters(std::map< TString, Double_t > tuneParameters)
set the tuning parameters accoding to the argument
Definition: MethodBDT.cxx:1112
#define ClassDef(name, id)
Definition: Rtypes.h:254
const std::vector< TMVA::DecisionTree * > & GetForest() const
Definition: MethodBDT.h:310
void MakeClassSpecificHeader(std::ostream &, const TString &) const
specific class header
Definition: MethodBDT.cxx:2697
Double_t fSigToBkgFraction
Definition: MethodBDT.h:213
virtual Bool_t HasAnalysisType(Types::EAnalysisType type, UInt_t numberClasses, UInt_t numberTargets)
BDT can handle classification with multiple classes and regression with one regression-target.
Definition: MethodBDT.cxx:275
Bool_t fDoPreselection
Definition: MethodBDT.h:273
void Reset(void)
reset the method, as if it had just been instantiated (forget all training etc.)
Definition: MethodBDT.cxx:718
void SetMinNodeSize(Double_t sizeInPercent)
Definition: MethodBDT.cxx:652
TString fAdaBoostR2Loss
Definition: MethodBDT.h:216
Int_t fMinNodeEvents
Definition: MethodBDT.h:230
Double_t AdaBoostR2(std::vector< const TMVA::Event *> &, DecisionTree *dt)
adaption of the AdaBoost to regression problems (see H.Drucker 1997)
Definition: MethodBDT.cxx:2067
void SetNTrees(Int_t d)
Definition: MethodBDT.h:137
std::vector< Double_t > fHighSigCut
Definition: MethodBDT.h:285
Double_t fCtb_ss
Definition: MethodBDT.h:270
const std::vector< Float_t > & GetMulticlassValues()
get the multiclass MVA response for the BDT classifier
Definition: MethodBDT.cxx:2368
Float_t fMinNodeSize
Definition: MethodBDT.h:231
Bool_t fNoNegWeightsInTraining
Definition: MethodBDT.h:255
const std::vector< Float_t > & GetRegressionValues()
get the regression value generated by the BDTs
Definition: MethodBDT.cxx:2403
void InitEventSample()
initialize the event sample (i.e. reset the boost-weights... etc)
Definition: MethodBDT.cxx:755
std::vector< Bool_t > fIsLowBkgCut
Definition: MethodBDT.h:289
void WriteMonitoringHistosToFile(void) const
Here we could write some histograms created during the processing to the output file.
Definition: MethodBDT.cxx:2489
Double_t fErrorFraction
Definition: MethodBDT.h:266
std::vector< Double_t > fLowBkgCut
Definition: MethodBDT.h:284
Double_t fNodePurityLimit
Definition: MethodBDT.h:239
Bool_t fHistoricBool
Definition: MethodBDT.h:293
void SetBaggedSampleFraction(Double_t f)
Definition: MethodBDT.h:142
void BoostMonitor(Int_t iTree)
fills the ROCIntegral vs Itree from the testSample for the monitoring plots during the training ...
Definition: MethodBDT.cxx:1620
Bool_t fTrainWithNegWeights
Definition: MethodBDT.h:258
Bool_t fSkipNormalization
Definition: MethodBDT.h:275
virtual ~MethodBDT(void)
destructor Note: fEventSample and ValidationSample are already deleted at the end of TRAIN When they ...
Definition: MethodBDT.cxx:747
void SetNodePurityLimit(Double_t l)
Definition: MethodBDT.h:139
Double_t PrivateGetMvaValue(const TMVA::Event *ev, Double_t *err=0, Double_t *errUpper=0, UInt_t useNTrees=0)
Return the MVA value (range [-1;1]) that classifies the event according to the majority vote from the...
Definition: MethodBDT.cxx:2341
unsigned int UInt_t
Definition: RtypesCore.h:42
Double_t GradBoost(std::vector< const TMVA::Event *> &, DecisionTree *dt, UInt_t cls=0)
Calculate the desired response value for each region.
Definition: MethodBDT.cxx:1471
TLine * l
Definition: textangle.C:4
Double_t fCbb
Definition: MethodBDT.h:271
SeparationBase * fSepType
Definition: MethodBDT.h:228
void Init(void)
common initialisation with defaults for the BDT-Method
Definition: MethodBDT.cxx:680
void ReadWeightsFromXML(void *parent)
reads the BDT from the xml file
Definition: MethodBDT.cxx:2215
Bool_t fUseExclusiveVars
Definition: MethodBDT.h:237
Double_t TestTreeQuality(DecisionTree *dt)
test the tree quality.. in terms of Miscalssification
Definition: MethodBDT.cxx:1565
DecisionTree::EPruneMethod fPruneMethod
Definition: MethodBDT.h:243
UInt_t fUseNvars
Definition: MethodBDT.h:249
UInt_t fMaxDepth
Definition: MethodBDT.h:241
double f(double x)
void ReadWeightsFromStream(std::istream &istr)
read the weights (BDT coefficients)
Definition: MethodBDT.cxx:2282
double Double_t
Definition: RtypesCore.h:55
Double_t ApplyPreselectionCuts(const Event *ev)
aply the preselection cuts before even bothing about any Decision Trees in the GetMVA ...
Definition: MethodBDT.cxx:2954
void UpdateTargets(std::vector< const TMVA::Event *> &, UInt_t cls=0)
Calculate residua for all events;.
Definition: MethodBDT.cxx:1424
void SetMaxDepth(Int_t d)
Definition: MethodBDT.h:133
int type
Definition: TGX11.cxx:120
void AddWeightsXMLTo(void *parent) const
write weights to XML
Definition: MethodBDT.cxx:2184
Double_t fShrinkage
Definition: MethodBDT.h:218
Bool_t fUsePoissonNvars
Definition: MethodBDT.h:250
void SetAdaBoostBeta(Double_t b)
Definition: MethodBDT.h:138
std::vector< const TMVA::Event * > * fTrainSample
Definition: MethodBDT.h:208
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
Double_t Boost(std::vector< const TMVA::Event *> &, DecisionTree *dt, UInt_t cls=0)
apply the boosting alogrithim (the algorithm is selecte via the the "option" given in the constructor...
Definition: MethodBDT.cxx:1586
LossFunctionBDT * fRegressionLossFunctionBDTG
Definition: MethodBDT.h:298
TString fBoostType
Definition: MethodBDT.h:214
Bool_t fUseFisherCuts
Definition: MethodBDT.h:235
UInt_t fUseNTrainEvents
Definition: MethodBDT.h:251
virtual std::map< TString, Double_t > OptimizeTuningParameters(TString fomType="ROCIntegral", TString fitType="FitGA")
call the Optimzier with the set of paremeters and ranges that are meant to be tuned.
Definition: MethodBDT.cxx:1059
TString fNegWeightTreatment
Definition: MethodBDT.h:254
Abstract ClassifierFactory template that handles arbitrary types.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
std::vector< const TMVA::Event * > fValidationSample
Definition: MethodBDT.h:206
std::vector< DecisionTree * > fForest
Definition: MethodBDT.h:211
std::vector< Double_t > GetVariableImportance()
Return the relative variable importance, normalized to all variables together having the importance 1...
Definition: MethodBDT.cxx:2504
Double_t fFValidationEvents
Definition: MethodBDT.h:246
const std::vector< double > & GetBoostWeights() const
Definition: MethodBDT.h:312
std::vector< Double_t > fLowSigCut
Definition: MethodBDT.h:283
void UpdateTargetsRegression(std::vector< const TMVA::Event *> &, Bool_t first=kFALSE)
Calculate current residuals for all events and update targets for next iteration. ...
Definition: MethodBDT.cxx:1457
A TTree object has a header with a name and a title.
Definition: TTree.h:98
std::map< const TMVA::Event *, std::vector< double > > fResiduals
Definition: MethodBDT.h:225
Definition: first.py:1
static const Int_t fgDebugLevel
Definition: MethodBDT.h:301
virtual void ReadWeightsFromStream(std::istream &)=0
Double_t fCss
Definition: MethodBDT.h:268
const Int_t n
Definition: legend1.C:16
std::vector< const TMVA::Event * > fSubSample
Definition: MethodBDT.h:207
UInt_t GetNTrees() const
Definition: MethodBDT.h:112
Double_t fHuberQuantile
Definition: MethodBDT.h:296
UInt_t fNNodesMax
Definition: MethodBDT.h:240
void InitGradBoost(std::vector< const TMVA::Event *> &)
initialize targets for first tree
Definition: MethodBDT.cxx:1527
std::vector< double > fBoostWeights
Definition: MethodBDT.h:212
MethodBDT(const TString &jobName, const TString &methodTitle, DataSetInfo &theData, const TString &theOption="")
the standard constructor for the "boosted decision trees"
Definition: MethodBDT.cxx:160