Logo ROOT   6.08/07
Reference Guide
TMVAClassificationApplication.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_tmva
3 /// \notebook -nodraw
4 /// This macro provides a simple example on how to use the trained classifiers
5 /// within an analysis module
6 /// - Project : TMVA - a Root-integrated toolkit for multivariate data analysis
7 /// - Package : TMVA
8 /// - Exectuable: TMVAClassificationApplication
9 ///
10 /// \macro_output
11 /// \macro_code
12 /// \author Andreas Hoecker
13 
14 #include <cstdlib>
15 #include <vector>
16 #include <iostream>
17 #include <map>
18 #include <string>
19 
20 #include "TFile.h"
21 #include "TTree.h"
22 #include "TString.h"
23 #include "TSystem.h"
24 #include "TROOT.h"
25 #include "TStopwatch.h"
26 
27 #include "TMVA/Tools.h"
28 #include "TMVA/Reader.h"
29 #include "TMVA/MethodCuts.h"
30 
31 using namespace TMVA;
32 
33 void TMVAClassificationApplication( TString myMethodList = "" )
34 {
35 
36  //---------------------------------------------------------------
37  // This loads the library
39 
40  // Default MVA methods to be trained + tested
41  std::map<std::string,int> Use;
42 
43  // Cut optimisation
44  Use["Cuts"] = 1;
45  Use["CutsD"] = 1;
46  Use["CutsPCA"] = 0;
47  Use["CutsGA"] = 0;
48  Use["CutsSA"] = 0;
49  //
50  // 1-dimensional likelihood ("naive Bayes estimator")
51  Use["Likelihood"] = 1;
52  Use["LikelihoodD"] = 0; // the "D" extension indicates decorrelated input variables (see option strings)
53  Use["LikelihoodPCA"] = 1; // the "PCA" extension indicates PCA-transformed input variables (see option strings)
54  Use["LikelihoodKDE"] = 0;
55  Use["LikelihoodMIX"] = 0;
56  //
57  // Mutidimensional likelihood and Nearest-Neighbour methods
58  Use["PDERS"] = 1;
59  Use["PDERSD"] = 0;
60  Use["PDERSPCA"] = 0;
61  Use["PDEFoam"] = 1;
62  Use["PDEFoamBoost"] = 0; // uses generalised MVA method boosting
63  Use["KNN"] = 1; // k-nearest neighbour method
64  //
65  // Linear Discriminant Analysis
66  Use["LD"] = 1; // Linear Discriminant identical to Fisher
67  Use["Fisher"] = 0;
68  Use["FisherG"] = 0;
69  Use["BoostedFisher"] = 0; // uses generalised MVA method boosting
70  Use["HMatrix"] = 0;
71  //
72  // Function Discriminant analysis
73  Use["FDA_GA"] = 1; // minimisation of user-defined function using Genetics Algorithm
74  Use["FDA_SA"] = 0;
75  Use["FDA_MC"] = 0;
76  Use["FDA_MT"] = 0;
77  Use["FDA_GAMT"] = 0;
78  Use["FDA_MCMT"] = 0;
79  //
80  // Neural Networks (all are feed-forward Multilayer Perceptrons)
81  Use["MLP"] = 0; // Recommended ANN
82  Use["MLPBFGS"] = 0; // Recommended ANN with optional training method
83  Use["MLPBNN"] = 1; // Recommended ANN with BFGS training method and bayesian regulator
84  Use["CFMlpANN"] = 0; // Depreciated ANN from ALEPH
85  Use["TMlpANN"] = 0; // ROOT's own ANN
86  Use["DNN"] = 0; // improved implementation of a NN
87  //
88  // Support Vector Machine
89  Use["SVM"] = 1;
90  //
91  // Boosted Decision Trees
92  Use["BDT"] = 1; // uses Adaptive Boost
93  Use["BDTG"] = 0; // uses Gradient Boost
94  Use["BDTB"] = 0; // uses Bagging
95  Use["BDTD"] = 0; // decorrelation + Adaptive Boost
96  Use["BDTF"] = 0; // allow usage of fisher discriminant for node splitting
97  //
98  // Friedman's RuleFit method, ie, an optimised series of cuts ("rules")
99  Use["RuleFit"] = 1;
100  // ---------------------------------------------------------------
101  Use["Plugin"] = 0;
102  Use["Category"] = 0;
103  Use["SVM_Gauss"] = 0;
104  Use["SVM_Poly"] = 0;
105  Use["SVM_Lin"] = 0;
106 
107  std::cout << std::endl;
108  std::cout << "==> Start TMVAClassificationApplication" << std::endl;
109 
110  // Select methods (don't look at this code - not of interest)
111  if (myMethodList != "") {
112  for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) it->second = 0;
113 
114  std::vector<TString> mlist = gTools().SplitString( myMethodList, ',' );
115  for (UInt_t i=0; i<mlist.size(); i++) {
116  std::string regMethod(mlist[i]);
117 
118  if (Use.find(regMethod) == Use.end()) {
119  std::cout << "Method \"" << regMethod
120  << "\" not known in TMVA under this name. Choose among the following:" << std::endl;
121  for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) {
122  std::cout << it->first << " ";
123  }
124  std::cout << std::endl;
125  return;
126  }
127  Use[regMethod] = 1;
128  }
129  }
130 
131  // --------------------------------------------------------------------------------------------------
132 
133  // Create the Reader object
134 
135  TMVA::Reader *reader = new TMVA::Reader( "!Color:!Silent" );
136 
137  // Create a set of variables and declare them to the reader
138  // - the variable names MUST corresponds in name and type to those given in the weight file(s) used
139  Float_t var1, var2;
140  Float_t var3, var4;
141  reader->AddVariable( "myvar1 := var1+var2", &var1 );
142  reader->AddVariable( "myvar2 := var1-var2", &var2 );
143  reader->AddVariable( "var3", &var3 );
144  reader->AddVariable( "var4", &var4 );
145 
146  // Spectator variables declared in the training have to be added to the reader, too
147  Float_t spec1,spec2;
148  reader->AddSpectator( "spec1 := var1*2", &spec1 );
149  reader->AddSpectator( "spec2 := var1*3", &spec2 );
150 
151  Float_t Category_cat1, Category_cat2, Category_cat3;
152  if (Use["Category"]){
153  // Add artificial spectators for distinguishing categories
154  reader->AddSpectator( "Category_cat1 := var3<=0", &Category_cat1 );
155  reader->AddSpectator( "Category_cat2 := (var3>0)&&(var4<0)", &Category_cat2 );
156  reader->AddSpectator( "Category_cat3 := (var3>0)&&(var4>=0)", &Category_cat3 );
157  }
158 
159  // Book the MVA methods
160 
161  TString dir = "dataset/weights/";
162  TString prefix = "TMVAClassification";
163 
164  // Book method(s)
165  for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) {
166  if (it->second) {
167  TString methodName = TString(it->first) + TString(" method");
168  TString weightfile = dir + prefix + TString("_") + TString(it->first) + TString(".weights.xml");
169  reader->BookMVA( methodName, weightfile );
170  }
171  }
172 
173  // Book output histograms
174  UInt_t nbin = 100;
175  TH1F *histLk(0), *histLkD(0), *histLkPCA(0), *histLkKDE(0), *histLkMIX(0), *histPD(0), *histPDD(0);
176  TH1F *histPDPCA(0), *histPDEFoam(0), *histPDEFoamErr(0), *histPDEFoamSig(0), *histKNN(0), *histHm(0);
177  TH1F *histFi(0), *histFiG(0), *histFiB(0), *histLD(0), *histNn(0),*histNnbfgs(0),*histNnbnn(0);
178  TH1F *histNnC(0), *histNnT(0), *histNdn(0), *histBdt(0), *histBdtG(0), *histBdtB(0), *histBdtD(0);
179  TH1F *histBdtF(0), *histRf(0), *histSVMG(0), *histSVMP(0), *histSVML(0), *histFDAMT(0), *histFDAGA(0);
180  TH1F *histCat(0), *histPBdt(0);
181 
182  if (Use["Likelihood"]) histLk = new TH1F( "MVA_Likelihood", "MVA_Likelihood", nbin, -1, 1 );
183  if (Use["LikelihoodD"]) histLkD = new TH1F( "MVA_LikelihoodD", "MVA_LikelihoodD", nbin, -1, 0.9999 );
184  if (Use["LikelihoodPCA"]) histLkPCA = new TH1F( "MVA_LikelihoodPCA", "MVA_LikelihoodPCA", nbin, -1, 1 );
185  if (Use["LikelihoodKDE"]) histLkKDE = new TH1F( "MVA_LikelihoodKDE", "MVA_LikelihoodKDE", nbin, -0.00001, 0.99999 );
186  if (Use["LikelihoodMIX"]) histLkMIX = new TH1F( "MVA_LikelihoodMIX", "MVA_LikelihoodMIX", nbin, 0, 1 );
187  if (Use["PDERS"]) histPD = new TH1F( "MVA_PDERS", "MVA_PDERS", nbin, 0, 1 );
188  if (Use["PDERSD"]) histPDD = new TH1F( "MVA_PDERSD", "MVA_PDERSD", nbin, 0, 1 );
189  if (Use["PDERSPCA"]) histPDPCA = new TH1F( "MVA_PDERSPCA", "MVA_PDERSPCA", nbin, 0, 1 );
190  if (Use["KNN"]) histKNN = new TH1F( "MVA_KNN", "MVA_KNN", nbin, 0, 1 );
191  if (Use["HMatrix"]) histHm = new TH1F( "MVA_HMatrix", "MVA_HMatrix", nbin, -0.95, 1.55 );
192  if (Use["Fisher"]) histFi = new TH1F( "MVA_Fisher", "MVA_Fisher", nbin, -4, 4 );
193  if (Use["FisherG"]) histFiG = new TH1F( "MVA_FisherG", "MVA_FisherG", nbin, -1, 1 );
194  if (Use["BoostedFisher"]) histFiB = new TH1F( "MVA_BoostedFisher", "MVA_BoostedFisher", nbin, -2, 2 );
195  if (Use["LD"]) histLD = new TH1F( "MVA_LD", "MVA_LD", nbin, -2, 2 );
196  if (Use["MLP"]) histNn = new TH1F( "MVA_MLP", "MVA_MLP", nbin, -1.25, 1.5 );
197  if (Use["MLPBFGS"]) histNnbfgs = new TH1F( "MVA_MLPBFGS", "MVA_MLPBFGS", nbin, -1.25, 1.5 );
198  if (Use["MLPBNN"]) histNnbnn = new TH1F( "MVA_MLPBNN", "MVA_MLPBNN", nbin, -1.25, 1.5 );
199  if (Use["CFMlpANN"]) histNnC = new TH1F( "MVA_CFMlpANN", "MVA_CFMlpANN", nbin, 0, 1 );
200  if (Use["TMlpANN"]) histNnT = new TH1F( "MVA_TMlpANN", "MVA_TMlpANN", nbin, -1.3, 1.3 );
201  if (Use["DNN"]) histNdn = new TH1F( "MVA_DNN", "MVA_DNN", nbin, -0.1, 1.1 );
202  if (Use["BDT"]) histBdt = new TH1F( "MVA_BDT", "MVA_BDT", nbin, -0.8, 0.8 );
203  if (Use["BDTG"]) histBdtG = new TH1F( "MVA_BDTG", "MVA_BDTG", nbin, -1.0, 1.0 );
204  if (Use["BDTB"]) histBdtB = new TH1F( "MVA_BDTB", "MVA_BDTB", nbin, -1.0, 1.0 );
205  if (Use["BDTD"]) histBdtD = new TH1F( "MVA_BDTD", "MVA_BDTD", nbin, -0.8, 0.8 );
206  if (Use["BDTF"]) histBdtF = new TH1F( "MVA_BDTF", "MVA_BDTF", nbin, -1.0, 1.0 );
207  if (Use["RuleFit"]) histRf = new TH1F( "MVA_RuleFit", "MVA_RuleFit", nbin, -2.0, 2.0 );
208  if (Use["SVM_Gauss"]) histSVMG = new TH1F( "MVA_SVM_Gauss", "MVA_SVM_Gauss", nbin, 0.0, 1.0 );
209  if (Use["SVM_Poly"]) histSVMP = new TH1F( "MVA_SVM_Poly", "MVA_SVM_Poly", nbin, 0.0, 1.0 );
210  if (Use["SVM_Lin"]) histSVML = new TH1F( "MVA_SVM_Lin", "MVA_SVM_Lin", nbin, 0.0, 1.0 );
211  if (Use["FDA_MT"]) histFDAMT = new TH1F( "MVA_FDA_MT", "MVA_FDA_MT", nbin, -2.0, 3.0 );
212  if (Use["FDA_GA"]) histFDAGA = new TH1F( "MVA_FDA_GA", "MVA_FDA_GA", nbin, -2.0, 3.0 );
213  if (Use["Category"]) histCat = new TH1F( "MVA_Category", "MVA_Category", nbin, -2., 2. );
214  if (Use["Plugin"]) histPBdt = new TH1F( "MVA_PBDT", "MVA_BDT", nbin, -0.8, 0.8 );
215 
216  // PDEFoam also returns per-event error, fill in histogram, and also fill significance
217  if (Use["PDEFoam"]) {
218  histPDEFoam = new TH1F( "MVA_PDEFoam", "MVA_PDEFoam", nbin, 0, 1 );
219  histPDEFoamErr = new TH1F( "MVA_PDEFoamErr", "MVA_PDEFoam error", nbin, 0, 1 );
220  histPDEFoamSig = new TH1F( "MVA_PDEFoamSig", "MVA_PDEFoam significance", nbin, 0, 10 );
221  }
222 
223  // Book example histogram for probability (the other methods are done similarly)
224  TH1F *probHistFi(0), *rarityHistFi(0);
225  if (Use["Fisher"]) {
226  probHistFi = new TH1F( "MVA_Fisher_Proba", "MVA_Fisher_Proba", nbin, 0, 1 );
227  rarityHistFi = new TH1F( "MVA_Fisher_Rarity", "MVA_Fisher_Rarity", nbin, 0, 1 );
228  }
229 
230  // Prepare input tree (this must be replaced by your data source)
231  // in this example, there is a toy tree with signal and one with background events
232  // we'll later on use only the "signal" events for the test in this example.
233  //
234  TFile *input(0);
235  TString fname = "./tmva_example.root";
236  if (!gSystem->AccessPathName( fname ))
237  input = TFile::Open( fname ); // check if file in local directory exists
238  else
239  input = TFile::Open( "http://root.cern.ch/files/tmva_class_example.root" ); // if not: download from ROOT server
240 
241  if (!input) {
242  std::cout << "ERROR: could not open data file" << std::endl;
243  exit(1);
244  }
245  std::cout << "--- TMVAClassificationApp : Using input file: " << input->GetName() << std::endl;
246 
247  // Event loop
248 
249  // Prepare the event tree
250  // - Here the variable names have to corresponds to your tree
251  // - You can use the same variables as above which is slightly faster,
252  // but of course you can use different ones and copy the values inside the event loop
253  //
254  std::cout << "--- Select signal sample" << std::endl;
255  TTree* theTree = (TTree*)input->Get("TreeS");
256  Float_t userVar1, userVar2;
257  theTree->SetBranchAddress( "var1", &userVar1 );
258  theTree->SetBranchAddress( "var2", &userVar2 );
259  theTree->SetBranchAddress( "var3", &var3 );
260  theTree->SetBranchAddress( "var4", &var4 );
261 
262  // Efficiency calculator for cut method
263  Int_t nSelCutsGA = 0;
264  Double_t effS = 0.7;
265 
266  std::vector<Float_t> vecVar(4); // vector for EvaluateMVA tests
267 
268  std::cout << "--- Processing: " << theTree->GetEntries() << " events" << std::endl;
269  TStopwatch sw;
270  sw.Start();
271  for (Long64_t ievt=0; ievt<theTree->GetEntries();ievt++) {
272 
273  if (ievt%1000 == 0) std::cout << "--- ... Processing event: " << ievt << std::endl;
274 
275  theTree->GetEntry(ievt);
276 
277  var1 = userVar1 + userVar2;
278  var2 = userVar1 - userVar2;
279 
280  // Return the MVA outputs and fill into histograms
281 
282  if (Use["CutsGA"]) {
283  // Cuts is a special case: give the desired signal efficienciy
284  Bool_t passed = reader->EvaluateMVA( "CutsGA method", effS );
285  if (passed) nSelCutsGA++;
286  }
287 
288  if (Use["Likelihood" ]) histLk ->Fill( reader->EvaluateMVA( "Likelihood method" ) );
289  if (Use["LikelihoodD" ]) histLkD ->Fill( reader->EvaluateMVA( "LikelihoodD method" ) );
290  if (Use["LikelihoodPCA"]) histLkPCA ->Fill( reader->EvaluateMVA( "LikelihoodPCA method" ) );
291  if (Use["LikelihoodKDE"]) histLkKDE ->Fill( reader->EvaluateMVA( "LikelihoodKDE method" ) );
292  if (Use["LikelihoodMIX"]) histLkMIX ->Fill( reader->EvaluateMVA( "LikelihoodMIX method" ) );
293  if (Use["PDERS" ]) histPD ->Fill( reader->EvaluateMVA( "PDERS method" ) );
294  if (Use["PDERSD" ]) histPDD ->Fill( reader->EvaluateMVA( "PDERSD method" ) );
295  if (Use["PDERSPCA" ]) histPDPCA ->Fill( reader->EvaluateMVA( "PDERSPCA method" ) );
296  if (Use["KNN" ]) histKNN ->Fill( reader->EvaluateMVA( "KNN method" ) );
297  if (Use["HMatrix" ]) histHm ->Fill( reader->EvaluateMVA( "HMatrix method" ) );
298  if (Use["Fisher" ]) histFi ->Fill( reader->EvaluateMVA( "Fisher method" ) );
299  if (Use["FisherG" ]) histFiG ->Fill( reader->EvaluateMVA( "FisherG method" ) );
300  if (Use["BoostedFisher"]) histFiB ->Fill( reader->EvaluateMVA( "BoostedFisher method" ) );
301  if (Use["LD" ]) histLD ->Fill( reader->EvaluateMVA( "LD method" ) );
302  if (Use["MLP" ]) histNn ->Fill( reader->EvaluateMVA( "MLP method" ) );
303  if (Use["MLPBFGS" ]) histNnbfgs ->Fill( reader->EvaluateMVA( "MLPBFGS method" ) );
304  if (Use["MLPBNN" ]) histNnbnn ->Fill( reader->EvaluateMVA( "MLPBNN method" ) );
305  if (Use["CFMlpANN" ]) histNnC ->Fill( reader->EvaluateMVA( "CFMlpANN method" ) );
306  if (Use["TMlpANN" ]) histNnT ->Fill( reader->EvaluateMVA( "TMlpANN method" ) );
307  if (Use["DNN" ]) histNdn ->Fill( reader->EvaluateMVA( "DNN method" ) );
308  if (Use["BDT" ]) histBdt ->Fill( reader->EvaluateMVA( "BDT method" ) );
309  if (Use["BDTG" ]) histBdtG ->Fill( reader->EvaluateMVA( "BDTG method" ) );
310  if (Use["BDTB" ]) histBdtB ->Fill( reader->EvaluateMVA( "BDTB method" ) );
311  if (Use["BDTD" ]) histBdtD ->Fill( reader->EvaluateMVA( "BDTD method" ) );
312  if (Use["BDTF" ]) histBdtF ->Fill( reader->EvaluateMVA( "BDTF method" ) );
313  if (Use["RuleFit" ]) histRf ->Fill( reader->EvaluateMVA( "RuleFit method" ) );
314  if (Use["SVM_Gauss" ]) histSVMG ->Fill( reader->EvaluateMVA( "SVM_Gauss method" ) );
315  if (Use["SVM_Poly" ]) histSVMP ->Fill( reader->EvaluateMVA( "SVM_Poly method" ) );
316  if (Use["SVM_Lin" ]) histSVML ->Fill( reader->EvaluateMVA( "SVM_Lin method" ) );
317  if (Use["FDA_MT" ]) histFDAMT ->Fill( reader->EvaluateMVA( "FDA_MT method" ) );
318  if (Use["FDA_GA" ]) histFDAGA ->Fill( reader->EvaluateMVA( "FDA_GA method" ) );
319  if (Use["Category" ]) histCat ->Fill( reader->EvaluateMVA( "Category method" ) );
320  if (Use["Plugin" ]) histPBdt ->Fill( reader->EvaluateMVA( "P_BDT method" ) );
321 
322  // Retrieve also per-event error
323  if (Use["PDEFoam"]) {
324  Double_t val = reader->EvaluateMVA( "PDEFoam method" );
325  Double_t err = reader->GetMVAError();
326  histPDEFoam ->Fill( val );
327  histPDEFoamErr->Fill( err );
328  if (err>1.e-50) histPDEFoamSig->Fill( val/err );
329  }
330 
331  // Retrieve probability instead of MVA output
332  if (Use["Fisher"]) {
333  probHistFi ->Fill( reader->GetProba ( "Fisher method" ) );
334  rarityHistFi->Fill( reader->GetRarity( "Fisher method" ) );
335  }
336  }
337 
338  // Get elapsed time
339  sw.Stop();
340  std::cout << "--- End of event loop: "; sw.Print();
341 
342  // Get efficiency for cuts classifier
343  if (Use["CutsGA"]) std::cout << "--- Efficiency for CutsGA method: " << double(nSelCutsGA)/theTree->GetEntries()
344  << " (for a required signal efficiency of " << effS << ")" << std::endl;
345 
346  if (Use["CutsGA"]) {
347 
348  // test: retrieve cuts for particular signal efficiency
349  // CINT ignores dynamic_casts so we have to use a cuts-secific Reader function to acces the pointer
350  TMVA::MethodCuts* mcuts = reader->FindCutsMVA( "CutsGA method" ) ;
351 
352  if (mcuts) {
353  std::vector<Double_t> cutsMin;
354  std::vector<Double_t> cutsMax;
355  mcuts->GetCuts( 0.7, cutsMin, cutsMax );
356  std::cout << "--- -------------------------------------------------------------" << std::endl;
357  std::cout << "--- Retrieve cut values for signal efficiency of 0.7 from Reader" << std::endl;
358  for (UInt_t ivar=0; ivar<cutsMin.size(); ivar++) {
359  std::cout << "... Cut: "
360  << cutsMin[ivar]
361  << " < \""
362  << mcuts->GetInputVar(ivar)
363  << "\" <= "
364  << cutsMax[ivar] << std::endl;
365  }
366  std::cout << "--- -------------------------------------------------------------" << std::endl;
367  }
368  }
369 
370  // Write histograms
371 
372  TFile *target = new TFile( "TMVApp.root","RECREATE" );
373  if (Use["Likelihood" ]) histLk ->Write();
374  if (Use["LikelihoodD" ]) histLkD ->Write();
375  if (Use["LikelihoodPCA"]) histLkPCA ->Write();
376  if (Use["LikelihoodKDE"]) histLkKDE ->Write();
377  if (Use["LikelihoodMIX"]) histLkMIX ->Write();
378  if (Use["PDERS" ]) histPD ->Write();
379  if (Use["PDERSD" ]) histPDD ->Write();
380  if (Use["PDERSPCA" ]) histPDPCA ->Write();
381  if (Use["KNN" ]) histKNN ->Write();
382  if (Use["HMatrix" ]) histHm ->Write();
383  if (Use["Fisher" ]) histFi ->Write();
384  if (Use["FisherG" ]) histFiG ->Write();
385  if (Use["BoostedFisher"]) histFiB ->Write();
386  if (Use["LD" ]) histLD ->Write();
387  if (Use["MLP" ]) histNn ->Write();
388  if (Use["MLPBFGS" ]) histNnbfgs ->Write();
389  if (Use["MLPBNN" ]) histNnbnn ->Write();
390  if (Use["CFMlpANN" ]) histNnC ->Write();
391  if (Use["TMlpANN" ]) histNnT ->Write();
392  if (Use["DNN" ]) histNdn ->Write();
393  if (Use["BDT" ]) histBdt ->Write();
394  if (Use["BDTG" ]) histBdtG ->Write();
395  if (Use["BDTB" ]) histBdtB ->Write();
396  if (Use["BDTD" ]) histBdtD ->Write();
397  if (Use["BDTF" ]) histBdtF ->Write();
398  if (Use["RuleFit" ]) histRf ->Write();
399  if (Use["SVM_Gauss" ]) histSVMG ->Write();
400  if (Use["SVM_Poly" ]) histSVMP ->Write();
401  if (Use["SVM_Lin" ]) histSVML ->Write();
402  if (Use["FDA_MT" ]) histFDAMT ->Write();
403  if (Use["FDA_GA" ]) histFDAGA ->Write();
404  if (Use["Category" ]) histCat ->Write();
405  if (Use["Plugin" ]) histPBdt ->Write();
406 
407  // Write also error and significance histos
408  if (Use["PDEFoam"]) { histPDEFoam->Write(); histPDEFoamErr->Write(); histPDEFoamSig->Write(); }
409 
410  // Write also probability hists
411  if (Use["Fisher"]) { if (probHistFi != 0) probHistFi->Write(); if (rarityHistFi != 0) rarityHistFi->Write(); }
412  target->Close();
413 
414  std::cout << "--- Created root file: \"TMVApp.root\" containing the MVA output histograms" << std::endl;
415 
416  delete reader;
417 
418  std::cout << "==> TMVAClassificationApplication is done!" << std::endl << std::endl;
419 }
420 
421 int main( int argc, char** argv )
422 {
423  TString methodList;
424  for (int i=1; i<argc; i++) {
425  TString regMethod(argv[i]);
426  if(regMethod=="-b" || regMethod=="--batch") continue;
427  if (!methodList.IsNull()) methodList += TString(",");
428  methodList += regMethod;
429  }
430  TMVAClassificationApplication(methodList);
431  return 0;
432 }
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1266
static Tools & Instance()
Definition: Tools.cxx:80
long long Long64_t
Definition: RtypesCore.h:69
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:58
void Print(Option_t *option="") const
Print the real and cpu time passed between the start and stop events.
Definition: TStopwatch.cxx:219
float Float_t
Definition: RtypesCore.h:53
void AddVariable(const TString &expression, Float_t *)
Add a float variable or expression to the reader.
Definition: Reader.cxx:309
THist< 1, float, THistStatContent, THistStatUncertainty > TH1F
Definition: THist.hxx:302
Double_t GetRarity(const TString &methodTag, Double_t mvaVal=-9999999)
evaluates the MVA&#39;s rarity
Definition: Reader.cxx:764
tomato 1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:575
Double_t GetCuts(Double_t effS, std::vector< Double_t > &cutMin, std::vector< Double_t > &cutMax) const
retrieve cut values for given signal efficiency
Definition: MethodCuts.cxx:563
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Double_t GetMVAError() const
Definition: Reader.h:107
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=1, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3907
const TString & GetInputVar(Int_t i) const
Definition: MethodBase.h:345
Tools & gTools()
Definition: Tools.cxx:79
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:77
IMethod * BookMVA(const TString &methodTag, const TString &weightfile)
read method name from weight file
Definition: Reader.cxx:378
R__EXTERN TSystem * gSystem
Definition: TSystem.h:549
unsigned int UInt_t
Definition: RtypesCore.h:42
void AddSpectator(const TString &expression, Float_t *)
Add a float spectator or expression to the reader.
Definition: Reader.cxx:327
double Double_t
Definition: RtypesCore.h:55
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
Double_t EvaluateMVA(const std::vector< Float_t > &, const TString &methodTag, Double_t aux=0)
Evaluate a std::vector<float> of input data for a given method The parameter aux is obligatory for th...
Definition: Reader.cxx:486
Abstract ClassifierFactory template that handles arbitrary types.
std::vector< TString > SplitString(const TString &theOpt, const char separator) const
splits the option string at &#39;separator&#39; and fills the list &#39;splitV&#39; with the primitive strings ...
Definition: Tools.cxx:1207
int main(int argc, char **argv)
Double_t GetProba(const TString &methodTag, Double_t ap_sig=0.5, Double_t mvaVal=-9999999)
evaluates probability of MVA for given set of input variables
Definition: Reader.cxx:733
MethodCuts * FindCutsMVA(const TString &methodTag)
special function for Cuts to avoid dynamic_casts in ROOT macros, which are not properly handled by CI...
Definition: Reader.cxx:725
Stopwatch class.
Definition: TStopwatch.h:30