Logo ROOT   6.10/09
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 #include "TH2.h"
44 #include "TTree.h"
45 #include "TMVA/MethodBase.h"
46 #include "TMVA/DecisionTree.h"
47 #include "TMVA/Event.h"
48 #include "TMVA/LossFunction.h"
49 
50 namespace TMVA {
51 
52  class SeparationBase;
53 
54  class MethodBDT : public MethodBase {
55 
56  public:
57  // constructor for training and reading
58  MethodBDT( const TString& jobName,
59  const TString& methodTitle,
60  DataSetInfo& theData,
61  const TString& theOption = "");
62 
63  // constructor for calculating BDT-MVA using previously generatad decision trees
64  MethodBDT( DataSetInfo& theData,
65  const TString& theWeightFile);
66 
67  virtual ~MethodBDT( void );
68 
69  virtual Bool_t HasAnalysisType( Types::EAnalysisType type, UInt_t numberClasses, UInt_t numberTargets );
70 
71 
72  // write all Events from the Tree into a vector of Events, that are
73  // more easily manipulated
74  void InitEventSample();
75 
76  // optimize tuning parameters
77  virtual std::map<TString,Double_t> OptimizeTuningParameters(TString fomType="ROCIntegral", TString fitType="FitGA");
78  virtual void SetTuneParameters(std::map<TString,Double_t> tuneParameters);
79 
80  // training method
81  void Train( void );
82 
83  // revoke training
84  void Reset( void );
85 
87 
88  // write weights to file
89  void AddWeightsXMLTo( void* parent ) const;
90 
91  // read weights from file
92  void ReadWeightsFromStream( std::istream& istr );
93  void ReadWeightsFromXML(void* parent);
94 
95  // write method specific histos to target file
96  void WriteMonitoringHistosToFile( void ) const;
97 
98  // calculate the MVA value
99  Double_t GetMvaValue( Double_t* err = 0, Double_t* errUpper = 0);
100 
101  // get the actual forest size (might be less than fNTrees, the requested one, if boosting is stopped early
102  UInt_t GetNTrees() const {return fForest.size();}
103  private:
104  Double_t GetMvaValue( Double_t* err, Double_t* errUpper, UInt_t useNTrees );
105  Double_t PrivateGetMvaValue( const TMVA::Event *ev, Double_t* err=0, Double_t* errUpper=0, UInt_t useNTrees=0 );
106  void BoostMonitor(Int_t iTree);
107 
108  public:
109  const std::vector<Float_t>& GetMulticlassValues();
110 
111  // regression response
112  const std::vector<Float_t>& GetRegressionValues();
113 
114  // apply the boost algorithm to a tree in the collection
115  Double_t Boost( std::vector<const TMVA::Event*>&, DecisionTree *dt, UInt_t cls = 0);
116 
117  // ranking of input variables
118  const Ranking* CreateRanking();
119 
120  // the option handling methods
121  void DeclareOptions();
122  void ProcessOptions();
124  void SetMinNodeSize(Double_t sizeInPercent);
125  void SetMinNodeSize(TString sizeInPercent);
126 
127  void SetNTrees(Int_t d){fNTrees = d;}
133 
134 
135  // get the forest
136  inline const std::vector<TMVA::DecisionTree*> & GetForest() const;
137 
138  // get the forest
139  inline const std::vector<const TMVA::Event*> & GetTrainingEvents() const;
140 
141  inline const std::vector<double> & GetBoostWeights() const;
142 
143  //return the individual relative variable importance
144  std::vector<Double_t> GetVariableImportance();
146 
148 
149  // make ROOT-independent C++ class for classifier response (classifier-specific implementation)
150  void MakeClassSpecific( std::ostream&, const TString& ) const;
151 
152  // header and auxiliary classes
153  void MakeClassSpecificHeader( std::ostream&, const TString& ) const;
154 
155  void MakeClassInstantiateNode( DecisionTreeNode *n, std::ostream& fout,
156  const TString& className ) const;
157 
158  void GetHelpMessage() const;
159 
160  protected:
162 
163  private:
164  // Init used in the various constructors
165  void Init( void );
166 
168 
169  // boosting algorithm (adaptive boosting)
170  Double_t AdaBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
171 
172  // boosting algorithm (adaptive boosting with cost matrix)
173  Double_t AdaCost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
174 
175  // boosting as a random re-weighting
176  Double_t Bagging( );
177 
178  // boosting special for regression
179  Double_t RegBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
180 
181  // adaboost adapted to regression
182  Double_t AdaBoostR2( std::vector<const TMVA::Event*>&, DecisionTree *dt );
183 
184  // binomial likelihood gradient boost for classification
185  // (see Friedman: "Greedy Function Approximation: a Gradient Boosting Machine"
186  // Technical report, Dept. of Statistics, Stanford University)
187  Double_t GradBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt, UInt_t cls = 0);
188  Double_t GradBoostRegression(std::vector<const TMVA::Event*>&, DecisionTree *dt );
189  void InitGradBoost( std::vector<const TMVA::Event*>&);
190  void UpdateTargets( std::vector<const TMVA::Event*>&, UInt_t cls = 0);
191  void UpdateTargetsRegression( std::vector<const TMVA::Event*>&,Bool_t first=kFALSE);
192  Double_t GetGradBoostMVA(const TMVA::Event *e, UInt_t nTrees);
193  void GetBaggedSubSample(std::vector<const TMVA::Event*>&);
194 
195  std::vector<const TMVA::Event*> fEventSample; // the training events
196  std::vector<const TMVA::Event*> fValidationSample;// the Validation events
197  std::vector<const TMVA::Event*> fSubSample; // subsample for bagged grad boost
198  std::vector<const TMVA::Event*> *fTrainSample; // pointer to sample actually used in training (fEventSample or fSubSample) for example
199 
200  Int_t fNTrees; // number of decision trees requested
201  std::vector<DecisionTree*> fForest; // the collection of decision trees
202  std::vector<double> fBoostWeights; // the weights applied in the individual boosts
203  Double_t fSigToBkgFraction;// Signal to Background fraction assumed during training
204  TString fBoostType; // string specifying the boost type
205  Double_t fAdaBoostBeta; // beta parameter for AdaBoost algorithm
206  TString fAdaBoostR2Loss; // loss type used in AdaBoostR2 (Linear,Quadratic or Exponential)
207  //Double_t fTransitionPoint; // break-down point for gradient regression
208  Double_t fShrinkage; // learning rate for gradient boost;
209  Bool_t fBaggedBoost; // turn bagging in combination with boost on/off
210  Bool_t fBaggedGradBoost; // turn bagging in combination with grad boost on/off
211  //Double_t fSumOfWeights; // sum of all event weights
212  //std::map< const TMVA::Event*, std::pair<Double_t, Double_t> > fWeightedResiduals; // weighted regression residuals
213  std::map< const TMVA::Event*, LossFunctionEventInfo> fLossFunctionEventInfo; // map event to true value, predicted value, and weight
214  // used by different loss functions for BDT regression
215  std::map< const TMVA::Event*,std::vector<double> > fResiduals; // individual event residuals for gradient boost
216 
217  //options for the decision Tree
218  SeparationBase *fSepType; // the separation used in node splitting
219  TString fSepTypeS; // the separation (option string) used in node splitting
220  Int_t fMinNodeEvents; // min number of events in node
221  Float_t fMinNodeSize; // min percentage of training events in node
222  TString fMinNodeSizeS; // string containing min percentage of training events in node
223 
224  Int_t fNCuts; // grid used in cut applied in node splitting
225  Bool_t fUseFisherCuts; // use multivariate splits using the Fisher criterium
226  Double_t fMinLinCorrForFisher; // the minimum linear correlation between two variables demanded for use in fisher criterium in node splitting
227  Bool_t fUseExclusiveVars; // individual variables already used in fisher criterium are not anymore analysed individually for node splitting
228  Bool_t fUseYesNoLeaf; // use sig or bkg classification in leave nodes or sig/bkg
229  Double_t fNodePurityLimit; // purity limit for sig/bkg nodes
230  UInt_t fNNodesMax; // max # of nodes
231  UInt_t fMaxDepth; // max depth
232 
233  DecisionTree::EPruneMethod fPruneMethod; // method used for prunig
234  TString fPruneMethodS; // prune method option String
235  Double_t fPruneStrength; // a parameter to set the "amount" of pruning..needs to be adjusted
236  Double_t fFValidationEvents; // fraction of events to use for pruning
237  Bool_t fAutomatic; // use user given prune strength or automatically determined one using a validation sample
238  Bool_t fRandomisedTrees; // choose a random subset of possible cut variables at each node during training
239  UInt_t fUseNvars; // the number of variables used in the randomised tree splitting
240  Bool_t fUsePoissonNvars; // use "fUseNvars" not as fixed number but as mean of a possion distr. in each split
241  UInt_t fUseNTrainEvents; // number of randomly picked training events used in randomised (and bagged) trees
242 
243  Double_t fBaggedSampleFraction; // relative size of bagged event sample to original sample size
244  TString fNegWeightTreatment; // variable that holds the option of how to treat negative event weights in training
245  Bool_t fNoNegWeightsInTraining; // ignore negative event weights in the training
246  Bool_t fInverseBoostNegWeights; // boost ev. with neg. weights with 1/boostweight rathre than boostweight
247  Bool_t fPairNegWeightsGlobal; // pair ev. with neg. and pos. weights in traning sample and "annihilate" them
248  Bool_t fTrainWithNegWeights; // yes there are negative event weights and we don't ignore them
249  Bool_t fDoBoostMonitor; //create control plot with ROC integral vs tree number
250 
251 
252  //some histograms for monitoring
253  TTree* fMonitorNtuple; // monitoring ntuple
254  Int_t fITree; // ntuple var: ith tree
255  Double_t fBoostWeight; // ntuple var: boost weight
256  Double_t fErrorFraction; // ntuple var: misclassification error fraction
257 
258  Double_t fCss; // Cost factor
259  Double_t fCts_sb; // Cost factor
260  Double_t fCtb_ss; // Cost factor
261  Double_t fCbb; // Cost factor
262 
263  Bool_t fDoPreselection; // do or do not perform automatic pre-selection of 100% eff. cuts
264 
265  Bool_t fSkipNormalization; // true for skipping normalization at initialization of trees
266 
267  std::vector<Double_t> fVariableImportance; // the relative importance of the different variables
268 
269 
270  void DeterminePreselectionCuts(const std::vector<const TMVA::Event*>& eventSample);
272 
273  std::vector<Double_t> fLowSigCut;
274  std::vector<Double_t> fLowBkgCut;
275  std::vector<Double_t> fHighSigCut;
276  std::vector<Double_t> fHighBkgCut;
277 
278  std::vector<Bool_t> fIsLowSigCut;
279  std::vector<Bool_t> fIsLowBkgCut;
280  std::vector<Bool_t> fIsHighSigCut;
281  std::vector<Bool_t> fIsHighBkgCut;
282 
283  Bool_t fHistoricBool; //historic variable, only needed for "CompatibilityOptions"
284 
285  TString fRegressionLossFunctionBDTGS; // the option string determining the loss function for BDT regression
286  Double_t fHuberQuantile; // the option string determining the quantile for the Huber Loss Function
287  // in BDT regression.
289 
290  // debugging flags
291  static const Int_t fgDebugLevel; // debug level determining some printout/control plots etc.
292 
293  // for backward compatibility
294 
295  ClassDef(MethodBDT,0); // Analysis of Boosted Decision Trees
296  };
297 
298 } // namespace TMVA
299 
300 const std::vector<TMVA::DecisionTree*>& TMVA::MethodBDT::GetForest() const { return fForest; }
301 const std::vector<const TMVA::Event*> & TMVA::MethodBDT::GetTrainingEvents() const { return fEventSample; }
302 const std::vector<double>& TMVA::MethodBDT::GetBoostWeights() const { return fBoostWeights; }
303 
304 #endif
Bool_t fUseYesNoLeaf
Definition: MethodBDT.h:228
void Train(void)
BDT training.
Definition: MethodBDT.cxx:1135
void PreProcessNegativeEventWeights()
O.k.
Definition: MethodBDT.cxx:926
void GetBaggedSubSample(std::vector< const TMVA::Event *> &)
Fills fEventSample with fBaggedSampleFraction*NEvents random training events.
Definition: MethodBDT.cxx:2031
std::vector< Bool_t > fIsLowSigCut
Definition: MethodBDT.h:278
Double_t RegBoost(std::vector< const TMVA::Event *> &, DecisionTree *dt)
A special boosting only for Regression (not implemented).
Definition: MethodBDT.cxx:2065
void DeclareCompatibilityOptions()
Options that are used ONLY for the READER to ensure backward compatibility.
Definition: MethodBDT.cxx:455
std::map< const TMVA::Event *, LossFunctionEventInfo > fLossFunctionEventInfo
Definition: MethodBDT.h:213
Bool_t fPairNegWeightsGlobal
Definition: MethodBDT.h:247
void SetUseNvars(Int_t n)
Definition: MethodBDT.h:131
Bool_t fRandomisedTrees
Definition: MethodBDT.h:238
TString fMinNodeSizeS
Definition: MethodBDT.h:222
const Ranking * CreateRanking()
Compute ranking of input variables.
Definition: MethodBDT.cxx:2552
float Float_t
Definition: RtypesCore.h:53
TString fPruneMethodS
Definition: MethodBDT.h:234
TString fSepTypeS
Definition: MethodBDT.h:219
TTree * fMonitorNtuple
Definition: MethodBDT.h:253
Double_t fAdaBoostBeta
Definition: MethodBDT.h:205
std::vector< Bool_t > fIsHighSigCut
Definition: MethodBDT.h:280
void DeclareOptions()
Define the options (their key words).
Definition: MethodBDT.cxx:334
std::vector< Double_t > fVariableImportance
Definition: MethodBDT.h:267
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:2858
EAnalysisType
Definition: Types.h:125
void MakeClassInstantiateNode(DecisionTreeNode *n, std::ostream &fout, const TString &className) const
Recursively descends a tree and writes the node instance to the output stream.
Definition: MethodBDT.cxx:2813
Double_t fMinLinCorrForFisher
Definition: MethodBDT.h:226
Virtual base Class for all MVA method.
Definition: MethodBase.h:106
std::vector< const TMVA::Event * > fEventSample
Definition: MethodBDT.h:195
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:2020
Bool_t fBaggedGradBoost
Definition: MethodBDT.h:210
Bool_t fDoBoostMonitor
Definition: MethodBDT.h:249
Basic string class.
Definition: TString.h:129
Ranking for variables in method (implementation)
Definition: Ranking.h:48
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Double_t AdaBoost(std::vector< const TMVA::Event *> &, DecisionTree *dt)
The AdaBoost implementation.
Definition: MethodBDT.cxx:1726
void ProcessOptions()
The option string is decoded, for available options see "DeclareOptions".
Definition: MethodBDT.cxx:471
Bool_t fBaggedBoost
Definition: MethodBDT.h:209
void GetHelpMessage() const
Get help message text.
Definition: MethodBDT.cxx:2569
std::vector< Bool_t > fIsHighBkgCut
Definition: MethodBDT.h:281
void SetShrinkage(Double_t s)
Definition: MethodBDT.h:130
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:1904
Bool_t fAutomatic
Definition: MethodBDT.h:237
void MakeClassSpecific(std::ostream &, const TString &) const
Make ROOT-independent C++ class for classifier response (classifier-specific implementation).
Definition: MethodBDT.cxx:2626
const std::vector< const TMVA::Event * > & GetTrainingEvents() const
Definition: MethodBDT.h:301
Double_t fCts_sb
Definition: MethodBDT.h:259
TString fRegressionLossFunctionBDTGS
Definition: MethodBDT.h:285
Double_t fBoostWeight
Definition: MethodBDT.h:255
Double_t GetMvaValue(Double_t *err=0, Double_t *errUpper=0)
Definition: MethodBDT.cxx:2323
Double_t fPruneStrength
Definition: MethodBDT.h:235
std::vector< Double_t > fHighBkgCut
Definition: MethodBDT.h:276
Double_t GetGradBoostMVA(const TMVA::Event *e, UInt_t nTrees)
Returns MVA value: -1 for background, 1 for signal.
Definition: MethodBDT.cxx:1411
Double_t fBaggedSampleFraction
Definition: MethodBDT.h:243
Bool_t fInverseBoostNegWeights
Definition: MethodBDT.h:246
Double_t GradBoostRegression(std::vector< const TMVA::Event *> &, DecisionTree *dt)
Implementation of M_TreeBoost using any loss function as described by Friedman 1999.
Definition: MethodBDT.cxx:1514
virtual void SetTuneParameters(std::map< TString, Double_t > tuneParameters)
Set the tuning parameters according to the argument.
Definition: MethodBDT.cxx:1115
#define ClassDef(name, id)
Definition: Rtypes.h:297
const std::vector< TMVA::DecisionTree * > & GetForest() const
Definition: MethodBDT.h:300
void MakeClassSpecificHeader(std::ostream &, const TString &) const
Specific class header.
Definition: MethodBDT.cxx:2702
Double_t fSigToBkgFraction
Definition: MethodBDT.h:203
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:281
Bool_t fDoPreselection
Definition: MethodBDT.h:263
void Reset(void)
Reset the method, as if it had just been instantiated (forget all training etc.). ...
Definition: MethodBDT.cxx:724
void SetMinNodeSize(Double_t sizeInPercent)
Definition: MethodBDT.cxx:659
TString fAdaBoostR2Loss
Definition: MethodBDT.h:206
Int_t fMinNodeEvents
Definition: MethodBDT.h:220
Double_t AdaBoostR2(std::vector< const TMVA::Event *> &, DecisionTree *dt)
Adaption of the AdaBoost to regression problems (see H.Drucker 1997).
Definition: MethodBDT.cxx:2073
void SetNTrees(Int_t d)
Definition: MethodBDT.h:127
std::vector< Double_t > fHighSigCut
Definition: MethodBDT.h:275
Class that contains all the data information.
Definition: DataSetInfo.h:60
Double_t fCtb_ss
Definition: MethodBDT.h:260
const std::vector< Float_t > & GetMulticlassValues()
Get the multiclass MVA response for the BDT classifier.
Definition: MethodBDT.cxx:2375
Float_t fMinNodeSize
Definition: MethodBDT.h:221
Bool_t fNoNegWeightsInTraining
Definition: MethodBDT.h:245
const std::vector< Float_t > & GetRegressionValues()
Get the regression value generated by the BDTs.
Definition: MethodBDT.cxx:2411
void InitEventSample()
Initialize the event sample (i.e. reset the boost-weights... etc).
Definition: MethodBDT.cxx:760
std::vector< Bool_t > fIsLowBkgCut
Definition: MethodBDT.h:279
void WriteMonitoringHistosToFile(void) const
Here we could write some histograms created during the processing to the output file.
Definition: MethodBDT.cxx:2497
Double_t fErrorFraction
Definition: MethodBDT.h:256
std::vector< Double_t > fLowBkgCut
Definition: MethodBDT.h:274
Double_t fNodePurityLimit
Definition: MethodBDT.h:229
Bool_t fHistoricBool
Definition: MethodBDT.h:283
void SetBaggedSampleFraction(Double_t f)
Definition: MethodBDT.h:132
void BoostMonitor(Int_t iTree)
Fills the ROCIntegral vs Itree from the testSample for the monitoring plots during the training ...
Definition: MethodBDT.cxx:1632
Bool_t fTrainWithNegWeights
Definition: MethodBDT.h:248
Bool_t fSkipNormalization
Definition: MethodBDT.h:265
virtual ~MethodBDT(void)
Destructor.
Definition: MethodBDT.cxx:752
void SetNodePurityLimit(Double_t l)
Definition: MethodBDT.h:129
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:2348
Implementation of a Decision Tree.
Definition: DecisionTree.h:59
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:1481
TLine * l
Definition: textangle.C:4
Double_t fCbb
Definition: MethodBDT.h:261
SeparationBase * fSepType
Definition: MethodBDT.h:218
void Init(void)
Common initialisation with defaults for the BDT-Method.
Definition: MethodBDT.cxx:686
void ReadWeightsFromXML(void *parent)
Reads the BDT from the xml file.
Definition: MethodBDT.cxx:2221
An interface to calculate the "SeparationGain" for different separation criteria used in various trai...
Bool_t fUseExclusiveVars
Definition: MethodBDT.h:227
Double_t TestTreeQuality(DecisionTree *dt)
Test the tree quality.. in terms of Misclassification.
Definition: MethodBDT.cxx:1577
const Bool_t kFALSE
Definition: RtypesCore.h:92
DecisionTree::EPruneMethod fPruneMethod
Definition: MethodBDT.h:233
UInt_t fUseNvars
Definition: MethodBDT.h:239
UInt_t fMaxDepth
Definition: MethodBDT.h:231
double f(double x)
void ReadWeightsFromStream(std::istream &istr)
Read the weights (BDT coefficients).
Definition: MethodBDT.cxx:2288
double Double_t
Definition: RtypesCore.h:55
Double_t ApplyPreselectionCuts(const Event *ev)
Apply the preselection cuts before even bothering about any Decision Trees in the GetMVA ...
Definition: MethodBDT.cxx:2958
void UpdateTargets(std::vector< const TMVA::Event *> &, UInt_t cls=0)
Calculate residual for all events.
Definition: MethodBDT.cxx:1425
void SetMaxDepth(Int_t d)
Definition: MethodBDT.h:123
int type
Definition: TGX11.cxx:120
void AddWeightsXMLTo(void *parent) const
Write weights to XML.
Definition: MethodBDT.cxx:2190
Double_t fShrinkage
Definition: MethodBDT.h:208
Bool_t fUsePoissonNvars
Definition: MethodBDT.h:240
void SetAdaBoostBeta(Double_t b)
Definition: MethodBDT.h:128
std::vector< const TMVA::Event * > * fTrainSample
Definition: MethodBDT.h:198
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 algorithm (the algorithm is selecte via the the "option" given in the constructor...
Definition: MethodBDT.cxx:1598
LossFunctionBDT * fRegressionLossFunctionBDTG
Definition: MethodBDT.h:288
TString fBoostType
Definition: MethodBDT.h:204
Bool_t fUseFisherCuts
Definition: MethodBDT.h:225
UInt_t fUseNTrainEvents
Definition: MethodBDT.h:241
virtual std::map< TString, Double_t > OptimizeTuningParameters(TString fomType="ROCIntegral", TString fitType="FitGA")
Call the Optimizer with the set of parameters and ranges that are meant to be tuned.
Definition: MethodBDT.cxx:1062
TString fNegWeightTreatment
Definition: MethodBDT.h:244
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:196
std::vector< DecisionTree * > fForest
Definition: MethodBDT.h:201
std::vector< Double_t > GetVariableImportance()
Return the relative variable importance, normalized to all variables together having the importance 1...
Definition: MethodBDT.cxx:2512
Double_t fFValidationEvents
Definition: MethodBDT.h:236
const std::vector< double > & GetBoostWeights() const
Definition: MethodBDT.h:302
std::vector< Double_t > fLowSigCut
Definition: MethodBDT.h:273
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:1467
A TTree object has a header with a name and a title.
Definition: TTree.h:78
std::map< const TMVA::Event *, std::vector< double > > fResiduals
Definition: MethodBDT.h:215
Definition: first.py:1
static const Int_t fgDebugLevel
Definition: MethodBDT.h:291
virtual void ReadWeightsFromStream(std::istream &)=0
Double_t fCss
Definition: MethodBDT.h:258
const Int_t n
Definition: legend1.C:16
std::vector< const TMVA::Event * > fSubSample
Definition: MethodBDT.h:197
UInt_t GetNTrees() const
Definition: MethodBDT.h:102
Analysis of Boosted Decision Trees.
Definition: MethodBDT.h:54
Double_t fHuberQuantile
Definition: MethodBDT.h:286
UInt_t fNNodesMax
Definition: MethodBDT.h:230
void InitGradBoost(std::vector< const TMVA::Event *> &)
Initialize targets for first tree.
Definition: MethodBDT.cxx:1539
std::vector< double > fBoostWeights
Definition: MethodBDT.h:202
MethodBDT(const TString &jobName, const TString &methodTitle, DataSetInfo &theData, const TString &theOption="")
The standard constructor for the "boosted decision trees".
Definition: MethodBDT.cxx:164