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