31 #ifndef ROOT_TMVA_MethodBDT
32 #define ROOT_TMVA_MethodBDT
61 class MethodBDT :
public MethodBase {
66 MethodBDT(
const TString& jobName,
67 const TString& methodTitle,
69 const TString& theOption =
"");
72 MethodBDT( DataSetInfo& theData,
73 const TString& theWeightFile);
75 virtual ~MethodBDT(
void );
77 virtual Bool_t HasAnalysisType( Types::EAnalysisType type, UInt_t numberClasses, UInt_t numberTargets );
82 void InitEventSample();
85 virtual std::map<TString,Double_t> OptimizeTuningParameters(TString fomType=
"ROCIntegral", TString fitType=
"FitGA");
86 virtual void SetTuneParameters(std::map<TString,Double_t> tuneParameters);
94 using MethodBase::ReadWeightsFromStream;
97 void AddWeightsXMLTo(
void* parent )
const;
100 void ReadWeightsFromStream( std::istream& istr );
101 void ReadWeightsFromXML(
void* parent);
104 void WriteMonitoringHistosToFile(
void )
const;
107 Double_t GetMvaValue( Double_t* err = 0, Double_t* errUpper = 0);
110 UInt_t GetNTrees()
const {
return fForest.size();}
113 Double_t GetMvaValue( Double_t* err, Double_t* errUpper, UInt_t useNTrees );
114 Double_t PrivateGetMvaValue(
const TMVA::Event *ev, Double_t* err=0, Double_t* errUpper=0, UInt_t useNTrees=0 );
115 void BoostMonitor(Int_t iTree);
118 const std::vector<Float_t>& GetMulticlassValues();
121 const std::vector<Float_t>& GetRegressionValues();
124 Double_t Boost( std::vector<const TMVA::Event*>&, DecisionTree *dt, UInt_t cls = 0);
127 const Ranking* CreateRanking();
130 void DeclareOptions();
131 void ProcessOptions();
132 void SetMaxDepth(Int_t d){fMaxDepth = d;}
133 void SetMinNodeSize(Double_t sizeInPercent);
134 void SetMinNodeSize(TString sizeInPercent);
136 void SetNTrees(Int_t d){fNTrees = d;}
137 void SetAdaBoostBeta(Double_t b){fAdaBoostBeta = b;}
138 void SetNodePurityLimit(Double_t l){fNodePurityLimit = l;}
139 void SetShrinkage(Double_t s){fShrinkage = s;}
140 void SetUseNvars(Int_t n){fUseNvars = n;}
141 void SetBaggedSampleFraction(Double_t f){fBaggedSampleFraction = f;}
145 inline const std::vector<TMVA::DecisionTree*> & GetForest()
const;
148 inline const std::vector<const TMVA::Event*> & GetTrainingEvents()
const;
150 inline const std::vector<double> & GetBoostWeights()
const;
153 std::vector<Double_t> GetVariableImportance();
154 Double_t GetVariableImportance(UInt_t ivar);
156 Double_t TestTreeQuality( DecisionTree *dt );
159 void MakeClassSpecific( std::ostream&,
const TString& )
const;
162 void MakeClassSpecificHeader( std::ostream&,
const TString& )
const;
164 void MakeClassInstantiateNode( DecisionTreeNode *n, std::ostream& fout,
165 const TString& className )
const;
167 void GetHelpMessage()
const;
170 void DeclareCompatibilityOptions();
176 void PreProcessNegativeEventWeights();
179 Double_t AdaBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
182 Double_t AdaCost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
188 Double_t RegBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
191 Double_t AdaBoostR2( std::vector<const TMVA::Event*>&, DecisionTree *dt );
196 Double_t GradBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt, UInt_t cls = 0);
197 Double_t GradBoostRegression(std::vector<const TMVA::Event*>&, DecisionTree *dt );
198 void InitGradBoost( std::vector<const TMVA::Event*>&);
199 void UpdateTargets( std::vector<const TMVA::Event*>&, UInt_t cls = 0);
200 void UpdateTargetsRegression( std::vector<const TMVA::Event*>&,Bool_t first=kFALSE);
201 Double_t GetGradBoostMVA(
const TMVA::Event *e, UInt_t nTrees);
202 void GetBaggedSubSample(std::vector<const TMVA::Event*>&);
204 std::vector<const TMVA::Event*> fEventSample;
205 std::vector<const TMVA::Event*> fValidationSample;
206 std::vector<const TMVA::Event*> fSubSample;
207 std::vector<const TMVA::Event*> *fTrainSample;
210 std::vector<DecisionTree*> fForest;
211 std::vector<double> fBoostWeights;
212 Double_t fSigToBkgFraction;
214 Double_t fAdaBoostBeta;
215 TString fAdaBoostR2Loss;
219 Bool_t fBaggedGradBoost;
222 std::map< const TMVA::Event*, LossFunctionEventInfo> fLossFunctionEventInfo;
224 std::map< const TMVA::Event*,std::vector<double> > fResiduals;
227 SeparationBase *fSepType;
229 Int_t fMinNodeEvents;
230 Float_t fMinNodeSize;
231 TString fMinNodeSizeS;
234 Bool_t fUseFisherCuts;
235 Double_t fMinLinCorrForFisher;
236 Bool_t fUseExclusiveVars;
237 Bool_t fUseYesNoLeaf;
238 Double_t fNodePurityLimit;
242 DecisionTree::EPruneMethod fPruneMethod;
243 TString fPruneMethodS;
244 Double_t fPruneStrength;
245 Double_t fFValidationEvents;
247 Bool_t fRandomisedTrees;
249 Bool_t fUsePoissonNvars;
250 UInt_t fUseNTrainEvents;
252 Double_t fBaggedSampleFraction;
253 TString fNegWeightTreatment;
254 Bool_t fNoNegWeightsInTraining;
255 Bool_t fInverseBoostNegWeights;
256 Bool_t fPairNegWeightsGlobal;
257 Bool_t fTrainWithNegWeights;
258 Bool_t fDoBoostMonitor;
262 TTree* fMonitorNtuple;
264 Double_t fBoostWeight;
265 Double_t fErrorFraction;
272 Bool_t fDoPreselection;
274 Bool_t fSkipNormalization;
276 std::vector<Double_t> fVariableImportance;
279 void DeterminePreselectionCuts(
const std::vector<const TMVA::Event*>& eventSample);
280 Double_t ApplyPreselectionCuts(
const Event* ev);
282 std::vector<Double_t> fLowSigCut;
283 std::vector<Double_t> fLowBkgCut;
284 std::vector<Double_t> fHighSigCut;
285 std::vector<Double_t> fHighBkgCut;
287 std::vector<Bool_t> fIsLowSigCut;
288 std::vector<Bool_t> fIsLowBkgCut;
289 std::vector<Bool_t> fIsHighSigCut;
290 std::vector<Bool_t> fIsHighBkgCut;
292 Bool_t fHistoricBool;
294 TString fRegressionLossFunctionBDTGS;
295 Double_t fHuberQuantile;
297 LossFunctionBDT* fRegressionLossFunctionBDTG;
300 static const Int_t fgDebugLevel;
303 ClassDef(MethodBDT,0);
308 const std::vector<TMVA::DecisionTree*>& TMVA::MethodBDT::GetForest()
const {
return fForest; }
309 const std::vector<const TMVA::Event*> & TMVA::MethodBDT::GetTrainingEvents()
const {
return fEventSample; }
310 const std::vector<double>& TMVA::MethodBDT::GetBoostWeights()
const {
return fBoostWeights; }