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