Logo ROOT   6.08/07
Reference Guide
OptimizeConfigParameters.cxx
Go to the documentation of this file.
1 /**********************************************************************************
2  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
3  * Package: TMVA *
4  * Class : OptimizeConfigParameters *
5  * Web : http://tmva.sourceforge.net *
6  * *
7  * Description: The OptimizeConfigParameters takes care of "scanning/fitting" *
8  * different tuning parameters in order to find the best set of *
9  * tuning paraemters which will be used in the end *
10  * *
11  * Authors (alphabetical): *
12  * Helge Voss <Helge.Voss@cern.ch> - MPI-K Heidelberg, Germany *
13  * *
14  * Copyright (c) 2005: *
15  * CERN, Switzerland *
16  * MPI-K Heidelberg, Germany *
17  * *
18  * Redistribution and use in source and binary forms, with or without *
19  * modification, are permitted according to the terms listed in LICENSE *
20  * (http://ttmva.sourceforge.net/LICENSE) *
21  **********************************************************************************/
22 
24 
25 #include "TMVA/DataSet.h"
26 #include "TMVA/DataSetInfo.h"
27 #include "TMVA/Event.h"
28 #include "TMVA/IFitterTarget.h"
29 #include "TMVA/FitterBase.h"
30 #include "TMVA/GeneticFitter.h"
31 #include "TMVA/IMethod.h"
32 #include "TMVA/Interval.h"
33 #include "TMVA/MethodBase.h"
34 #include "TMVA/MethodFDA.h"
35 #include "TMVA/MsgLogger.h"
36 #include "TMVA/MinuitFitter.h"
37 #include "TMVA/PDF.h"
38 #include "TMVA/Tools.h"
39 #include "TMVA/Types.h"
40 
41 #include "TDirectory.h"
42 #include "TGraph.h"
43 #include "TH1.h"
44 #include "TH2.h"
45 #include "TMath.h"
46 
47 #include <cstdlib>
48 #include <limits>
49 
50 
52 
53 ////////////////////////////////////////////////////////////////////////////////
54 /// Constructor which sets either "Classification or Regression"
55 
56 TMVA::OptimizeConfigParameters::OptimizeConfigParameters(MethodBase * const method, std::map<TString,TMVA::Interval*> tuneParameters, TString fomType, TString optimizationFitType)
57 : fMethod(method),
58  fTuneParameters(tuneParameters),
59  fFOMType(fomType),
60  fOptimizationFitType(optimizationFitType),
61  fMvaSig(NULL),
62  fMvaBkg(NULL),
63  fMvaSigFineBin(NULL),
64  fMvaBkgFineBin(NULL),
65  fNotDoneYet(kFALSE)
66 {
67  std::string name = "OptimizeConfigParameters_";
68  name += std::string(GetMethod()->GetName());
69  fLogger = new MsgLogger(name);
70  if (fMethod->DoRegression()){
71  Log() << kFATAL << " ERROR: Sorry, Regression is not yet implement for automatic parameter optimization"
72  << " --> exit" << Endl;
73  }
74 
75  Log() << kINFO << "Automatic optimisation of tuning parameters in "
76  << GetMethod()->GetName() << " uses:" << Endl;
77 
78  std::map<TString,TMVA::Interval*>::iterator it;
79  for (it=fTuneParameters.begin(); it!=fTuneParameters.end();it++) {
80  Log() << kINFO << it->first
81  << " in range from: " << it->second->GetMin()
82  << " to: " << it->second->GetMax()
83  << " in : " << it->second->GetNbins() << " steps"
84  << Endl;
85  }
86  Log() << kINFO << " using the options: " << fFOMType << " and " << fOptimizationFitType << Endl;
87 
88 
89 
90 }
91 
92 ////////////////////////////////////////////////////////////////////////////////
93 /// the destructor (delete the OptimizeConfigParameters, store the graph and .. delete it)
94 
96 {
97  if(!GetMethod()->IsSilentFile()) GetMethod()->BaseDir()->cd();
98  Int_t n=Int_t(fFOMvsIter.size());
99  Float_t *x = new Float_t[n];
100  Float_t *y = new Float_t[n];
101  Float_t ymin=+999999999;
102  Float_t ymax=-999999999;
103 
104  for (Int_t i=0;i<n;i++){
105  x[i] = Float_t(i);
106  y[i] = fFOMvsIter[i];
107  if (ymin>y[i]) ymin=y[i];
108  if (ymax<y[i]) ymax=y[i];
109  }
110 
111  TH2D *h=new TH2D(TString(GetMethod()->GetName())+"_FOMvsIterFrame","",2,0,n,2,ymin*0.95,ymax*1.05);
112  h->SetXTitle("#iteration "+fOptimizationFitType);
113  h->SetYTitle(fFOMType);
114  TGraph *gFOMvsIter = new TGraph(n,x,y);
115  gFOMvsIter->SetName((TString(GetMethod()->GetName())+"_FOMvsIter").Data());
116  if(!GetMethod()->IsSilentFile()) gFOMvsIter->Write();
117  if(!GetMethod()->IsSilentFile()) h->Write();
118 
119  delete [] x;
120  delete [] y;
121  // delete fFOMvsIter;
122 }
123 ////////////////////////////////////////////////////////////////////////////////
124 
125 std::map<TString,Double_t> TMVA::OptimizeConfigParameters::optimize()
126 {
127  if (fOptimizationFitType == "Scan" ) this->optimizeScan();
128  else if (fOptimizationFitType == "FitGA" || fOptimizationFitType == "Minuit" ) this->optimizeFit();
129  else {
130  Log() << kFATAL << "You have chosen as optimization type " << fOptimizationFitType
131  << " that is not (yet) coded --> exit()" << Endl;
132  }
133 
134  Log() << kINFO << "For " << GetMethod()->GetName() << " the optimized Parameters are: " << Endl;
135  std::map<TString,Double_t>::iterator it;
136  for(it=fTunedParameters.begin(); it!= fTunedParameters.end(); it++){
137  Log() << kINFO << it->first << " = " << it->second << Endl;
138  }
139  return fTunedParameters;
140 
141 }
142 
143 ////////////////////////////////////////////////////////////////////////////////
144 /// helper function to scan through the all the combinations in the
145 /// parameter space
146 
147 std::vector< int > TMVA::OptimizeConfigParameters::GetScanIndices( int val, std::vector<int> base){
148  std::vector < int > indices;
149  for (UInt_t i=0; i< base.size(); i++){
150  indices.push_back(val % base[i] );
151  val = int( floor( float(val) / float(base[i]) ) );
152  }
153  return indices;
154 }
155 
156 ////////////////////////////////////////////////////////////////////////////////
157 /// do the actual optimization using a simple scan method,
158 /// i.e. calcualte the FOM for
159 /// different tuning paraemters and remember which one is
160 /// gave the best FOM
161 
163 {
164 
165  Double_t bestFOM=-1000000, currentFOM;
166 
167  std::map<TString,Double_t> currentParameters;
168  std::map<TString,TMVA::Interval*>::iterator it;
169 
170  // for the scan, start at the lower end of the interval and then "move upwards"
171  // initialize all parameters in currentParameter
172  currentParameters.clear();
173  fTunedParameters.clear();
174 
175  for (it=fTuneParameters.begin(); it!=fTuneParameters.end(); it++){
176  currentParameters.insert(std::pair<TString,Double_t>(it->first,it->second->GetMin()));
177  fTunedParameters.insert(std::pair<TString,Double_t>(it->first,it->second->GetMin()));
178  }
179  // now loop over all the parameters and get for each combination the figure of merit
180 
181  // in order to loop over all the parameters, I first create an "array" (tune parameters)
182  // of arrays (the different values of the tune parameter)
183 
184  std::vector< std::vector <Double_t> > v;
185  for (it=fTuneParameters.begin(); it!=fTuneParameters.end(); it++){
186  std::vector< Double_t > tmp;
187  for (Int_t k=0; k<it->second->GetNbins(); k++){
188  tmp.push_back(it->second->GetElement(k));
189  }
190  v.push_back(tmp);
191  }
192  Int_t Ntot = 1;
193  std::vector< int > Nindividual;
194  for (UInt_t i=0; i<v.size(); i++) {
195  Ntot *= v[i].size();
196  Nindividual.push_back(v[i].size());
197  }
198  //loop on the total number of differnt combinations
199 
200  for (int i=0; i<Ntot; i++){
201  UInt_t index=0;
202  std::vector<int> indices = GetScanIndices(i, Nindividual );
203  for (it=fTuneParameters.begin(), index=0; index< indices.size(); index++, it++){
204  currentParameters[it->first] = v[index][indices[index]];
205  }
206  Log() << kINFO << "--------------------------" << Endl;
207  Log() << kINFO <<"Settings being evaluated:" << Endl;
208  for (std::map<TString,Double_t>::iterator it_print=currentParameters.begin();
209  it_print!=currentParameters.end(); it_print++){
210  Log() << kINFO << " " << it_print->first << " = " << it_print->second << Endl;
211  }
212 
213  GetMethod()->Reset();
214  GetMethod()->SetTuneParameters(currentParameters);
215  // now do the training for the current parameters:
216  if(!GetMethod()->IsSilentFile()) GetMethod()->BaseDir()->cd();
218  GetMethod()->Data()->GetEventCollection());
220  GetMethod()->Train();
222  currentFOM = GetFOM();
223  Log() << kINFO << "FOM was found : " << currentFOM << "; current best is " << bestFOM << Endl;
224 
225  if (currentFOM > bestFOM) {
226  bestFOM = currentFOM;
227  for (std::map<TString,Double_t>::iterator iter=currentParameters.begin();
228  iter != currentParameters.end(); iter++){
229  fTunedParameters[iter->first]=iter->second;
230  }
231  }
232  }
233 
234  GetMethod()->Reset();
236 }
237 
239 {
240  // ranges (intervals) in which the fit varies the parameters
241  std::vector<TMVA::Interval*> ranges; // intervals of the fit ranges
242  std::map<TString, TMVA::Interval*>::iterator it;
243  std::vector<Double_t> pars; // current (starting) fit parameters
244 
245  for (it=fTuneParameters.begin(); it != fTuneParameters.end(); it++){
246  ranges.push_back(new TMVA::Interval(*(it->second)));
247  pars.push_back( (it->second)->GetMean() ); // like this the order is "right". Always keep the
248  // order in the vector "pars" the same as the iterator
249  // iterates through the tuneParameters !!!!
250  }
251 
252  // added to allow for transformation on input variables i.e. norm
253  GetMethod()->GetTransformationHandler().CalcTransformations(GetMethod()->Data()->GetEventCollection());
254 
255  // create the fitter
256 
257  FitterBase* fitter = NULL;
258 
259  if ( fOptimizationFitType == "Minuit" ) {
260  TString opt="";
261  fitter = new MinuitFitter( *this,
262  "FitterMinuit_BDTOptimize",
263  ranges, opt );
264  }else if ( fOptimizationFitType == "FitGA" ) {
265  TString opt="PopSize=20:Steps=30:Cycles=3:ConvCrit=0.01:SaveBestCycle=5";
266  fitter = new GeneticFitter( *this,
267  "FitterGA_BDTOptimize",
268  ranges, opt );
269  } else {
270  Log() << kWARNING << " you did not specify a valid OptimizationFitType "
271  << " will use the default (FitGA) " << Endl;
272  TString opt="PopSize=20:Steps=30:Cycles=3:ConvCrit=0.01:SaveBestCycle=5";
273  fitter = new GeneticFitter( *this,
274  "FitterGA_BDTOptimize",
275  ranges, opt );
276  }
277 
278  fitter->CheckForUnusedOptions();
279 
280  // perform the fit
281  fitter->Run(pars);
282 
283  // clean up
284  for (UInt_t ipar=0; ipar<ranges.size(); ipar++) delete ranges[ipar];
285 
286 
287  GetMethod()->Reset();
288 
289  fTunedParameters.clear();
290  Int_t jcount=0;
291  for (it=fTuneParameters.begin(); it!=fTuneParameters.end(); it++){
292  fTunedParameters.insert(std::pair<TString,Double_t>(it->first,pars[jcount++]));
293  }
294 
296 
297 }
298 
299 ////////////////////////////////////////////////////////////////////////////////
300 /// return the estimator (from current FOM) for the fitting interface
301 
303 {
304  std::map< std::vector<Double_t> , Double_t>::const_iterator iter;
305  iter = fAlreadyTrainedParCombination.find(pars);
306 
307  if (iter != fAlreadyTrainedParCombination.end()) {
308  // std::cout << "I had trained Depth=" <<Int_t(pars[0])
309  // <<" MinEv=" <<Int_t(pars[1])
310  // <<" already --> FOM="<< iter->second <<std::endl;
311  return iter->second;
312  }else{
313  std::map<TString,Double_t> currentParameters;
314  Int_t icount =0; // map "pars" to the map of Tuneparameter, make sure
315  // you never screw up this order!!
316  std::map<TString, TMVA::Interval*>::iterator it;
317  for (it=fTuneParameters.begin(); it!=fTuneParameters.end(); it++){
318  currentParameters[it->first] = pars[icount++];
319  }
320  GetMethod()->Reset();
321  GetMethod()->SetTuneParameters(currentParameters);
322  if(!GetMethod()->IsSilentFile()) GetMethod()->BaseDir()->cd();
323 
324  if (fNotDoneYet){
326  CalcTransformations(GetMethod()->Data()->GetEventCollection());
328  }
330  GetMethod()->Train();
332 
333 
334  Double_t currentFOM = GetFOM();
335 
336  fAlreadyTrainedParCombination.insert(std::make_pair(pars,-currentFOM));
337  return -currentFOM;
338  }
339 }
340 
341 ////////////////////////////////////////////////////////////////////////////////
342 /// Return the Figure of Merit (FOM) used in the parameter
343 /// optimization process
344 
346 {
347  Double_t fom=0;
348  if (fMethod->DoRegression()){
349  std::cout << " ERROR: Sorry, Regression is not yet implement for automatic parameter optimisation"
350  << " --> exit" << std::endl;
351  std::exit(1);
352  }else{
353  if (fFOMType == "Separation") fom = GetSeparation();
354  else if (fFOMType == "ROCIntegral") fom = GetROCIntegral();
355  else if (fFOMType == "SigEffAtBkgEff01") fom = GetSigEffAtBkgEff(0.1);
356  else if (fFOMType == "SigEffAtBkgEff001") fom = GetSigEffAtBkgEff(0.01);
357  else if (fFOMType == "SigEffAtBkgEff002") fom = GetSigEffAtBkgEff(0.02);
358  else if (fFOMType == "BkgRejAtSigEff05") fom = GetBkgRejAtSigEff(0.5);
359  else if (fFOMType == "BkgEffAtSigEff05") fom = GetBkgEffAtSigEff(0.5);
360  else {
361  Log()<<kFATAL << " ERROR, you've specified as Figure of Merit in the "
362  << " parameter optimisation " << fFOMType << " which has not"
363  << " been implemented yet!! ---> exit " << Endl;
364  }
365  }
366  fFOMvsIter.push_back(fom);
367  // std::cout << "fom="<<fom<<std::endl; // should write that into a debug log (as option)
368  return fom;
369 }
370 
371 ////////////////////////////////////////////////////////////////////////////////
372 /// fill the private histograms with the mva distributinos for sig/bkg
373 
375 {
376  if (fMvaSig) fMvaSig->Delete();
377  if (fMvaBkg) fMvaBkg->Delete();
380 
381  // maybe later on this should be done a bit more clever (time consuming) by
382  // first determining proper ranges, removing outliers, as we do in the
383  // MVA output calculation in MethodBase::TestClassifier...
384  // --> then it might be possible also to use the splined PDF's which currently
385  // doesn't seem to work
386 
387  fMvaSig = new TH1D("fMvaSig","",100,-1.5,1.5); //used for spline fit
388  fMvaBkg = new TH1D("fMvaBkg","",100,-1.5,1.5); //used for spline fit
389  fMvaSigFineBin = new TH1D("fMvaSigFineBin","",100000,-1.5,1.5);
390  fMvaBkgFineBin = new TH1D("fMvaBkgFineBin","",100000,-1.5,1.5);
391 
392  const std::vector< Event*> events=fMethod->Data()->GetEventCollection(Types::kTesting);
393 
394  UInt_t signalClassNr = fMethod->DataInfo().GetClassInfo("Signal")->GetNumber();
395 
396  // fMethod->GetTransformationHandler().CalcTransformations(fMethod->Data()->GetEventCollection(Types::kTesting));
397 
398  for (UInt_t iev=0; iev < events.size() ; iev++){
399  // std::cout << " GetMVADists event " << iev << std::endl;
400  // std::cout << " Class = " << events[iev]->GetClass() << std::endl;
401  // std::cout << " MVA Value = " << fMethod->GetMvaValue(events[iev]) << std::endl;
402  if (events[iev]->GetClass() == signalClassNr) {
403  fMvaSig->Fill(fMethod->GetMvaValue(events[iev]),events[iev]->GetWeight());
404  fMvaSigFineBin->Fill(fMethod->GetMvaValue(events[iev]),events[iev]->GetWeight());
405  } else {
406  fMvaBkg->Fill(fMethod->GetMvaValue(events[iev]),events[iev]->GetWeight());
407  fMvaBkgFineBin->Fill(fMethod->GetMvaValue(events[iev]),events[iev]->GetWeight());
408  }
409  }
410 }
411 ////////////////////////////////////////////////////////////////////////////////
412 /// return the searation between the signal and background
413 /// MVA ouput distribution
414 
416 {
417  GetMVADists();
418  if (1){
419  PDF *splS = new PDF( " PDF Sig", fMvaSig, PDF::kSpline2 );
420  PDF *splB = new PDF( " PDF Bkg", fMvaBkg, PDF::kSpline2 );
421  return gTools().GetSeparation(*splS,*splB);
422  }else{
423  std::cout << "Separation caclulcaton via histograms (not PDFs) seems to give still strange results!! Don't do that, check!!"<<std::endl;
424  return gTools().GetSeparation(fMvaSigFineBin,fMvaBkgFineBin); // somehow sitll gives strange results!!!! Check!!!
425  }
426 }
427 
428 
429 ////////////////////////////////////////////////////////////////////////////////
430 /// calculate the area (integral) under the ROC curve as a
431 /// overall quality measure of the classification
432 ///
433 /// makeing pdfs out of the MVA-ouput distributions doesn't work
434 /// reliably for cases where the MVA-ouput isn't a smooth distribution.
435 /// this happens "frequently" in BDTs for example when the number of
436 /// trees is small resulting in only some discrete possible MVA ouput values.
437 /// (I still leave the code here, but use this with care!!! The default
438 /// however is to use the distributions!!!
439 
441 {
442  GetMVADists();
443 
444  Double_t integral = 0;
445  if (0){
446  PDF *pdfS = new PDF( " PDF Sig", fMvaSig, PDF::kSpline2 );
447  PDF *pdfB = new PDF( " PDF Bkg", fMvaBkg, PDF::kSpline2 );
448 
449  Double_t xmin = TMath::Min(pdfS->GetXmin(), pdfB->GetXmin());
450  Double_t xmax = TMath::Max(pdfS->GetXmax(), pdfB->GetXmax());
451 
452  UInt_t nsteps = 1000;
453  Double_t step = (xmax-xmin)/Double_t(nsteps);
454  Double_t cut = xmin;
455  for (UInt_t i=0; i<nsteps; i++){
456  integral += (1-pdfB->GetIntegral(cut,xmax)) * pdfS->GetVal(cut);
457  cut+=step;
458  }
459  integral*=step;
460  }else{
461  // sanity checks
464  std::cout << " Error in OptimizeConfigParameters GetROCIntegral, unequal histograms for sig and bkg.." << std::endl;
465  std::exit(1);
466  }else{
467 
468  Double_t *cumulator = fMvaBkgFineBin->GetIntegral();
470  // get the true signal integral (CompuetIntegral just return 1 as they
471  // automatically normalize. IN ADDITION, they do not account for variable
472  // bin sizes (which you migh perhaps use later on for the fMvaSig/Bkg histograms)
473  Double_t sigIntegral = 0;
474  for (Int_t ibin=1; ibin<=nbins; ibin++){
475  sigIntegral += fMvaSigFineBin->GetBinContent(ibin) * fMvaSigFineBin->GetBinWidth(ibin);
476  }
477  //gTools().NormHist( fMvaSigFineBin ); // also doesn't use variable bin width. And callse TH1::Scale, which oddly enough does not change the SumOfWeights !!!
478 
479  for (Int_t ibin=1; ibin <= nbins; ibin++){ // don't include under- and overflow bin
480  integral += (cumulator[ibin]) * fMvaSigFineBin->GetBinContent(ibin)/sigIntegral * fMvaSigFineBin->GetBinWidth(ibin) ;
481  }
482  }
483  }
484 
485  return integral;
486 }
487 
488 
489 ////////////////////////////////////////////////////////////////////////////////
490 /// calculate the signal efficiency for a given background efficiency
491 
493 {
494  GetMVADists();
495  Double_t sigEff=0;
496 
497  // sanity checks
500  std::cout << " Error in OptimizeConfigParameters GetSigEffAt, unequal histograms for sig and bkg.." << std::endl;
501  std::exit(1);
502  }else{
503  Double_t *bkgCumulator = fMvaBkgFineBin->GetIntegral();
504  Double_t *sigCumulator = fMvaSigFineBin->GetIntegral();
505 
507  Int_t ibin=0;
508 
509  // std::cout << " bkgIntegral="<<bkgIntegral
510  // << " sigIntegral="<<sigIntegral
511  // << " bkgCumulator[nbins]="<<bkgCumulator[nbins]
512  // << " sigCumulator[nbins]="<<sigCumulator[nbins]
513  // << std::endl;
514 
515  while (bkgCumulator[nbins-ibin] > (1-bkgEff)) {
516  sigEff = sigCumulator[nbins]-sigCumulator[nbins-ibin];
517  ibin++;
518  }
519  }
520  return sigEff;
521 }
522 
523 
524 //__adaptated_by_marc-olivier.bettler@cern.ch__________________________
525 ////////////////////////////////////////////////////////////////////////////////
526 /// calculate the background efficiency for a given signal efficiency
527 
529 {
530  GetMVADists();
531  Double_t bkgEff=0;
532 
533  // sanity checks
536  std::cout << " Error in OptimizeConfigParameters GetBkgEffAt, unequal histograms for sig and bkg.." << std::endl;
537  std::exit(1);
538  }else{
539 
540  Double_t *bkgCumulator = fMvaBkgFineBin->GetIntegral();
541  Double_t *sigCumulator = fMvaSigFineBin->GetIntegral();
542 
544  Int_t ibin=0;
545 
546  // std::cout << " bkgIntegral="<<bkgIntegral
547  // << " sigIntegral="<<sigIntegral
548  // << " bkgCumulator[nbins]="<<bkgCumulator[nbins]
549  // << " sigCumulator[nbins]="<<sigCumulator[nbins]
550  // << std::endl;
551 
552  while ( sigCumulator[nbins]-sigCumulator[nbins-ibin] < sigEff) {
553  bkgEff = bkgCumulator[nbins]-bkgCumulator[nbins-ibin];
554  ibin++;
555  }
556  }
557  return bkgEff;
558 }
559 
560 //__adaptated_by_marc-olivier.bettler@cern.ch__________________________
561 ////////////////////////////////////////////////////////////////////////////////
562 /// calculate the background rejection for a given signal efficiency
563 
565 {
566  GetMVADists();
567  Double_t bkgRej=0;
568 
569  // sanity checks
572  std::cout << " Error in OptimizeConfigParameters GetBkgEffAt, unequal histograms for sig and bkg.." << std::endl;
573  std::exit(1);
574  }else{
575 
576  Double_t *bkgCumulator = fMvaBkgFineBin->GetIntegral();
577  Double_t *sigCumulator = fMvaSigFineBin->GetIntegral();
578 
580  Int_t ibin=0;
581 
582  // std::cout << " bkgIntegral="<<bkgIntegral
583  // << " sigIntegral="<<sigIntegral
584  // << " bkgCumulator[nbins]="<<bkgCumulator[nbins]
585  // << " sigCumulator[nbins]="<<sigCumulator[nbins]
586  // << std::endl;
587 
588  while ( sigCumulator[nbins]-sigCumulator[nbins-ibin] < sigEff) {
589  bkgRej = bkgCumulator[nbins-ibin];
590  ibin++;
591  }
592  }
593  return bkgRej;
594 }
std::map< TString, Double_t > fTunedParameters
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TObject.cxx:830
std::string GetName(const std::string &scope_name)
Definition: Cppyy.cxx:140
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3125
virtual void SetTuneParameters(std::map< TString, Double_t > tuneParameters)
set the tuning parameters accoding to the argument This is just a dummy .
Definition: MethodBase.cxx:638
std::map< TString, TMVA::Interval * > fTuneParameters
float xmin
Definition: THbookFile.cxx:93
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
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
float ymin
Definition: THbookFile.cxx:93
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition: TNamed.cxx:131
TH1 * h
Definition: legend2.C:5
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4638
Basic string class.
Definition: TString.h:137
Double_t GetSeparation()
return the searation between the signal and background MVA ouput distribution
TransformationHandler & GetTransformationHandler(Bool_t takeReroutedIfAvailable=true)
Definition: MethodBase.h:390
void optimizeScan()
do the actual optimization using a simple scan method, i.e.
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
int Int_t
Definition: RtypesCore.h:41
virtual void SetYTitle(const char *title)
Definition: TH1.h:414
const Bool_t kFALSE
Definition: Rtypes.h:92
int nbins[3]
const std::vector< Event * > & GetEventCollection(Types::ETreeType type=Types::kMaxTreeType) const
Definition: DataSet.h:239
STL namespace.
Double_t GetFOM()
Return the Figure of Merit (FOM) used in the parameter optimization process.
static void SetIsTraining(Bool_t)
when this static function is called, it sets the flag whether events with negative event weight shoul...
Definition: Event.cxx:388
TClass * GetClass(T *)
Definition: TClass.h:555
Tools & gTools()
Definition: Tools.cxx:79
Double_t GetXmin() const
Definition: TAxis.h:139
void GetMVADists()
fill the private histograms with the mva distributinos for sig/bkg
Double_t x[n]
Definition: legend1.C:17
Double_t Run()
estimator function interface for fitting
Definition: FitterBase.cxx:80
DataSet * Data() const
Definition: MethodBase.h:405
std::vector< std::vector< double > > Data
DataSetInfo & DataInfo() const
Definition: MethodBase.h:406
Bool_t DoRegression() const
Definition: MethodBase.h:434
Definition: PDF.h:71
TCppMethod_t GetMethod(TCppScope_t scope, TCppIndex_t imeth)
Definition: Cppyy.cxx:722
const std::vector< Event * > * CalcTransformations(const std::vector< Event *> &, Bool_t createNewVector=kFALSE)
computation of transformation
virtual void Delete(Option_t *option="")
Delete this object.
Definition: TObject.cxx:229
virtual Double_t * GetIntegral()
Return a pointer to the array of bins integral.
Definition: TH1.cxx:2403
float ymax
Definition: THbookFile.cxx:93
virtual void Train()=0
Double_t GetROCIntegral()
calculate the area (integral) under the ROC curve as a overall quality measure of the classification ...
SVector< double, 2 > v
Definition: Dict.h:5
const char * GetName() const
Definition: MethodBase.h:330
ClassInfo * GetClassInfo(Int_t clNum) const
std::map< TString, Double_t > optimize()
virtual ~OptimizeConfigParameters()
the destructor (delete the OptimizeConfigParameters, store the graph and .. delete it) ...
unsigned int UInt_t
Definition: RtypesCore.h:42
Double_t GetSigEffAtBkgEff(Double_t bkgEff=0.1)
calculate the signal efficiency for a given background efficiency
double floor(double)
Double_t GetXmin() const
Definition: PDF.h:112
float xmax
Definition: THbookFile.cxx:93
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width for 1D histogram.
Definition: TH1.cxx:8273
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
Double_t GetXmax() const
Definition: PDF.h:113
Double_t y[n]
Definition: legend1.C:17
virtual void Reset()
Definition: MethodBase.h:197
Double_t GetBkgEffAtSigEff(Double_t sigEff=0.5)
calculate the background efficiency for a given signal efficiency
THist< 2, double, THistStatContent, THistStatUncertainty > TH2D
Definition: THist.hxx:307
UInt_t GetNumber() const
Definition: ClassInfo.h:73
Abstract ClassifierFactory template that handles arbitrary types.
virtual void SetXTitle(const char *title)
Definition: TH1.h:413
virtual Bool_t cd(const char *path=0)
Change current directory to "this" directory.
Definition: TDirectory.cxx:435
std::vector< int > GetScanIndices(int val, std::vector< int > base)
helper function to scan through the all the combinations in the parameter space
Double_t GetBkgRejAtSigEff(Double_t sigEff=0.5)
calculate the background rejection for a given signal efficiency
TDirectory * BaseDir() const
returns the ROOT directory where info/histograms etc of the corresponding MVA method instance are sto...
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:53
Double_t EstimatorFunction(std::vector< Double_t > &)
return the estimator (from current FOM) for the fitting interface
#define NULL
Definition: Rtypes.h:82
std::map< std::vector< Double_t >, Double_t > fAlreadyTrainedParCombination
THist< 1, double, THistStatContent, THistStatUncertainty > TH1D
Definition: THist.hxx:301
Double_t GetIntegral(Double_t xmin, Double_t xmax)
computes PDF integral within given ranges
Definition: PDF.cxx:652
virtual Int_t GetNbinsX() const
Definition: TH1.h:301
const Bool_t kTRUE
Definition: Rtypes.h:91
void CheckForUnusedOptions() const
checks for unused options in option string
const Int_t n
Definition: legend1.C:16
char name[80]
Definition: TGX11.cxx:109
TAxis * GetXaxis()
Definition: TH1.h:324
Double_t GetVal(Double_t x) const
returns value PDF(x)
Definition: PDF.cxx:699
Double_t GetSeparation(TH1 *S, TH1 *B) const
compute "separation" defined as <s2> = (1/2) Int_-oo..+oo { (S(x) - B(x))^2/(S(x) + B(x)) dx } ...
Definition: Tools.cxx:136
tomato 2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:296