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