ROOT  6.06/09
Reference Guide
MethodBoost.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Andreas Hoecker, Joerg Stelzer, Helge Voss, Kai Voss,Or Cohen, Jan Therhaag, Eckhard von Toerne
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : MethodCompositeBase *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Virtual base class for all MVA method *
12  * *
13  * Authors (alphabetical): *
14  * Andreas Hoecker <Andreas.Hocker@cern.ch> - CERN, Switzerland *
15  * Peter Speckmayer <Peter.Speckmazer@cern.ch> - CERN, Switzerland *
16  * Joerg Stelzer <Joerg.Stelzer@cern.ch> - CERN, Switzerland *
17  * Helge Voss <Helge.Voss@cern.ch> - MPI-K Heidelberg, Germany *
18  * Jan Therhaag <Jan.Therhaag@cern.ch> - U of Bonn, Germany *
19  * Eckhard v. Toerne <evt@uni-bonn.de> - U of Bonn, Germany *
20  * *
21  * Copyright (c) 2005-2011: *
22  * CERN, Switzerland *
23  * U. of Victoria, Canada *
24  * MPI-K Heidelberg, Germany *
25  * U. of Bonn, Germany *
26  * *
27  * Redistribution and use in source and binary forms, with or without *
28  * modification, are permitted according to the terms listed in LICENSE *
29  * (http://tmva.sourceforge.net/LICENSE) *
30  **********************************************************************************/
31 
32 //_______________________________________________________________________
33 //
34 // This class is meant to boost a single classifier. Boosting means //
35 // training the classifier a few times. Everytime the wieghts of the //
36 // events are modified according to how well the classifier performed //
37 // on the test sample. //
38 ////////////////////////////////////////////////////////////////////////////////
39 
40 #include <algorithm>
41 #include <iomanip>
42 #include <vector>
43 #include <cmath>
44 
45 #include "Riostream.h"
46 #include "TRandom3.h"
47 #include "TMath.h"
48 #include "TObjString.h"
49 #include "TH1F.h"
50 #include "TGraph.h"
51 #include "TSpline.h"
52 #include "TDirectory.h"
53 
55 #include "TMVA/MethodBase.h"
56 #include "TMVA/MethodBoost.h"
57 #include "TMVA/MethodCategory.h"
58 #include "TMVA/MethodDT.h"
59 #include "TMVA/MethodFisher.h"
60 #include "TMVA/Tools.h"
61 #include "TMVA/ClassifierFactory.h"
62 #include "TMVA/Timer.h"
63 #include "TMVA/Types.h"
64 #include "TMVA/PDF.h"
65 #include "TMVA/Results.h"
66 #include "TMVA/Config.h"
67 
68 #include "TMVA/SeparationBase.h"
70 #include "TMVA/GiniIndex.h"
71 #include "TMVA/CrossEntropy.h"
74 
75 REGISTER_METHOD(Boost)
76 
77 ClassImp(TMVA::MethodBoost)
78 
79 ////////////////////////////////////////////////////////////////////////////////
80 
81 TMVA::MethodBoost::MethodBoost( const TString& jobName,
82  const TString& methodTitle,
83  DataSetInfo& theData,
84  const TString& theOption,
85  TDirectory* theTargetDir ) :
86  TMVA::MethodCompositeBase( jobName, Types::kBoost, methodTitle, theData, theOption, theTargetDir )
87  , fBoostNum(0)
88  , fDetailedMonitoring(kFALSE)
89  , fAdaBoostBeta(0)
90  , fRandomSeed(0)
91  , fBaggedSampleFraction(0)
92  , fBoostedMethodTitle(methodTitle)
93  , fBoostedMethodOptions(theOption)
94  , fMonitorBoostedMethod(kFALSE)
95  , fMonitorTree(0)
96  , fBoostWeight(0)
97  , fMethodError(0)
98  , fROC_training(0.0)
99  , fOverlap_integral(0.0)
100  , fMVAvalues(0)
101 {
102  fMVAvalues = new std::vector<Float_t>;
103 }
104 
105 ////////////////////////////////////////////////////////////////////////////////
106 
108  const TString& theWeightFile,
109  TDirectory* theTargetDir )
110  : TMVA::MethodCompositeBase( Types::kBoost, dsi, theWeightFile, theTargetDir )
111  , fBoostNum(0)
112  , fDetailedMonitoring(kFALSE)
113  , fAdaBoostBeta(0)
114  , fRandomSeed(0)
115  , fBaggedSampleFraction(0)
116  , fBoostedMethodTitle("")
117  , fBoostedMethodOptions("")
118  , fMonitorBoostedMethod(kFALSE)
119  , fMonitorTree(0)
120  , fBoostWeight(0)
121  , fMethodError(0)
122  , fROC_training(0.0)
123  , fOverlap_integral(0.0)
124  , fMVAvalues(0)
125 {
126  fMVAvalues = new std::vector<Float_t>;
127 }
128 
129 ////////////////////////////////////////////////////////////////////////////////
130 /// destructor
131 
133 {
134  fMethodWeight.clear();
135 
136  // the histogram themselves are deleted when the file is closed
137 
138  fTrainSigMVAHist.clear();
139  fTrainBgdMVAHist.clear();
140  fBTrainSigMVAHist.clear();
141  fBTrainBgdMVAHist.clear();
142  fTestSigMVAHist.clear();
143  fTestBgdMVAHist.clear();
144 
145  if (fMVAvalues) {
146  delete fMVAvalues;
147  fMVAvalues = 0;
148  }
149 }
150 
151 
152 ////////////////////////////////////////////////////////////////////////////////
153 /// Boost can handle classification with 2 classes and regression with one regression-target
154 
156 {
157  if (type == Types::kClassification && numberClasses == 2) return kTRUE;
158  // if (type == Types::kRegression && numberTargets == 1) return kTRUE;
159  return kFALSE;
160 }
161 
162 
163 ////////////////////////////////////////////////////////////////////////////////
164 
166 {
167  DeclareOptionRef( fBoostNum = 1, "Boost_Num",
168  "Number of times the classifier is boosted" );
169 
170  DeclareOptionRef( fMonitorBoostedMethod = kTRUE, "Boost_MonitorMethod",
171  "Write monitoring histograms for each boosted classifier" );
172 
173  DeclareOptionRef( fDetailedMonitoring = kFALSE, "Boost_DetailedMonitoring",
174  "Produce histograms for detailed boost monitoring" );
175 
176  DeclareOptionRef( fBoostType = "AdaBoost", "Boost_Type", "Boosting type for the classifiers" );
177  AddPreDefVal(TString("RealAdaBoost"));
178  AddPreDefVal(TString("AdaBoost"));
179  AddPreDefVal(TString("Bagging"));
180 
181  DeclareOptionRef(fBaggedSampleFraction=.6,"Boost_BaggedSampleFraction","Relative size of bagged event sample to original size of the data sample (used whenever bagging is used)" );
182 
183  DeclareOptionRef( fAdaBoostBeta = 1.0, "Boost_AdaBoostBeta",
184  "The ADA boost parameter that sets the effect of every boost step on the events' weights" );
185 
186  DeclareOptionRef( fTransformString = "step", "Boost_Transform",
187  "Type of transform applied to every boosted method linear, log, step" );
188  AddPreDefVal(TString("step"));
189  AddPreDefVal(TString("linear"));
190  AddPreDefVal(TString("log"));
191  AddPreDefVal(TString("gauss"));
192 
193  DeclareOptionRef( fRandomSeed = 0, "Boost_RandomSeed",
194  "Seed for random number generator used for bagging" );
195 
196  TMVA::MethodCompositeBase::fMethods.reserve(fBoostNum);
197 }
198 
199 ////////////////////////////////////////////////////////////////////////////////
200 /// options that are used ONLY for the READER to ensure backward compatibility
201 /// they are hence without any effect (the reader is only reading the training
202 /// options that HAD been used at the training of the .xml weightfile at hand
203 
205 {
206 
208 
209  DeclareOptionRef( fHistoricOption = "ByError", "Boost_MethodWeightType",
210  "How to set the final weight of the boosted classifiers" );
211  AddPreDefVal(TString("ByError"));
212  AddPreDefVal(TString("Average"));
213  AddPreDefVal(TString("ByROC"));
214  AddPreDefVal(TString("ByOverlap"));
215  AddPreDefVal(TString("LastMethod"));
216 
217  DeclareOptionRef( fHistoricOption = "step", "Boost_Transform",
218  "Type of transform applied to every boosted method linear, log, step" );
219  AddPreDefVal(TString("step"));
220  AddPreDefVal(TString("linear"));
221  AddPreDefVal(TString("log"));
222  AddPreDefVal(TString("gauss"));
223 
224  // this option here
225  //DeclareOptionRef( fBoostType = "AdaBoost", "Boost_Type", "Boosting type for the classifiers" );
226  // still exists, but these two possible values
227  AddPreDefVal(TString("HighEdgeGauss"));
228  AddPreDefVal(TString("HighEdgeCoPara"));
229  // have been deleted .. hope that works :)
230 
231  DeclareOptionRef( fHistoricBoolOption, "Boost_RecalculateMVACut",
232  "Recalculate the classifier MVA Signallike cut at every boost iteration" );
233 
234 }
235 ////////////////////////////////////////////////////////////////////////////////
236 /// just registering the string from which the boosted classifier will be created
237 
238 Bool_t TMVA::MethodBoost::BookMethod( Types::EMVA theMethod, TString methodTitle, TString theOption )
239 {
240  fBoostedMethodName = Types::Instance().GetMethodName( theMethod );
241  fBoostedMethodTitle = methodTitle;
242  fBoostedMethodOptions = theOption;
243  TString opts=theOption;
244  opts.ToLower();
245 // if (opts.Contains("vartransform")) Log() << kFATAL << "It is not possible to use boost in conjunction with variable transform. Please remove either Boost_Num or VarTransform from the option string"<< methodTitle<<Endl;
246 
247  return kTRUE;
248 }
249 
250 ////////////////////////////////////////////////////////////////////////////////
251 
253 {
254 }
255 
256 ////////////////////////////////////////////////////////////////////////////////
257 /// initialisation routine
258 
260 {
261 
262  Results* results = Data()->GetResults(GetMethodName(), Types::kTraining, GetAnalysisType());
263 
264  results->Store(new TH1F("MethodWeight","Normalized Classifier Weight",fBoostNum,0,fBoostNum),"ClassifierWeight");
265  results->Store(new TH1F("BoostWeight","Boost Weight",fBoostNum,0,fBoostNum),"BoostWeight");
266  results->Store(new TH1F("ErrFraction","Error Fraction (by boosted event weights)",fBoostNum,0,fBoostNum),"ErrorFraction");
267  if (fDetailedMonitoring){
268  results->Store(new TH1F("ROCIntegral_test","ROC integral of single classifier (testing sample)",fBoostNum,0,fBoostNum),"ROCIntegral_test");
269  results->Store(new TH1F("ROCIntegralBoosted_test","ROC integral of boosted method (testing sample)",fBoostNum,0,fBoostNum),"ROCIntegralBoosted_test");
270  results->Store(new TH1F("ROCIntegral_train","ROC integral of single classifier (training sample)",fBoostNum,0,fBoostNum),"ROCIntegral_train");
271  results->Store(new TH1F("ROCIntegralBoosted_train","ROC integral of boosted method (training sample)",fBoostNum,0,fBoostNum),"ROCIntegralBoosted_train");
272  results->Store(new TH1F("OverlapIntegal_train","Overlap integral (training sample)",fBoostNum,0,fBoostNum),"Overlap");
273  }
274 
275 
276  results->GetHist("ClassifierWeight")->GetXaxis()->SetTitle("Index of boosted classifier");
277  results->GetHist("ClassifierWeight")->GetYaxis()->SetTitle("Classifier Weight");
278  results->GetHist("BoostWeight")->GetXaxis()->SetTitle("Index of boosted classifier");
279  results->GetHist("BoostWeight")->GetYaxis()->SetTitle("Boost Weight");
280  results->GetHist("ErrorFraction")->GetXaxis()->SetTitle("Index of boosted classifier");
281  results->GetHist("ErrorFraction")->GetYaxis()->SetTitle("Error Fraction");
282  if (fDetailedMonitoring){
283  results->GetHist("ROCIntegral_test")->GetXaxis()->SetTitle("Index of boosted classifier");
284  results->GetHist("ROCIntegral_test")->GetYaxis()->SetTitle("ROC integral of single classifier");
285  results->GetHist("ROCIntegralBoosted_test")->GetXaxis()->SetTitle("Number of boosts");
286  results->GetHist("ROCIntegralBoosted_test")->GetYaxis()->SetTitle("ROC integral boosted");
287  results->GetHist("ROCIntegral_train")->GetXaxis()->SetTitle("Index of boosted classifier");
288  results->GetHist("ROCIntegral_train")->GetYaxis()->SetTitle("ROC integral of single classifier");
289  results->GetHist("ROCIntegralBoosted_train")->GetXaxis()->SetTitle("Number of boosts");
290  results->GetHist("ROCIntegralBoosted_train")->GetYaxis()->SetTitle("ROC integral boosted");
291  results->GetHist("Overlap")->GetXaxis()->SetTitle("Index of boosted classifier");
292  results->GetHist("Overlap")->GetYaxis()->SetTitle("Overlap integral");
293  }
294 
295  results->Store(new TH1F("SoverBtotal","S/B in reweighted training sample",fBoostNum,0,fBoostNum),"SoverBtotal");
296  results->GetHist("SoverBtotal")->GetYaxis()->SetTitle("S/B (boosted sample)");
297  results->GetHist("SoverBtotal")->GetXaxis()->SetTitle("Index of boosted classifier");
298 
299  results->Store(new TH1F("SeparationGain","SeparationGain",fBoostNum,0,fBoostNum),"SeparationGain");
300  results->GetHist("SeparationGain")->GetYaxis()->SetTitle("SeparationGain");
301  results->GetHist("SeparationGain")->GetXaxis()->SetTitle("Index of boosted classifier");
302 
303 
304 
305  fMonitorTree= new TTree("MonitorBoost","Boost variables");
306  fMonitorTree->Branch("iMethod",&fCurrentMethodIdx,"iMethod/I");
307  fMonitorTree->Branch("boostWeight",&fBoostWeight,"boostWeight/D");
308  fMonitorTree->Branch("errorFraction",&fMethodError,"errorFraction/D");
309  fMonitorBoostedMethod = kTRUE;
310 
311 }
312 
313 
314 ////////////////////////////////////////////////////////////////////////////////
315 
317 {
318  Log() << kDEBUG << "CheckSetup: fBoostType="<<fBoostType << Endl;
319  Log() << kDEBUG << "CheckSetup: fAdaBoostBeta="<<fAdaBoostBeta<<Endl;
320  Log() << kDEBUG << "CheckSetup: fBoostWeight="<<fBoostWeight<<Endl;
321  Log() << kDEBUG << "CheckSetup: fMethodError="<<fMethodError<<Endl;
322  Log() << kDEBUG << "CheckSetup: fBoostNum="<<fBoostNum << Endl;
323  Log() << kDEBUG << "CheckSetup: fRandomSeed=" << fRandomSeed<< Endl;
324  Log() << kDEBUG << "CheckSetup: fTrainSigMVAHist.size()="<<fTrainSigMVAHist.size()<<Endl;
325  Log() << kDEBUG << "CheckSetup: fTestSigMVAHist.size()="<<fTestSigMVAHist.size()<<Endl;
326  Log() << kDEBUG << "CheckSetup: fMonitorBoostedMethod=" << (fMonitorBoostedMethod? "true" : "false") << Endl;
327  Log() << kDEBUG << "CheckSetup: MName=" << fBoostedMethodName << " Title="<< fBoostedMethodTitle<< Endl;
328  Log() << kDEBUG << "CheckSetup: MOptions="<< fBoostedMethodOptions << Endl;
329  Log() << kDEBUG << "CheckSetup: fMonitorTree=" << fMonitorTree <<Endl;
330  Log() << kDEBUG << "CheckSetup: fCurrentMethodIdx=" <<fCurrentMethodIdx << Endl;
331  if (fMethods.size()>0) Log() << kDEBUG << "CheckSetup: fMethods[0]" <<fMethods[0]<<Endl;
332  Log() << kDEBUG << "CheckSetup: fMethodWeight.size()" << fMethodWeight.size() << Endl;
333  if (fMethodWeight.size()>0) Log() << kDEBUG << "CheckSetup: fMethodWeight[0]="<<fMethodWeight[0]<<Endl;
334  Log() << kDEBUG << "CheckSetup: trying to repair things" << Endl;
335 
336 }
337 ////////////////////////////////////////////////////////////////////////////////
338 
340 {
341  TDirectory* methodDir( 0 );
342  TString dirName,dirTitle;
343  Int_t StopCounter=0;
344  Results* results = Data()->GetResults(GetMethodName(), Types::kTraining, GetAnalysisType());
345 
346 
347  InitHistos();
348 
349  if (Data()->GetNTrainingEvents()==0) Log() << kFATAL << "<Train> Data() has zero events" << Endl;
350  Data()->SetCurrentType(Types::kTraining);
351 
352  if (fMethods.size() > 0) fMethods.clear();
353  fMVAvalues->resize(Data()->GetNTrainingEvents(), 0.0);
354 
355  Log() << kINFO << "Training "<< fBoostNum << " " << fBoostedMethodName << " with title " << fBoostedMethodTitle << " Classifiers ... patience please" << Endl;
356  Timer timer( fBoostNum, GetName() );
357 
358  ResetBoostWeights();
359 
360  // clean boosted method options
361  CleanBoostOptions();
362 
363 
364  // remove transformations for individual boosting steps
365  // the transformation of the main method will be rerouted to each of the boost steps
366  Ssiz_t varTrafoStart=fBoostedMethodOptions.Index("~VarTransform=");
367  if (varTrafoStart >0) {
368  Ssiz_t varTrafoEnd =fBoostedMethodOptions.Index(":",varTrafoStart);
369  if (varTrafoEnd<varTrafoStart)
370  varTrafoEnd=fBoostedMethodOptions.Length();
371  fBoostedMethodOptions.Remove(varTrafoStart,varTrafoEnd-varTrafoStart);
372  }
373 
374  //
375  // training and boosting the classifiers
376  for (fCurrentMethodIdx=0;fCurrentMethodIdx<fBoostNum;fCurrentMethodIdx++) {
377  // the first classifier shows the option string output, the rest not
378  if (fCurrentMethodIdx>0) TMVA::MsgLogger::InhibitOutput();
379 
380  IMethod* method = ClassifierFactory::Instance().Create(std::string(fBoostedMethodName),
381  GetJobName(),
382  Form("%s_B%04i", fBoostedMethodTitle.Data(),fCurrentMethodIdx),
383  DataInfo(),
384  fBoostedMethodOptions);
386 
387  // supressing the rest of the classifier output the right way
388  fCurrentMethod = (dynamic_cast<MethodBase*>(method));
389 
390  if (fCurrentMethod==0) {
391  Log() << kFATAL << "uups.. guess the booking of the " << fCurrentMethodIdx << "-th classifier somehow failed" << Endl;
392  return; // hope that makes coverity happy (as if fears I migh use the pointer later on, not knowing that FATAL exits
393  }
394 
395  // set fDataSetManager if MethodCategory (to enable Category to create datasetinfo objects) // DSMTEST
396  if (fCurrentMethod->GetMethodType() == Types::kCategory) { // DSMTEST
397  MethodCategory *methCat = (dynamic_cast<MethodCategory*>(fCurrentMethod)); // DSMTEST
398  if (!methCat) // DSMTEST
399  Log() << kFATAL << "Method with type kCategory cannot be casted to MethodCategory. /MethodBoost" << Endl; // DSMTEST
400  methCat->fDataSetManager = fDataSetManager; // DSMTEST
401  } // DSMTEST
402 
403  fCurrentMethod->SetMsgType(kWARNING);
404  fCurrentMethod->SetupMethod();
405  fCurrentMethod->ParseOptions();
406  // put SetAnalysisType here for the needs of MLP
407  fCurrentMethod->SetAnalysisType( GetAnalysisType() );
408  fCurrentMethod->ProcessSetup();
409  fCurrentMethod->CheckSetup();
410 
411 
412  // reroute transformationhandler
413  fCurrentMethod->RerouteTransformationHandler (&(this->GetTransformationHandler()));
414 
415 
416  // creating the directory of the classifier
417  if (fMonitorBoostedMethod) {
418  methodDir=MethodBaseDir()->GetDirectory(dirName=Form("%s_B%04i",fBoostedMethodName.Data(),fCurrentMethodIdx));
419  if (methodDir==0) {
420  methodDir=BaseDir()->mkdir(dirName,dirTitle=Form("Directory Boosted %s #%04i", fBoostedMethodName.Data(),fCurrentMethodIdx));
421  }
422  fCurrentMethod->SetMethodDir(methodDir);
423  fCurrentMethod->BaseDir()->cd();
424  }
425 
426  // training
427  TMVA::MethodCompositeBase::fMethods.push_back(method);
428  timer.DrawProgressBar( fCurrentMethodIdx );
429  if (fCurrentMethodIdx==0) MonitorBoost(Types::kBoostProcBegin,fCurrentMethodIdx);
430  MonitorBoost(Types::kBeforeTraining,fCurrentMethodIdx);
431  TMVA::MsgLogger::InhibitOutput(); //supressing Logger outside the method
432  if (fBoostType=="Bagging") Bagging(); // you want also to train the first classifier on a bagged sample
433  SingleTrain();
435  fCurrentMethod->WriteMonitoringHistosToFile();
436 
437  // calculate MVA values of current method for all events in training sample
438  // (used later on to get 'misclassified events' etc for the boosting
439  CalcMVAValues();
440 
441  if (fCurrentMethodIdx==0 && fMonitorBoostedMethod) CreateMVAHistorgrams();
442 
443  // get ROC integral and overlap integral for single method on
444  // training sample if fMethodWeightType == "ByROC" or the user
445  // wants detailed monitoring
446 
447  // boosting (reweight training sample)
448  MonitorBoost(Types::kBeforeBoosting,fCurrentMethodIdx);
449  SingleBoost(fCurrentMethod);
450 
451  MonitorBoost(Types::kAfterBoosting,fCurrentMethodIdx);
452  results->GetHist("BoostWeight")->SetBinContent(fCurrentMethodIdx+1,fBoostWeight);
453  results->GetHist("ErrorFraction")->SetBinContent(fCurrentMethodIdx+1,fMethodError);
454 
455  if (fDetailedMonitoring) {
456  fROC_training = GetBoostROCIntegral(kTRUE, Types::kTraining, kTRUE);
457  results->GetHist("ROCIntegral_test")->SetBinContent(fCurrentMethodIdx+1, GetBoostROCIntegral(kTRUE, Types::kTesting));
458  results->GetHist("ROCIntegralBoosted_test")->SetBinContent(fCurrentMethodIdx+1, GetBoostROCIntegral(kFALSE, Types::kTesting));
459  results->GetHist("ROCIntegral_train")->SetBinContent(fCurrentMethodIdx+1, fROC_training);
460  results->GetHist("ROCIntegralBoosted_train")->SetBinContent(fCurrentMethodIdx+1, GetBoostROCIntegral(kFALSE, Types::kTraining));
461  results->GetHist("Overlap")->SetBinContent(fCurrentMethodIdx+1, fOverlap_integral);
462  }
463 
464 
465 
466  fMonitorTree->Fill();
467 
468  // stop boosting if needed when error has reached 0.5
469  // thought of counting a few steps, but it doesn't seem to be necessary
470  Log() << kDEBUG << "AdaBoost (methodErr) err = " << fMethodError << Endl;
471  if (fMethodError > 0.49999) StopCounter++;
472  if (StopCounter > 0 && fBoostType != "Bagging") {
473  timer.DrawProgressBar( fBoostNum );
474  fBoostNum = fCurrentMethodIdx+1;
475  Log() << kINFO << "Error rate has reached 0.5 ("<< fMethodError<<"), boosting process stopped at #" << fBoostNum << " classifier" << Endl;
476  if (fBoostNum < 5)
477  Log() << kINFO << "The classifier might be too strong to boost with Beta = " << fAdaBoostBeta << ", try reducing it." <<Endl;
478  break;
479  }
480  }
481 
482  //as MethodBoost acts not on a private event sample (like MethodBDT does), we need to remember not
483  // to leave "boosted" events to the next classifier in the factory
484 
485  ResetBoostWeights();
486 
487  Timer* timer1= new Timer( fBoostNum, GetName() );
488  // normalizing the weights of the classifiers
489  for (fCurrentMethodIdx=0;fCurrentMethodIdx<fBoostNum;fCurrentMethodIdx++) {
490  // pefroming post-boosting actions
491 
492  timer1->DrawProgressBar( fCurrentMethodIdx );
493 
494  if (fCurrentMethodIdx==fBoostNum) {
495  Log() << kINFO << "Elapsed time: " << timer1->GetElapsedTime()
496  << " " << Endl;
497  }
498 
499  TH1F* tmp = dynamic_cast<TH1F*>( results->GetHist("ClassifierWeight") );
500  if (tmp) tmp->SetBinContent(fCurrentMethodIdx+1,fMethodWeight[fCurrentMethodIdx]);
501 
502  }
503 
504  // Ensure that in case of only 1 boost the method weight equals
505  // 1.0. This avoids unexpected behaviour in case of very bad
506  // classifiers which have fBoostWeight=1 or fMethodError=0.5,
507  // because their weight would be set to zero. This behaviour is
508  // not ok if one boosts just one time.
509  if (fMethods.size()==1) fMethodWeight[0] = 1.0;
510 
511  MonitorBoost(Types::kBoostProcEnd);
512 
513  delete timer1;
514 }
515 
516 ////////////////////////////////////////////////////////////////////////////////
517 
519 {
520  fBoostedMethodOptions=GetOptions();
521 }
522 
523 ////////////////////////////////////////////////////////////////////////////////
524 
526 {
527  if (fBoostNum <=0) Log() << kFATAL << "CreateHistorgrams called before fBoostNum is initialized" << Endl;
528  // calculating histograms boundries and creating histograms..
529  // nrms = number of rms around the average to use for outline (of the 0 classifier)
530  Double_t meanS, meanB, rmsS, rmsB, xmin, xmax, nrms = 10;
531  Int_t signalClass = 0;
532  if (DataInfo().GetClassInfo("Signal") != 0) {
533  signalClass = DataInfo().GetClassInfo("Signal")->GetNumber();
534  }
535  gTools().ComputeStat( GetEventCollection( Types::kMaxTreeType ), fMVAvalues,
536  meanS, meanB, rmsS, rmsB, xmin, xmax, signalClass );
537 
539  xmin = TMath::Max( TMath::Min(meanS - nrms*rmsS, meanB - nrms*rmsB ), xmin );
540  xmax = TMath::Min( TMath::Max(meanS + nrms*rmsS, meanB + nrms*rmsB ), xmax ) + 0.00001;
541 
542  // creating all the historgrams
543  for (UInt_t imtd=0; imtd<fBoostNum; imtd++) {
544  fTrainSigMVAHist .push_back( new TH1F( Form("MVA_Train_S_%04i",imtd), "MVA_Train_S", fNbins, xmin, xmax ) );
545  fTrainBgdMVAHist .push_back( new TH1F( Form("MVA_Train_B%04i", imtd), "MVA_Train_B", fNbins, xmin, xmax ) );
546  fBTrainSigMVAHist.push_back( new TH1F( Form("MVA_BTrain_S%04i",imtd), "MVA_BoostedTrain_S", fNbins, xmin, xmax ) );
547  fBTrainBgdMVAHist.push_back( new TH1F( Form("MVA_BTrain_B%04i",imtd), "MVA_BoostedTrain_B", fNbins, xmin, xmax ) );
548  fTestSigMVAHist .push_back( new TH1F( Form("MVA_Test_S%04i", imtd), "MVA_Test_S", fNbins, xmin, xmax ) );
549  fTestBgdMVAHist .push_back( new TH1F( Form("MVA_Test_B%04i", imtd), "MVA_Test_B", fNbins, xmin, xmax ) );
550  }
551 }
552 
553 ////////////////////////////////////////////////////////////////////////////////
554 /// resetting back the boosted weights of the events to 1
555 
557 {
558  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) {
559  const Event *ev = Data()->GetEvent(ievt);
560  ev->SetBoostWeight( 1.0 );
561  }
562 }
563 
564 ////////////////////////////////////////////////////////////////////////////////
565 
567 {
568  TDirectory* dir=0;
569  if (fMonitorBoostedMethod) {
570  for (UInt_t imtd=0;imtd<fBoostNum;imtd++) {
571 
572  //writing the histograms in the specific classifier's directory
573  MethodBase* m = dynamic_cast<MethodBase*>(fMethods[imtd]);
574  if (!m) continue;
575  dir = m->BaseDir();
576  dir->cd();
577  fTrainSigMVAHist[imtd]->SetDirectory(dir);
578  fTrainSigMVAHist[imtd]->Write();
579  fTrainBgdMVAHist[imtd]->SetDirectory(dir);
580  fTrainBgdMVAHist[imtd]->Write();
581  fBTrainSigMVAHist[imtd]->SetDirectory(dir);
582  fBTrainSigMVAHist[imtd]->Write();
583  fBTrainBgdMVAHist[imtd]->SetDirectory(dir);
584  fBTrainBgdMVAHist[imtd]->Write();
585  }
586  }
587 
588  // going back to the original folder
589  BaseDir()->cd();
590 
591  fMonitorTree->Write();
592 }
593 
594 ////////////////////////////////////////////////////////////////////////////////
595 
597 {
599  if (fMonitorBoostedMethod) {
600  UInt_t nloop = fTestSigMVAHist.size();
601  if (fMethods.size()<nloop) nloop = fMethods.size();
602  //running over all the events and populating the test MVA histograms
603  Data()->SetCurrentType(Types::kTesting);
604  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) {
605  const Event* ev = GetEvent(ievt);
606  Float_t w = ev->GetWeight();
607  if (DataInfo().IsSignal(ev)) {
608  for (UInt_t imtd=0; imtd<nloop; imtd++) {
609  fTestSigMVAHist[imtd]->Fill(fMethods[imtd]->GetMvaValue(),w);
610  }
611  }
612  else {
613  for (UInt_t imtd=0; imtd<nloop; imtd++) {
614  fTestBgdMVAHist[imtd]->Fill(fMethods[imtd]->GetMvaValue(),w);
615  }
616  }
617  }
618  Data()->SetCurrentType(Types::kTraining);
619  }
620 }
621 
622 ////////////////////////////////////////////////////////////////////////////////
623 
625 {
627  if (treetype==Types::kTraining) return;
628  UInt_t nloop = fTestSigMVAHist.size();
629  if (fMethods.size()<nloop) nloop = fMethods.size();
630  if (fMonitorBoostedMethod) {
631  TDirectory* dir=0;
632  for (UInt_t imtd=0;imtd<nloop;imtd++) {
633  //writing the histograms in the specific classifier's directory
634  MethodBase* mva = dynamic_cast<MethodBase*>(fMethods[imtd]);
635  if (!mva) continue;
636  dir = mva->BaseDir();
637  if (dir==0) continue;
638  dir->cd();
639  fTestSigMVAHist[imtd]->SetDirectory(dir);
640  fTestSigMVAHist[imtd]->Write();
641  fTestBgdMVAHist[imtd]->SetDirectory(dir);
642  fTestBgdMVAHist[imtd]->Write();
643  }
644  }
645 }
646 
647 ////////////////////////////////////////////////////////////////////////////////
648 /// process user options
649 
651 {
652 }
653 
654 ////////////////////////////////////////////////////////////////////////////////
655 /// initialization
656 
658 {
659  Data()->SetCurrentType(Types::kTraining);
660  MethodBase* meth = dynamic_cast<MethodBase*>(GetLastMethod());
661  if (meth) meth->TrainMethod();
662 }
663 
664 ////////////////////////////////////////////////////////////////////////////////
665 /// find the CUT on the individual MVA that defines an event as
666 /// correct or misclassified (to be used in the boosting process)
667 
669 {
670  if (!method || method->GetMethodType() == Types::kDT ){ return;}
671 
672  // creating a fine histograms containing the error rate
673  const Int_t nBins=10001;
674  Double_t minMVA=150000;
675  Double_t maxMVA=-150000;
676  for (Long64_t ievt=0; ievt<Data()->GetNEvents(); ievt++) {
677  GetEvent(ievt);
678  Double_t val=method->GetMvaValue();
679  //Helge .. I think one could very well use fMVAValues for that ... -->to do
680  if (val>maxMVA) maxMVA=val;
681  if (val<minMVA) minMVA=val;
682  }
683  maxMVA = maxMVA+(maxMVA-minMVA)/nBins;
684 
685  Double_t sum = 0.;
686 
687  TH1D *mvaS = new TH1D(Form("MVAS_%d",fCurrentMethodIdx) ,"",nBins,minMVA,maxMVA);
688  TH1D *mvaB = new TH1D(Form("MVAB_%d",fCurrentMethodIdx) ,"",nBins,minMVA,maxMVA);
689  TH1D *mvaSC = new TH1D(Form("MVASC_%d",fCurrentMethodIdx),"",nBins,minMVA,maxMVA);
690  TH1D *mvaBC = new TH1D(Form("MVABC_%d",fCurrentMethodIdx),"",nBins,minMVA,maxMVA);
691 
692 
693  Results* results = Data()->GetResults(GetMethodName(), Types::kTraining, GetAnalysisType());
694  if (fDetailedMonitoring){
695  results->Store(mvaS, Form("MVAS_%d",fCurrentMethodIdx));
696  results->Store(mvaB, Form("MVAB_%d",fCurrentMethodIdx));
697  results->Store(mvaSC,Form("MVASC_%d",fCurrentMethodIdx));
698  results->Store(mvaBC,Form("MVABC_%d",fCurrentMethodIdx));
699  }
700 
701  for (Long64_t ievt=0; ievt<Data()->GetNEvents(); ievt++) {
702 
703  Double_t weight = GetEvent(ievt)->GetWeight();
704  Double_t mvaVal=method->GetMvaValue();
705  sum +=weight;
706  if (DataInfo().IsSignal(GetEvent(ievt))){
707  mvaS->Fill(mvaVal,weight);
708  }else {
709  mvaB->Fill(mvaVal,weight);
710  }
711  }
712  SeparationBase *sepGain;
713 
714 
715  // Boosting should use Miscalssification not Gini Index (changed, Helge 31.5.2013)
716  // ACHTUNG !! mit "Misclassification" geht es NUR wenn man die Signal zu Background bei jedem Boost schritt
717  // wieder hinbiegt. Es gibt aber komischerweise bessere Ergebnisse (genau wie bei BDT auch schon beobachtet) wenn
718  // man GiniIndex benutzt und akzeptiert dass jedes andere mal KEIN vernuenftiger Cut gefunden wird - d.h. der
719  // Cut liegt dann ausserhalb der MVA value range, alle events sind als Bkg classifiziert und dann wird entpsrehcend
720  // des Boost algorithmus 'automitisch' etwas renormiert .. sodass im naechsten Schritt dann wieder was vernuenftiges
721  // rauskommt. Komisch .. dass DAS richtig sein soll ??
722 
723  // SeparationBase *sepGain2 = new MisClassificationError();
724  //sepGain = new MisClassificationError();
725  sepGain = new GiniIndex();
726  //sepGain = new CrossEntropy();
727 
728  Double_t sTot = mvaS->GetSum();
729  Double_t bTot = mvaB->GetSum();
730 
731  mvaSC->SetBinContent(1,mvaS->GetBinContent(1));
732  mvaBC->SetBinContent(1,mvaB->GetBinContent(1));
733  Double_t sSel=0;
734  Double_t bSel=0;
735  Double_t separationGain=sepGain->GetSeparationGain(sSel,bSel,sTot,bTot);
736  Double_t mvaCut=mvaSC->GetBinLowEdge(1);
737  Double_t sSelCut=sSel;
738  Double_t bSelCut=bSel;
739  // std::cout << "minMVA =" << minMVA << " maxMVA = " << maxMVA << " width = " << mvaSC->GetBinWidth(1) << std::endl;
740 
741  // for (Int_t ibin=1;ibin<=nBins;ibin++) std::cout << " cutvalues[" << ibin<<"]="<<mvaSC->GetBinLowEdge(ibin) << " " << mvaSC->GetBinCenter(ibin) << std::endl;
742  Double_t mvaCutOrientation=1; // 1 if mva > mvaCut --> Signal and -1 if mva < mvaCut (i.e. mva*-1 > mvaCut*-1) --> Signal
743  for (Int_t ibin=1;ibin<=nBins;ibin++){
744  mvaSC->SetBinContent(ibin,mvaS->GetBinContent(ibin)+mvaSC->GetBinContent(ibin-1));
745  mvaBC->SetBinContent(ibin,mvaB->GetBinContent(ibin)+mvaBC->GetBinContent(ibin-1));
746 
747  sSel=mvaSC->GetBinContent(ibin);
748  bSel=mvaBC->GetBinContent(ibin);
749 
750  // if (ibin==nBins){
751  // std::cout << "Last bin s="<< sSel <<" b="<<bSel << " s="<< sTot-sSel <<" b="<<bTot-bSel << endl;
752  // }
753 
754  if (separationGain < sepGain->GetSeparationGain(sSel,bSel,sTot,bTot)
755  // && (mvaSC->GetBinCenter(ibin) >0 || (fCurrentMethodIdx+1)%2 )
756  ){
757  separationGain = sepGain->GetSeparationGain(sSel,bSel,sTot,bTot);
758  // mvaCut=mvaSC->GetBinCenter(ibin);
759  mvaCut=mvaSC->GetBinLowEdge(ibin+1);
760  // if (sSel/bSel > (sTot-sSel)/(bTot-bSel)) mvaCutOrientation=-1;
761  if (sSel*(bTot-bSel) > (sTot-sSel)*bSel) mvaCutOrientation=-1;
762  else mvaCutOrientation=1;
763  sSelCut=sSel;
764  bSelCut=bSel;
765  // std::cout << "new cut at " << mvaCut << "with s="<<sTot-sSel << " b="<<bTot-bSel << std::endl;
766  }
767  /*
768  Double_t ori;
769  if (sSel/bSel > (sTot-sSel)/(bTot-bSel)) ori=-1;
770  else ori=1;
771  std::cout << ibin << " mvacut="<<mvaCut
772  << " sTot=" << sTot
773  << " bTot=" << bTot
774  << " sSel=" << sSel
775  << " bSel=" << bSel
776  << " s/b(1)=" << sSel/bSel
777  << " s/b(2)=" << (sTot-sSel)/(bTot-bSel)
778  << " sepGain="<<sepGain->GetSeparationGain(sSel,bSel,sTot,bTot)
779  << " sepGain2="<<sepGain2->GetSeparationGain(sSel,bSel,sTot,bTot)
780  << " " <<ori
781  << std::endl;
782  */
783 
784  }
785 
786  if (0){
787  double parentIndex=sepGain->GetSeparationIndex(sTot,bTot);
788  double leftIndex =sepGain->GetSeparationIndex(sSelCut,bSelCut);
789  double rightIndex =sepGain->GetSeparationIndex(sTot-sSelCut,bTot-bSelCut);
790  std::cout
791  << " sTot=" << sTot
792  << " bTot=" << bTot
793  << " s="<<sSelCut
794  << " b="<<bSelCut
795  << " s2="<<(sTot-sSelCut)
796  << " b2="<<(bTot-bSelCut)
797  << " s/b(1)=" << sSelCut/bSelCut
798  << " s/b(2)=" << (sTot-sSelCut)/(bTot-bSelCut)
799  << " index before cut=" << parentIndex
800  << " after: left=" << leftIndex
801  << " after: right=" << rightIndex
802  << " sepGain=" << parentIndex-( (sSelCut+bSelCut) * leftIndex + (sTot-sSelCut+bTot-bSelCut) * rightIndex )/(sTot+bTot)
803  << " sepGain="<<separationGain
804  << " sepGain="<<sepGain->GetSeparationGain(sSelCut,bSelCut,sTot,bTot)
805  << " cut=" << mvaCut
806  << " idx="<<fCurrentMethodIdx
807  << " cutOrientation="<<mvaCutOrientation
808  << std::endl;
809  }
810  method->SetSignalReferenceCut(mvaCut);
811  method->SetSignalReferenceCutOrientation(mvaCutOrientation);
812 
813  results->GetHist("SeparationGain")->SetBinContent(fCurrentMethodIdx+1,separationGain);
814 
815 
816  Log() << kDEBUG << "(old step) Setting method cut to " <<method->GetSignalReferenceCut()<< Endl;
817 
818  // mvaS ->Delete();
819  // mvaB ->Delete();
820  // mvaSC->Delete();
821  // mvaBC->Delete();
822 }
823 
824 ////////////////////////////////////////////////////////////////////////////////
825 
827 {
828  Double_t returnVal=-1;
829 
830 
831  if (fBoostType=="AdaBoost") returnVal = this->AdaBoost (method,1);
832  else if (fBoostType=="RealAdaBoost") returnVal = this->AdaBoost (method,0);
833  else if (fBoostType=="Bagging") returnVal = this->Bagging ();
834  else{
835  Log() << kFATAL << "<Boost> unknown boost option " << fBoostType<< " called" << Endl;
836  }
837  fMethodWeight.push_back(returnVal);
838  return returnVal;
839 }
840 ////////////////////////////////////////////////////////////////////////////////
841 /// the standard (discrete or real) AdaBoost algorithm
842 
844 {
845  if (!method) {
846  Log() << kWARNING << " AdaBoost called without classifier reference - needed for calulating AdaBoost " << Endl;
847  return 0;
848  }
849 
850  Float_t w,v; Bool_t sig=kTRUE;
851  Double_t sumAll=0, sumWrong=0;
852  Bool_t* WrongDetection=new Bool_t[GetNEvents()];
853  QuickMVAProbEstimator *MVAProb=NULL;
854 
855  if (discreteAdaBoost) {
856  FindMVACut(method);
857  Log() << kDEBUG << " individual mva cut value = " << method->GetSignalReferenceCut() << Endl;
858  } else {
859  MVAProb=new TMVA::QuickMVAProbEstimator();
860  // the RealAdaBoost does use a simple "yes (signal)" or "no (background)"
861  // answer from your single MVA, but a "signal probability" instead (in the BDT case,
862  // that would be the 'purity' in the leaf node. For some MLP parameter, the MVA output
863  // can also interpreted as a probability, but here I try a genera aproach to get this
864  // probability from the MVA distributions...
865 
866  for (Long64_t evt=0; evt<GetNEvents(); evt++) {
867  const Event* ev = Data()->GetEvent(evt);
868  MVAProb->AddEvent(fMVAvalues->at(evt),ev->GetWeight(),ev->GetClass());
869  }
870  }
871 
872 
873  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) WrongDetection[ievt]=kTRUE;
874 
875  // finding the wrong events and calculating their total weights
876  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) {
877  const Event* ev = GetEvent(ievt);
878  sig=DataInfo().IsSignal(ev);
879  v = fMVAvalues->at(ievt);
880  w = ev->GetWeight();
881  sumAll += w;
882  if (fMonitorBoostedMethod) {
883  if (sig) {
884  fBTrainSigMVAHist[fCurrentMethodIdx]->Fill(v,w);
885  fTrainSigMVAHist[fCurrentMethodIdx]->Fill(v,ev->GetOriginalWeight());
886  }
887  else {
888  fBTrainBgdMVAHist[fCurrentMethodIdx]->Fill(v,w);
889  fTrainBgdMVAHist[fCurrentMethodIdx]->Fill(v,ev->GetOriginalWeight());
890  }
891  }
892 
893  if (discreteAdaBoost){
894  if (sig == method->IsSignalLike(fMVAvalues->at(ievt))){
895  WrongDetection[ievt]=kFALSE;
896  }else{
897  WrongDetection[ievt]=kTRUE;
898  sumWrong+=w;
899  }
900  }else{
901  Double_t mvaProb = MVAProb->GetMVAProbAt((Float_t)fMVAvalues->at(ievt));
902  mvaProb = 2*(mvaProb-0.5);
903  Int_t trueType;
904  if (DataInfo().IsSignal(ev)) trueType = 1;
905  else trueType = -1;
906  sumWrong+= w*trueType*mvaProb;
907  }
908  }
909 
910  fMethodError=sumWrong/sumAll;
911 
912  // calculating the fMethodError and the boostWeight out of it uses the formula
913  // w = ((1-err)/err)^beta
914 
915  Double_t boostWeight=0;
916 
917  if (fMethodError == 0) { //no misclassification made.. perfect, no boost ;)
918  Log() << kWARNING << "Your classifier worked perfectly on the training sample --> serious overtraining expected and no boosting done " << Endl;
919  }else{
920 
921  if (discreteAdaBoost)
922  boostWeight = TMath::Log((1.-fMethodError)/fMethodError)*fAdaBoostBeta;
923  else
924  boostWeight = TMath::Log((1.+fMethodError)/(1-fMethodError))*fAdaBoostBeta;
925 
926 
927  // std::cout << "boostweight = " << boostWeight << std::endl;
928 
929  // ADA boosting, rescaling the weight of the wrong events according to the error level
930  // over the entire test sample rescaling all the weights to have the same sum, but without
931  // touching the original weights (changing only the boosted weight of all the events)
932  // first reweight
933 
934  Double_t newSum=0., oldSum=0.;
935 
936 
937  Double_t boostfactor = TMath::Exp(boostWeight);
938 
939 
940  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) {
941  const Event* ev = Data()->GetEvent(ievt);
942  oldSum += ev->GetWeight();
943  if (discreteAdaBoost){
944  // events are classified as Signal OR background .. right or wrong
945  if (WrongDetection[ievt] && boostWeight != 0) {
946  if (ev->GetWeight() > 0) ev->ScaleBoostWeight(boostfactor);
947  else ev->ScaleBoostWeight(1./boostfactor);
948  }
949  // if (ievt<30) std::cout<<ievt<<" var0="<<ev->GetValue(0)<<" var1="<<ev->GetValue(1)<<" weight="<<ev->GetWeight() << " boostby:"<<boostfactor<<std::endl;
950 
951  }else{
952  // events are classified by their probability of being signal or background
953  // (eventually you should write this one - i.e. re-use the MVA value that were already
954  // calcualted and stroed.. however ,for the moement ..
955  Double_t mvaProb = MVAProb->GetMVAProbAt((Float_t)fMVAvalues->at(ievt));
956  mvaProb = 2*(mvaProb-0.5);
957  // mvaProb = (1-mvaProb);
958 
959  Int_t trueType=1;
960  if (DataInfo().IsSignal(ev)) trueType = 1;
961  else trueType = -1;
962 
963  boostfactor = TMath::Exp(-1*boostWeight*trueType*mvaProb);
964  if (ev->GetWeight() > 0) ev->ScaleBoostWeight(boostfactor);
965  else ev->ScaleBoostWeight(1./boostfactor);
966 
967  }
968  newSum += ev->GetWeight();
969  }
970 
971  Double_t normWeight = oldSum/newSum;
972  // next normalize the weights
973  Double_t normSig=0, normBkg=0;
974  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) {
975  const Event* ev = Data()->GetEvent(ievt);
976  ev->ScaleBoostWeight(normWeight);
977  if (ev->GetClass()) normSig+=ev->GetWeight();
978  else normBkg+=ev->GetWeight();
979  }
980 
981  Results* results = Data()->GetResults(GetMethodName(), Types::kTraining, GetAnalysisType());
982  results->GetHist("SoverBtotal")->SetBinContent(fCurrentMethodIdx+1, normSig/normBkg);
983 
984  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) {
985  const Event* ev = Data()->GetEvent(ievt);
986 
987  if (ev->GetClass()) ev->ScaleBoostWeight(oldSum/normSig/2);
988  else ev->ScaleBoostWeight(oldSum/normBkg/2);
989  }
990  }
991 
992  delete[] WrongDetection;
993  if (MVAProb) delete MVAProb;
994 
995  fBoostWeight = boostWeight; // used ONLY for the monitoring tree
996 
997  return boostWeight;
998 }
999 
1000 
1001 ////////////////////////////////////////////////////////////////////////////////
1002 /// Bagging or Bootstrap boosting, gives new random poisson weight for every event
1003 
1005 {
1006  TRandom3 *trandom = new TRandom3(fRandomSeed+fMethods.size());
1007  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) {
1008  const Event* ev = Data()->GetEvent(ievt);
1009  ev->SetBoostWeight(trandom->PoissonD(fBaggedSampleFraction));
1010  }
1011  fBoostWeight = 1; // used ONLY for the monitoring tree
1012  return 1.;
1013 }
1014 
1015 
1016 ////////////////////////////////////////////////////////////////////////////////
1017 /// Get help message text
1018 ///
1019 /// typical length of text line:
1020 /// "|--------------------------------------------------------------|"
1021 
1023 {
1024  Log() << Endl;
1025  Log() << gTools().Color("bold") << "--- Short description:" << gTools().Color("reset") << Endl;
1026  Log() << Endl;
1027  Log() << "This method combines several classifier of one species in a "<<Endl;
1028  Log() << "single multivariate quantity via the boost algorithm." << Endl;
1029  Log() << "the output is a weighted sum over all individual classifiers" <<Endl;
1030  Log() << "By default, the AdaBoost method is employed, which gives " << Endl;
1031  Log() << "events that were misclassified in the previous tree a larger " << Endl;
1032  Log() << "weight in the training of the following classifier."<<Endl;
1033  Log() << "Optionally, Bagged boosting can also be applied." << Endl;
1034  Log() << Endl;
1035  Log() << gTools().Color("bold") << "--- Performance tuning via configuration options:" << gTools().Color("reset") << Endl;
1036  Log() << Endl;
1037  Log() << "The most important parameter in the configuration is the "<<Endl;
1038  Log() << "number of boosts applied (Boost_Num) and the choice of boosting"<<Endl;
1039  Log() << "(Boost_Type), which can be set to either AdaBoost or Bagging." << Endl;
1040  Log() << "AdaBoosting: The most important parameters in this configuration" <<Endl;
1041  Log() << "is the beta parameter (Boost_AdaBoostBeta) " << Endl;
1042  Log() << "When boosting a linear classifier, it is sometimes advantageous"<<Endl;
1043  Log() << "to transform the MVA output non-linearly. The following options" <<Endl;
1044  Log() << "are available: step, log, and minmax, the default is no transform."<<Endl;
1045  Log() <<Endl;
1046  Log() << "Some classifiers are hard to boost and do not improve much in"<<Endl;
1047  Log() << "their performance by boosting them, some even slightly deteriorate"<< Endl;
1048  Log() << "due to the boosting." <<Endl;
1049  Log() << "The booking of the boost method is special since it requires"<<Endl;
1050  Log() << "the booing of the method to be boosted and the boost itself."<<Endl;
1051  Log() << "This is solved by booking the method to be boosted and to add"<<Endl;
1052  Log() << "all Boost parameters, which all begin with \"Boost_\" to the"<<Endl;
1053  Log() << "options string. The factory separates the options and initiates"<<Endl;
1054  Log() << "the boost process. The TMVA macro directory contains the example"<<Endl;
1055  Log() << "macro \"Boost.C\"" <<Endl;
1056 }
1057 
1058 ////////////////////////////////////////////////////////////////////////////////
1059 
1061 {
1062  return 0;
1063 }
1064 
1065 ////////////////////////////////////////////////////////////////////////////////
1066 /// return boosted MVA response
1067 
1069 {
1070  Double_t mvaValue = 0;
1071  Double_t norm = 0;
1072  Double_t epsilon = TMath::Exp(-1.);
1073  //Double_t fact = TMath::Exp(-1.)+TMath::Exp(1.);
1074  for (UInt_t i=0;i< fMethods.size(); i++){
1075  MethodBase* m = dynamic_cast<MethodBase*>(fMethods[i]);
1076  if (m==0) continue;
1077  Double_t val = fTmpEvent ? m->GetMvaValue(fTmpEvent) : m->GetMvaValue();
1078  Double_t sigcut = m->GetSignalReferenceCut();
1079 
1080  // default is no transform
1081  if (fTransformString == "linear"){
1082 
1083  }
1084  else if (fTransformString == "log"){
1085  if (val < sigcut) val = sigcut;
1086 
1087  val = TMath::Log((val-sigcut)+epsilon);
1088  }
1089  else if (fTransformString == "step" ){
1090  if (m->IsSignalLike(val)) val = 1.;
1091  else val = -1.;
1092  }
1093  else if (fTransformString == "gauss"){
1094  val = TMath::Gaus((val-sigcut),1);
1095  }
1096  else {
1097  Log() << kFATAL << "error unknown transformation " << fTransformString<<Endl;
1098  }
1099  mvaValue+=val*fMethodWeight[i];
1100  norm +=fMethodWeight[i];
1101  // std::cout << "mva("<<i<<") = "<<val<<" " << valx<< " " << mvaValue<<" and sigcut="<<sigcut << std::endl;
1102  }
1103  mvaValue/=norm;
1104  // cannot determine error
1105  NoErrorCalc(err, errUpper);
1106 
1107  return mvaValue;
1108 }
1109 
1110 ////////////////////////////////////////////////////////////////////////////////
1111 /// Calculate the ROC integral of a single classifier or even the
1112 /// whole boosted classifier. The tree type (training or testing
1113 /// sample) is specified by 'eTT'.
1114 ///
1115 /// If tree type kTraining is set, the original training sample is
1116 /// used to compute the ROC integral (original weights).
1117 ///
1118 /// - singleMethod - if kTRUE, return ROC integral of single (last
1119 /// trained) classifier; if kFALSE, return ROC
1120 /// integral of full classifier
1121 ///
1122 /// - eTT - tree type (Types::kTraining / Types::kTesting)
1123 ///
1124 /// - CalcOverlapIntergral - if kTRUE, the overlap integral of the
1125 /// signal/background MVA distributions
1126 /// is calculated and stored in
1127 /// 'fOverlap_integral'
1128 
1130 {
1131  // set data sample training / testing
1132  Data()->SetCurrentType(eTT);
1133 
1134  MethodBase* method = singleMethod ? dynamic_cast<MethodBase*>(fMethods.back()) : 0; // ToDo CoVerity flags this line as there is no prtection against a zero-pointer delivered by dynamic_cast
1135  // to make CoVerity happy (although, OF COURSE, the last method in the commitee
1136  // has to be also of type MethodBase as ANY method is... hence the dynamic_cast
1137  // will never by "zero" ...
1138  if (singleMethod && !method) {
1139  Log() << kFATAL << " What do you do? Your method:"
1140  << fMethods.back()->GetName()
1141  << " seems not to be a propper TMVA method"
1142  << Endl;
1143  std::exit(1);
1144  }
1145  Double_t err = 0.0;
1146 
1147  // temporary renormalize the method weights in case of evaluation
1148  // of full classifier.
1149  // save the old normalization of the methods
1150  std::vector<Double_t> OldMethodWeight(fMethodWeight);
1151  if (!singleMethod) {
1152  // calculate sum of weights of all methods
1153  Double_t AllMethodsWeight = 0;
1154  for (UInt_t i=0; i<=fCurrentMethodIdx; i++)
1155  AllMethodsWeight += fMethodWeight.at(i);
1156  // normalize the weights of the classifiers
1157  if (AllMethodsWeight != 0.0) {
1158  for (UInt_t i=0; i<=fCurrentMethodIdx; i++)
1159  fMethodWeight[i] /= AllMethodsWeight;
1160  }
1161  }
1162 
1163  // calculate MVA values
1164  Double_t meanS, meanB, rmsS, rmsB, xmin, xmax, nrms = 10;
1165  std::vector <Float_t>* mvaRes;
1166  if (singleMethod && eTT==Types::kTraining)
1167  mvaRes = fMVAvalues; // values already calculated
1168  else {
1169  mvaRes = new std::vector <Float_t>(GetNEvents());
1170  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) {
1171  GetEvent(ievt);
1172  (*mvaRes)[ievt] = singleMethod ? method->GetMvaValue(&err) : GetMvaValue(&err);
1173  }
1174  }
1175 
1176  // restore the method weights
1177  if (!singleMethod)
1178  fMethodWeight = OldMethodWeight;
1179 
1180  // now create histograms for calculation of the ROC integral
1181  Int_t signalClass = 0;
1182  if (DataInfo().GetClassInfo("Signal") != 0) {
1183  signalClass = DataInfo().GetClassInfo("Signal")->GetNumber();
1184  }
1185  gTools().ComputeStat( GetEventCollection(eTT), mvaRes,
1186  meanS, meanB, rmsS, rmsB, xmin, xmax, signalClass );
1187 
1189  xmin = TMath::Max( TMath::Min(meanS - nrms*rmsS, meanB - nrms*rmsB ), xmin );
1190  xmax = TMath::Min( TMath::Max(meanS + nrms*rmsS, meanB + nrms*rmsB ), xmax ) + 0.0001;
1191 
1192  // calculate ROC integral
1193  TH1* mva_s = new TH1F( "MVA_S", "MVA_S", fNbins, xmin, xmax );
1194  TH1* mva_b = new TH1F( "MVA_B", "MVA_B", fNbins, xmin, xmax );
1195  TH1 *mva_s_overlap=0, *mva_b_overlap=0;
1196  if (CalcOverlapIntergral) {
1197  mva_s_overlap = new TH1F( "MVA_S_OVERLAP", "MVA_S_OVERLAP", fNbins, xmin, xmax );
1198  mva_b_overlap = new TH1F( "MVA_B_OVERLAP", "MVA_B_OVERLAP", fNbins, xmin, xmax );
1199  }
1200  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) {
1201  const Event* ev = GetEvent(ievt);
1202  Float_t w = (eTT==Types::kTesting ? ev->GetWeight() : ev->GetOriginalWeight());
1203  if (DataInfo().IsSignal(ev)) mva_s->Fill( (*mvaRes)[ievt], w );
1204  else mva_b->Fill( (*mvaRes)[ievt], w );
1205 
1206  if (CalcOverlapIntergral) {
1207  Float_t w_ov = ev->GetWeight();
1208  if (DataInfo().IsSignal(ev))
1209  mva_s_overlap->Fill( (*mvaRes)[ievt], w_ov );
1210  else
1211  mva_b_overlap->Fill( (*mvaRes)[ievt], w_ov );
1212  }
1213  }
1214  gTools().NormHist( mva_s );
1215  gTools().NormHist( mva_b );
1216  PDF *fS = new PDF( "PDF Sig", mva_s, PDF::kSpline2 );
1217  PDF *fB = new PDF( "PDF Bkg", mva_b, PDF::kSpline2 );
1218 
1219  // calculate ROC integral from fS, fB
1220  Double_t ROC = MethodBase::GetROCIntegral(fS, fB);
1221 
1222  // calculate overlap integral
1223  if (CalcOverlapIntergral) {
1224  gTools().NormHist( mva_s_overlap );
1225  gTools().NormHist( mva_b_overlap );
1226 
1227  fOverlap_integral = 0.0;
1228  for (Int_t bin=1; bin<=mva_s_overlap->GetNbinsX(); bin++){
1229  Double_t bc_s = mva_s_overlap->GetBinContent(bin);
1230  Double_t bc_b = mva_b_overlap->GetBinContent(bin);
1231  if (bc_s > 0.0 && bc_b > 0.0)
1232  fOverlap_integral += TMath::Min(bc_s, bc_b);
1233  }
1234 
1235  delete mva_s_overlap;
1236  delete mva_b_overlap;
1237  }
1238 
1239  delete mva_s;
1240  delete mva_b;
1241  delete fS;
1242  delete fB;
1243  if (!(singleMethod && eTT==Types::kTraining)) delete mvaRes;
1244 
1245  Data()->SetCurrentType(Types::kTraining);
1246 
1247  return ROC;
1248 }
1249 
1251 {
1252  // Calculate MVA values of current method fMethods.back() on
1253  // training sample
1254 
1255  Data()->SetCurrentType(Types::kTraining);
1256  MethodBase* method = dynamic_cast<MethodBase*>(fMethods.back());
1257  if (!method) {
1258  Log() << kFATAL << "dynamic cast to MethodBase* failed" <<Endl;
1259  return;
1260  }
1261  // calculate MVA values
1262  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) {
1263  GetEvent(ievt);
1264  fMVAvalues->at(ievt) = method->GetMvaValue();
1265  }
1266 
1267  // fill cumulative mva distribution
1268 
1269 
1270 }
1271 
1272 
1273 ////////////////////////////////////////////////////////////////////////////////
1274 /// fill various monitoring histograms from information of the individual classifiers that
1275 /// have been boosted.
1276 /// of course.... this depends very much on the individual classifiers, and so far, only for
1277 /// Decision Trees, this monitoring is actually implemented
1278 
1280 {
1281  Results* results = Data()->GetResults(GetMethodName(), Types::kTraining, GetAnalysisType());
1282 
1283  if (GetCurrentMethod(methodIndex)->GetMethodType() == TMVA::Types::kDT) {
1284  TMVA::MethodDT* currentDT=dynamic_cast<TMVA::MethodDT*>(GetCurrentMethod(methodIndex));
1285  if (currentDT){
1286  if (stage == Types::kBoostProcBegin){
1287  results->Store(new TH1I("NodesBeforePruning","nodes before pruning",this->GetBoostNum(),0,this->GetBoostNum()),"NodesBeforePruning");
1288  results->Store(new TH1I("NodesAfterPruning","nodes after pruning",this->GetBoostNum(),0,this->GetBoostNum()),"NodesAfterPruning");
1289  }
1290 
1291  if (stage == Types::kBeforeTraining){
1292  }
1293  else if (stage == Types::kBeforeBoosting){
1294  results->GetHist("NodesBeforePruning")->SetBinContent(methodIndex+1,currentDT->GetNNodesBeforePruning());
1295  results->GetHist("NodesAfterPruning")->SetBinContent(methodIndex+1,currentDT->GetNNodes());
1296  }
1297  else if (stage == Types::kAfterBoosting){
1298 
1299  }
1300  else if (stage != Types::kBoostProcEnd){
1301  Log() << kINFO << "<Train> average number of nodes before/after pruning : "
1302  << results->GetHist("NodesBeforePruning")->GetMean() << " / "
1303  << results->GetHist("NodesAfterPruning")->GetMean()
1304  << Endl;
1305  }
1306  }
1307 
1308  }else if (GetCurrentMethod(methodIndex)->GetMethodType() == TMVA::Types::kFisher) {
1309  if (stage == Types::kAfterBoosting){
1311  }
1312  }else{
1313  if (methodIndex < 3){
1314  Log() << kINFO << "No detailed boost monitoring for "
1315  << GetCurrentMethod(methodIndex)->GetMethodName()
1316  << " yet available " << Endl;
1317  }
1318  }
1319 
1320  //boosting plots universal for all classifiers 'typically for debug purposes only as they are not general enough'
1321 
1322  if (stage == Types::kBeforeBoosting){
1323  // if you want to display the weighted events for 2D case at each boost step:
1324  if (fDetailedMonitoring){
1325  // the following code is useful only for 2D examples - mainly illustration for debug/educational purposes:
1326  if (DataInfo().GetNVariables() == 2) {
1327  results->Store(new TH2F(Form("EventDistSig_%d",methodIndex),Form("EventDistSig_%d",methodIndex),100,0,7,100,0,7));
1328  results->GetHist(Form("EventDistSig_%d",methodIndex))->SetMarkerColor(4);
1329  results->Store(new TH2F(Form("EventDistBkg_%d",methodIndex),Form("EventDistBkg_%d",methodIndex),100,0,7,100,0,7));
1330  results->GetHist(Form("EventDistBkg_%d",methodIndex))->SetMarkerColor(2);
1331 
1332  Data()->SetCurrentType(Types::kTraining);
1333  for (Long64_t ievt=0; ievt<GetNEvents(); ievt++) {
1334  const Event* ev = GetEvent(ievt);
1335  Float_t w = ev->GetWeight();
1336  Float_t v0= ev->GetValue(0);
1337  Float_t v1= ev->GetValue(1);
1338  // if (ievt<3) std::cout<<ievt<<" var0="<<v0<<" var1="<<v1<<" weight="<<w<<std::endl;
1339  TH2* h;
1340  if (DataInfo().IsSignal(ev)) h=results->GetHist2D(Form("EventDistSig_%d",methodIndex));
1341  else h=results->GetHist2D(Form("EventDistBkg_%d",methodIndex));
1342  if (h) h->Fill(v0,v1,w);
1343  }
1344  }
1345  }
1346  }
1347 
1348  return;
1349 }
1350 
1351 
IMethod * Create(const std::string &name, const TString &job, const TString &title, DataSetInfo &dsi, const TString &option)
creates the method if needed based on the method name using the creator function the factory has stor...
static ClassifierFactory & Instance()
access to the ClassifierFactory singleton creates the instance if needed
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
Double_t GetBoostROCIntegral(Bool_t, Types::ETreeType, Bool_t CalcOverlapIntergral=kFALSE)
Calculate the ROC integral of a single classifier or even the whole boosted classifier.
float xmin
Definition: THbookFile.cxx:93
Random number generator class based on M.
Definition: TRandom3.h:29
void MonitorBoost(Types::EBoostStage stage, UInt_t methodIdx=0)
fill various monitoring histograms from information of the individual classifiers that have been boos...
std::vector< Float_t > * fMVAvalues
Definition: MethodBoost.h:192
virtual Double_t PoissonD(Double_t mean)
Generates a random number according to a Poisson law.
Definition: TRandom.cxx:414
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
TH1 * GetHist(const TString &alias) const
Definition: Results.cxx:107
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4629
long long Long64_t
Definition: RtypesCore.h:69
const Double_t * v1
Definition: TArcBall.cxx:33
void SingleTrain()
initialization
Stat_t GetSum() const
Definition: TArrayD.h:48
Double_t Bagging()
Bagging or Bootstrap boosting, gives new random poisson weight for every event.
Double_t Log(Double_t x)
Definition: TMath.h:526
virtual Double_t GetMvaValue(Double_t *errLower=0, Double_t *errUpper=0)=0
float Float_t
Definition: RtypesCore.h:53
THist< 2, float > TH2F
Definition: THist.h:321
Double_t AdaBoost(MethodBase *method, Bool_t useYesNoLeaf)
the standard (discrete or real) AdaBoost algorithm
virtual void WriteEvaluationHistosToFile(Types::ETreeType treetype)
writes all MVA evaluation histograms to file
Config & gConfig()
TH1 * h
Definition: legend2.C:5
virtual Bool_t HasAnalysisType(Types::EAnalysisType type, UInt_t numberClasses, UInt_t)
Boost can handle classification with 2 classes and regression with one regression-target.
EAnalysisType
Definition: Types.h:124
THist< 1, float > TH1F
Definition: THist.h:315
void SetSignalReferenceCutOrientation(Double_t cutOrientation)
Definition: MethodBase.h:330
void SetBoostWeight(Double_t w) const
Definition: Event.h:110
THist< 1, int > TH1I
Definition: THist.h:317
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1088
int Int_t
Definition: RtypesCore.h:41
virtual TDirectory * mkdir(const char *name, const char *title="")
Create a sub-directory and return a pointer to the created directory.
Definition: TDirectory.cxx:955
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
const Ranking * CreateRanking()
virtual Int_t GetNbinsX() const
Definition: TH1.h:296
void ResetBoostWeights()
resetting back the boosted weights of the events to 1
Double_t GetWeight() const
return the event weight - depending on whether the flag IgnoreNegWeightsInTraining is or not...
Definition: Event.cxx:376
THist< 1, double > TH1D
Definition: THist.h:314
virtual Bool_t IsSignalLike()
uses a pre-set cut on the MVA output (SetSignalReferenceCut and SetSignalReferenceCutOrientation) for...
Definition: MethodBase.cxx:838
Float_t GetValue(UInt_t ivar) const
return value of i'th variable
Definition: Event.cxx:231
TH2 * GetHist2D(const TString &alias) const
Definition: Results.cxx:116
static Types & Instance()
the the single instance of "Types" if existin already, or create it (Signleton)
Definition: Types.cxx:61
static void InhibitOutput()
Definition: MsgLogger.cxx:68
Tools & gTools()
Definition: Tools.cxx:79
void FindMVACut(MethodBase *method)
find the CUT on the individual MVA that defines an event as correct or misclassified (to be used in t...
TStopwatch timer
Definition: pirndm.C:37
MethodBoost(const TString &jobName, const TString &methodTitle, DataSetInfo &theData, const TString &theOption="", TDirectory *theTargetDir=NULL)
virtual Double_t GetBinLowEdge(Int_t bin) const
return bin lower edge for 1D historam Better to use h1.GetXaxis().GetBinLowEdge(bin) ...
Definition: TH1.cxx:8481
void AddEvent(Double_t val, Double_t weight, Int_t type)
void ProcessOptions()
process user options
Double_t SingleBoost(MethodBase *method)
virtual ~MethodBoost(void)
destructor
void ScaleBoostWeight(Double_t s) const
Definition: Event.h:111
virtual void SetMarkerColor(Color_t mcolor=1)
Definition: TAttMarker.h:51
Int_t GetNNodes()
Definition: MethodDT.h:111
std::vector< std::vector< double > > Data
Double_t NormHist(TH1 *theHist, Double_t norm=1.0)
normalises histogram
Definition: Tools.cxx:395
void GetHelpMessage() const
Get help message text.
Types::EMVA GetMethodType() const
Definition: MethodBase.h:298
Definition: PDF.h:71
Double_t GetOriginalWeight() const
Definition: Event.h:84
TString GetElapsedTime(Bool_t Scientific=kTRUE)
Definition: Timer.cxx:131
Bool_t BookMethod(Types::EMVA theMethod, TString methodTitle, TString theOption)
just registering the string from which the boosted classifier will be created
RooCmdArg Timer(Bool_t flag=kTRUE)
virtual Double_t GetMean(Int_t axis=1) const
For axis = 1,2 or 3 returns the mean value of the histogram along X,Y or Z axis.
Definition: TH1.cxx:7014
std::string GetMethodName(TCppMethod_t)
Definition: Cppyy.cxx:707
Service class for 2-Dim histogram classes.
Definition: TH2.h:36
SVector< double, 2 > v
Definition: Dict.h:5
TString GetMethodName(Types::EMVA method) const
Definition: Types.cxx:127
class TMVA::Config::VariablePlotting fVariablePlotting
virtual void WriteEvaluationHistosToFile(Types::ETreeType treetype)
writes all MVA evaluation histograms to file
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition: TH1.cxx:8543
unsigned int UInt_t
Definition: RtypesCore.h:42
TMarker * m
Definition: textangle.C:8
char * Form(const char *fmt,...)
DataSetManager * fDataSetManager
TAxis * GetYaxis()
Definition: TH1.h:320
float xmax
Definition: THbookFile.cxx:93
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
REAL epsilon
Definition: triangle.c:617
void CreateMVAHistorgrams()
Double_t Gaus(Double_t x, Double_t mean=0, Double_t sigma=1, Bool_t norm=kFALSE)
Calculate a gaussian function with mean and sigma.
Definition: TMath.cxx:453
virtual Double_t GetSeparationGain(const Double_t &nSelS, const Double_t &nSelB, const Double_t &nTotS, const Double_t &nTotB)
Separation Gain: the measure of how the quality of separation of the sample increases by splitting th...
int Ssiz_t
Definition: RtypesCore.h:63
void DeclareCompatibilityOptions()
options that are used ONLY for the READER to ensure backward compatibility they are hence without any...
Double_t Exp(Double_t x)
Definition: TMath.h:495
void CheckSetup()
check may be overridden by derived class (sometimes, eg, fitters are used which can only be implement...
double Double_t
Definition: RtypesCore.h:55
Describe directory structure in memory.
Definition: TDirectory.h:41
virtual Double_t GetROCIntegral(TH1D *histS, TH1D *histB) const
calculate the area (integral) under the ROC curve as a overall quality measure of the classification ...
int type
Definition: TGX11.cxx:120
TDirectory * BaseDir() const
returns the ROOT directory where info/histograms etc of the corresponding MVA method instance are sto...
The TH1 histogram class.
Definition: TH1.h:80
UInt_t GetClass() const
Definition: Event.h:86
Int_t GetNNodesBeforePruning()
Definition: MethodDT.h:110
virtual void TestClassification()
initialization
const TString & Color(const TString &)
human readable color strings
Definition: Tools.cxx:837
void WriteMonitoringHistosToFile(void) const
write special monitoring histograms to file dummy implementation here --------------— ...
#define REGISTER_METHOD(CLASS)
for example
std::vector< IMethod * > fMethods
Abstract ClassifierFactory template that handles arbitrary types.
Double_t GetMVAProbAt(Double_t value)
virtual Bool_t cd(const char *path=0)
Change current directory to "this" directory.
Definition: TDirectory.cxx:433
virtual void DeclareCompatibilityOptions()
options that are used ONLY for the READER to ensure backward compatibility they are hence without any...
Definition: MethodBase.cxx:599
virtual Double_t GetSeparationIndex(const Double_t &s, const Double_t &b)=0
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
void InitHistos()
initialisation routine
virtual TDirectory * GetDirectory(const char *namecycle, Bool_t printError=false, const char *funcname="GetDirectory")
Find a directory using apath.
Definition: TDirectory.cxx:336
#define NULL
Definition: Rtypes.h:82
void DrawProgressBar(Int_t, const TString &comment="")
draws progress bar in color or B&W caution:
Definition: Timer.cxx:183
A TTree object has a header with a name and a title.
Definition: TTree.h:94
Double_t GetSignalReferenceCut() const
Definition: MethodBase.h:325
void Store(TObject *obj, const char *alias=0)
Definition: Results.cxx:63
static void EnableOutput()
Definition: MsgLogger.cxx:69
const Bool_t kTRUE
Definition: Rtypes.h:91
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:285
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
Definition: TNamed.cxx:152
Double_t GetMvaValue(Double_t *err=0, Double_t *errUpper=0)
return boosted MVA response
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
void ComputeStat(const std::vector< TMVA::Event * > &, std::vector< Float_t > *, Double_t &, Double_t &, Double_t &, Double_t &, Double_t &, Double_t &, Int_t signalClass, Bool_t norm=kFALSE)
sanity check
Definition: Tools.cxx:215
virtual void TestClassification()
initialization
Definition: math.cpp:60
TAxis * GetXaxis()
Definition: TH1.h:319
void SetSignalReferenceCut(Double_t cut)
Definition: MethodBase.h:329
ClassImp(TMVA::MethodBoost) TMVA
Definition: MethodBoost.cxx:77