Logo ROOT   6.08/07
Reference Guide
TMVAClassification.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_tmva
3 /// \notebook -nodraw
4 /// This macro provides examples for the training and testing of the
5 /// TMVA classifiers.
6 ///
7 /// As input data is used a toy-MC sample consisting of four Gaussian-distributed
8 /// and linearly correlated input variables.
9 /// The methods to be used can be switched on and off by means of booleans, or
10 /// via the prompt command, for example:
11 ///
12 /// root -l ./TMVAClassification.C\(\"Fisher,Likelihood\"\)
13 ///
14 /// (note that the backslashes are mandatory)
15 /// If no method given, a default set of classifiers is used.
16 /// The output file "TMVA.root" can be analysed with the use of dedicated
17 /// macros (simply say: root -l <macro.C>), which can be conveniently
18 /// invoked through a GUI that will appear at the end of the run of this macro.
19 /// Launch the GUI via the command:
20 ///
21 /// root -l ./TMVAGui.C
22 ///
23 /// You can also compile and run the example with the following commands
24 ///
25 /// make
26 /// ./TMVAClassification <Methods>
27 ///
28 /// where: `<Methods> = "method1 method2"` are the TMVA classifier names
29 /// example:
30 ///
31 /// ./TMVAClassification Fisher LikelihoodPCA BDT
32 ///
33 /// If no method given, a default set is of classifiers is used
34 ///
35 /// - Project : TMVA - a ROOT-integrated toolkit for multivariate data analysis
36 /// - Package : TMVA
37 /// - Root Macro: TMVAClassification
38 ///
39 /// \macro_output
40 /// \macro_code
41 /// \author Andreas Hoecker
42 
43 
44 #include <cstdlib>
45 #include <iostream>
46 #include <map>
47 #include <string>
48 
49 #include "TChain.h"
50 #include "TFile.h"
51 #include "TTree.h"
52 #include "TString.h"
53 #include "TObjString.h"
54 #include "TSystem.h"
55 #include "TROOT.h"
56 
57 #include "TMVA/Factory.h"
58 #include "TMVA/DataLoader.h"
59 #include "TMVA/Tools.h"
60 #include "TMVA/TMVAGui.h"
61 
62 int TMVAClassification( TString myMethodList = "" )
63 {
64  // The explicit loading of the shared libTMVA is done in TMVAlogon.C, defined in .rootrc
65  // if you use your private .rootrc, or run from a different directory, please copy the
66  // corresponding lines from .rootrc
67 
68  // Methods to be processed can be given as an argument; use format:
69  //
70  // mylinux~> root -l TMVAClassification.C\(\"myMethod1,myMethod2,myMethod3\"\)
71 
72  //---------------------------------------------------------------
73  // This loads the library
75 
76  // Default MVA methods to be trained + tested
77  std::map<std::string,int> Use;
78 
79  // Cut optimisation
80  Use["Cuts"] = 1;
81  Use["CutsD"] = 1;
82  Use["CutsPCA"] = 0;
83  Use["CutsGA"] = 0;
84  Use["CutsSA"] = 0;
85  //
86  // 1-dimensional likelihood ("naive Bayes estimator")
87  Use["Likelihood"] = 1;
88  Use["LikelihoodD"] = 0; // the "D" extension indicates decorrelated input variables (see option strings)
89  Use["LikelihoodPCA"] = 1; // the "PCA" extension indicates PCA-transformed input variables (see option strings)
90  Use["LikelihoodKDE"] = 0;
91  Use["LikelihoodMIX"] = 0;
92  //
93  // Mutidimensional likelihood and Nearest-Neighbour methods
94  Use["PDERS"] = 1;
95  Use["PDERSD"] = 0;
96  Use["PDERSPCA"] = 0;
97  Use["PDEFoam"] = 1;
98  Use["PDEFoamBoost"] = 0; // uses generalised MVA method boosting
99  Use["KNN"] = 1; // k-nearest neighbour method
100  //
101  // Linear Discriminant Analysis
102  Use["LD"] = 1; // Linear Discriminant identical to Fisher
103  Use["Fisher"] = 0;
104  Use["FisherG"] = 0;
105  Use["BoostedFisher"] = 0; // uses generalised MVA method boosting
106  Use["HMatrix"] = 0;
107  //
108  // Function Discriminant analysis
109  Use["FDA_GA"] = 1; // minimisation of user-defined function using Genetics Algorithm
110  Use["FDA_SA"] = 0;
111  Use["FDA_MC"] = 0;
112  Use["FDA_MT"] = 0;
113  Use["FDA_GAMT"] = 0;
114  Use["FDA_MCMT"] = 0;
115  //
116  // Neural Networks (all are feed-forward Multilayer Perceptrons)
117  Use["MLP"] = 0; // Recommended ANN
118  Use["MLPBFGS"] = 0; // Recommended ANN with optional training method
119  Use["MLPBNN"] = 1; // Recommended ANN with BFGS training method and bayesian regulator
120  Use["CFMlpANN"] = 0; // Depreciated ANN from ALEPH
121  Use["TMlpANN"] = 0; // ROOT's own ANN
122  Use["DNN"] = 0; // Deep Neural Network
123  Use["DNN_GPU"] = 0; // CUDA-accelerated DNN training.
124  Use["DNN_CPU"] = 0; // Multi-core accelerated DNN.
125  //
126  // Support Vector Machine
127  Use["SVM"] = 1;
128  //
129  // Boosted Decision Trees
130  Use["BDT"] = 1; // uses Adaptive Boost
131  Use["BDTG"] = 0; // uses Gradient Boost
132  Use["BDTB"] = 0; // uses Bagging
133  Use["BDTD"] = 0; // decorrelation + Adaptive Boost
134  Use["BDTF"] = 0; // allow usage of fisher discriminant for node splitting
135  //
136  // Friedman's RuleFit method, ie, an optimised series of cuts ("rules")
137  Use["RuleFit"] = 1;
138  // ---------------------------------------------------------------
139 
140  std::cout << std::endl;
141  std::cout << "==> Start TMVAClassification" << std::endl;
142 
143  // Select methods (don't look at this code - not of interest)
144  if (myMethodList != "") {
145  for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) it->second = 0;
146 
147  std::vector<TString> mlist = TMVA::gTools().SplitString( myMethodList, ',' );
148  for (UInt_t i=0; i<mlist.size(); i++) {
149  std::string regMethod(mlist[i]);
150 
151  if (Use.find(regMethod) == Use.end()) {
152  std::cout << "Method \"" << regMethod << "\" not known in TMVA under this name. Choose among the following:" << std::endl;
153  for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) std::cout << it->first << " ";
154  std::cout << std::endl;
155  return 1;
156  }
157  Use[regMethod] = 1;
158  }
159  }
160 
161  // --------------------------------------------------------------------------------------------------
162 
163  // Here the preparation phase begins
164 
165  // Read training and test data
166  // (it is also possible to use ASCII format as input -> see TMVA Users Guide)
167  TString fname = "./tmva_class_example.root";
168 
169  if (gSystem->AccessPathName( fname )) // file does not exist in local directory
170  gSystem->Exec("curl -O http://root.cern.ch/files/tmva_class_example.root");
171 
172  TFile *input = TFile::Open( fname );
173 
174  std::cout << "--- TMVAClassification : Using input file: " << input->GetName() << std::endl;
175 
176  // Register the training and test trees
177 
178  TTree *signalTree = (TTree*)input->Get("TreeS");
179  TTree *background = (TTree*)input->Get("TreeB");
180 
181  // Create a ROOT output file where TMVA will store ntuples, histograms, etc.
182  TString outfileName( "TMVA.root" );
183  TFile* outputFile = TFile::Open( outfileName, "RECREATE" );
184 
185  // Create the factory object. Later you can choose the methods
186  // whose performance you'd like to investigate. The factory is
187  // the only TMVA object you have to interact with
188  //
189  // The first argument is the base of the name of all the
190  // weightfiles in the directory weight/
191  //
192  // The second argument is the output file for the training results
193  // All TMVA output can be suppressed by removing the "!" (not) in
194  // front of the "Silent" argument in the option string
195  TMVA::Factory *factory = new TMVA::Factory( "TMVAClassification", outputFile,
196  "!V:!Silent:Color:DrawProgressBar:Transformations=I;D;P;G,D:AnalysisType=Classification" );
197 
199  // If you wish to modify default settings
200  // (please check "src/Config.h" to see all available global options)
201  //
202  // (TMVA::gConfig().GetVariablePlotting()).fTimesRMS = 8.0;
203  // (TMVA::gConfig().GetIONames()).fWeightFileDir = "myWeightDirectory";
204 
205  // Define the input variables that shall be used for the MVA training
206  // note that you may also use variable expressions, such as: "3*var1/var2*abs(var3)"
207  // [all types of expressions that can also be parsed by TTree::Draw( "expression" )]
208  dataloader->AddVariable( "myvar1 := var1+var2", 'F' );
209  dataloader->AddVariable( "myvar2 := var1-var2", "Expression 2", "", 'F' );
210  dataloader->AddVariable( "var3", "Variable 3", "units", 'F' );
211  dataloader->AddVariable( "var4", "Variable 4", "units", 'F' );
212 
213  // You can add so-called "Spectator variables", which are not used in the MVA training,
214  // but will appear in the final "TestTree" produced by TMVA. This TestTree will contain the
215  // input variables, the response values of all trained MVAs, and the spectator variables
216 
217  dataloader->AddSpectator( "spec1 := var1*2", "Spectator 1", "units", 'F' );
218  dataloader->AddSpectator( "spec2 := var1*3", "Spectator 2", "units", 'F' );
219 
220 
221  // global event weights per tree (see below for setting event-wise weights)
222  Double_t signalWeight = 1.0;
223  Double_t backgroundWeight = 1.0;
224 
225  // You can add an arbitrary number of signal or background trees
226  dataloader->AddSignalTree ( signalTree, signalWeight );
227  dataloader->AddBackgroundTree( background, backgroundWeight );
228 
229  // To give different trees for training and testing, do as follows:
230  //
231  // dataloader->AddSignalTree( signalTrainingTree, signalTrainWeight, "Training" );
232  // dataloader->AddSignalTree( signalTestTree, signalTestWeight, "Test" );
233 
234  // Use the following code instead of the above two or four lines to add signal and background
235  // training and test events "by hand"
236  // NOTE that in this case one should not give expressions (such as "var1+var2") in the input
237  // variable definition, but simply compute the expression before adding the event
238  // ```cpp
239  // // --- begin ----------------------------------------------------------
240  // std::vector<Double_t> vars( 4 ); // vector has size of number of input variables
241  // Float_t treevars[4], weight;
242  //
243  // // Signal
244  // for (UInt_t ivar=0; ivar<4; ivar++) signalTree->SetBranchAddress( Form( "var%i", ivar+1 ), &(treevars[ivar]) );
245  // for (UInt_t i=0; i<signalTree->GetEntries(); i++) {
246  // signalTree->GetEntry(i);
247  // for (UInt_t ivar=0; ivar<4; ivar++) vars[ivar] = treevars[ivar];
248  // // add training and test events; here: first half is training, second is testing
249  // // note that the weight can also be event-wise
250  // if (i < signalTree->GetEntries()/2.0) dataloader->AddSignalTrainingEvent( vars, signalWeight );
251  // else dataloader->AddSignalTestEvent ( vars, signalWeight );
252  // }
253  //
254  // // Background (has event weights)
255  // background->SetBranchAddress( "weight", &weight );
256  // for (UInt_t ivar=0; ivar<4; ivar++) background->SetBranchAddress( Form( "var%i", ivar+1 ), &(treevars[ivar]) );
257  // for (UInt_t i=0; i<background->GetEntries(); i++) {
258  // background->GetEntry(i);
259  // for (UInt_t ivar=0; ivar<4; ivar++) vars[ivar] = treevars[ivar];
260  // // add training and test events; here: first half is training, second is testing
261  // // note that the weight can also be event-wise
262  // if (i < background->GetEntries()/2) dataloader->AddBackgroundTrainingEvent( vars, backgroundWeight*weight );
263  // else dataloader->AddBackgroundTestEvent ( vars, backgroundWeight*weight );
264  // }
265  // // --- end ------------------------------------------------------------
266  // ```
267  // End of tree registration
268 
269  // Set individual event weights (the variables must exist in the original TTree)
270  // - for signal : `dataloader->SetSignalWeightExpression ("weight1*weight2");`
271  // - for background: `dataloader->SetBackgroundWeightExpression("weight1*weight2");`
272  dataloader->SetBackgroundWeightExpression( "weight" );
273 
274  // Apply additional cuts on the signal and background samples (can be different)
275  TCut mycuts = ""; // for example: TCut mycuts = "abs(var1)<0.5 && abs(var2-0.5)<1";
276  TCut mycutb = ""; // for example: TCut mycutb = "abs(var1)<0.5";
277 
278  // Tell the dataloader how to use the training and testing events
279  //
280  // If no numbers of events are given, half of the events in the tree are used
281  // for training, and the other half for testing:
282  //
283  // dataloader->PrepareTrainingAndTestTree( mycut, "SplitMode=random:!V" );
284  //
285  // To also specify the number of testing events, use:
286  //
287  // dataloader->PrepareTrainingAndTestTree( mycut,
288  // "NSigTrain=3000:NBkgTrain=3000:NSigTest=3000:NBkgTest=3000:SplitMode=Random:!V" );
289  dataloader->PrepareTrainingAndTestTree( mycuts, mycutb,
290  "nTrain_Signal=1000:nTrain_Background=1000:SplitMode=Random:NormMode=NumEvents:!V" );
291 
292  // ### Book MVA methods
293  //
294  // Please lookup the various method configuration options in the corresponding cxx files, eg:
295  // src/MethoCuts.cxx, etc, or here: http://tmva.sourceforge.net/optionRef.html
296  // it is possible to preset ranges in the option string in which the cut optimisation should be done:
297  // "...:CutRangeMin[2]=-1:CutRangeMax[2]=1"...", where [2] is the third input variable
298 
299  // Cut optimisation
300  if (Use["Cuts"])
301  factory->BookMethod( dataloader, TMVA::Types::kCuts, "Cuts",
302  "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart" );
303 
304  if (Use["CutsD"])
305  factory->BookMethod( dataloader, TMVA::Types::kCuts, "CutsD",
306  "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart:VarTransform=Decorrelate" );
307 
308  if (Use["CutsPCA"])
309  factory->BookMethod( dataloader, TMVA::Types::kCuts, "CutsPCA",
310  "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart:VarTransform=PCA" );
311 
312  if (Use["CutsGA"])
313  factory->BookMethod( dataloader, TMVA::Types::kCuts, "CutsGA",
314  "H:!V:FitMethod=GA:CutRangeMin[0]=-10:CutRangeMax[0]=10:VarProp[1]=FMax:EffSel:Steps=30:Cycles=3:PopSize=400:SC_steps=10:SC_rate=5:SC_factor=0.95" );
315 
316  if (Use["CutsSA"])
317  factory->BookMethod( dataloader, TMVA::Types::kCuts, "CutsSA",
318  "!H:!V:FitMethod=SA:EffSel:MaxCalls=150000:KernelTemp=IncAdaptive:InitialTemp=1e+6:MinTemp=1e-6:Eps=1e-10:UseDefaultScale" );
319 
320  // Likelihood ("naive Bayes estimator")
321  if (Use["Likelihood"])
322  factory->BookMethod( dataloader, TMVA::Types::kLikelihood, "Likelihood",
323  "H:!V:TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmoothBkg[1]=10:NSmooth=1:NAvEvtPerBin=50" );
324 
325  // Decorrelated likelihood
326  if (Use["LikelihoodD"])
327  factory->BookMethod( dataloader, TMVA::Types::kLikelihood, "LikelihoodD",
328  "!H:!V:TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmooth=5:NAvEvtPerBin=50:VarTransform=Decorrelate" );
329 
330  // PCA-transformed likelihood
331  if (Use["LikelihoodPCA"])
332  factory->BookMethod( dataloader, TMVA::Types::kLikelihood, "LikelihoodPCA",
333  "!H:!V:!TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmooth=5:NAvEvtPerBin=50:VarTransform=PCA" );
334 
335  // Use a kernel density estimator to approximate the PDFs
336  if (Use["LikelihoodKDE"])
337  factory->BookMethod( dataloader, TMVA::Types::kLikelihood, "LikelihoodKDE",
338  "!H:!V:!TransformOutput:PDFInterpol=KDE:KDEtype=Gauss:KDEiter=Adaptive:KDEFineFactor=0.3:KDEborder=None:NAvEvtPerBin=50" );
339 
340  // Use a variable-dependent mix of splines and kernel density estimator
341  if (Use["LikelihoodMIX"])
342  factory->BookMethod( dataloader, TMVA::Types::kLikelihood, "LikelihoodMIX",
343  "!H:!V:!TransformOutput:PDFInterpolSig[0]=KDE:PDFInterpolBkg[0]=KDE:PDFInterpolSig[1]=KDE:PDFInterpolBkg[1]=KDE:PDFInterpolSig[2]=Spline2:PDFInterpolBkg[2]=Spline2:PDFInterpolSig[3]=Spline2:PDFInterpolBkg[3]=Spline2:KDEtype=Gauss:KDEiter=Nonadaptive:KDEborder=None:NAvEvtPerBin=50" );
344 
345  // Test the multi-dimensional probability density estimator
346  // here are the options strings for the MinMax and RMS methods, respectively:
347  //
348  // "!H:!V:VolumeRangeMode=MinMax:DeltaFrac=0.2:KernelEstimator=Gauss:GaussSigma=0.3" );
349  // "!H:!V:VolumeRangeMode=RMS:DeltaFrac=3:KernelEstimator=Gauss:GaussSigma=0.3" );
350  if (Use["PDERS"])
351  factory->BookMethod( dataloader, TMVA::Types::kPDERS, "PDERS",
352  "!H:!V:NormTree=T:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600" );
353 
354  if (Use["PDERSD"])
355  factory->BookMethod( dataloader, TMVA::Types::kPDERS, "PDERSD",
356  "!H:!V:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600:VarTransform=Decorrelate" );
357 
358  if (Use["PDERSPCA"])
359  factory->BookMethod( dataloader, TMVA::Types::kPDERS, "PDERSPCA",
360  "!H:!V:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600:VarTransform=PCA" );
361 
362  // Multi-dimensional likelihood estimator using self-adapting phase-space binning
363  if (Use["PDEFoam"])
364  factory->BookMethod( dataloader, TMVA::Types::kPDEFoam, "PDEFoam",
365  "!H:!V:SigBgSeparate=F:TailCut=0.001:VolFrac=0.0666:nActiveCells=500:nSampl=2000:nBin=5:Nmin=100:Kernel=None:Compress=T" );
366 
367  if (Use["PDEFoamBoost"])
368  factory->BookMethod( dataloader, TMVA::Types::kPDEFoam, "PDEFoamBoost",
369  "!H:!V:Boost_Num=30:Boost_Transform=linear:SigBgSeparate=F:MaxDepth=4:UseYesNoCell=T:DTLogic=MisClassificationError:FillFoamWithOrigWeights=F:TailCut=0:nActiveCells=500:nBin=20:Nmin=400:Kernel=None:Compress=T" );
370 
371  // K-Nearest Neighbour classifier (KNN)
372  if (Use["KNN"])
373  factory->BookMethod( dataloader, TMVA::Types::kKNN, "KNN",
374  "H:nkNN=20:ScaleFrac=0.8:SigmaFact=1.0:Kernel=Gaus:UseKernel=F:UseWeight=T:!Trim" );
375 
376  // H-Matrix (chi2-squared) method
377  if (Use["HMatrix"])
378  factory->BookMethod( dataloader, TMVA::Types::kHMatrix, "HMatrix", "!H:!V:VarTransform=None" );
379 
380  // Linear discriminant (same as Fisher discriminant)
381  if (Use["LD"])
382  factory->BookMethod( dataloader, TMVA::Types::kLD, "LD", "H:!V:VarTransform=None:CreateMVAPdfs:PDFInterpolMVAPdf=Spline2:NbinsMVAPdf=50:NsmoothMVAPdf=10" );
383 
384  // Fisher discriminant (same as LD)
385  if (Use["Fisher"])
386  factory->BookMethod( dataloader, TMVA::Types::kFisher, "Fisher", "H:!V:Fisher:VarTransform=None:CreateMVAPdfs:PDFInterpolMVAPdf=Spline2:NbinsMVAPdf=50:NsmoothMVAPdf=10" );
387 
388  // Fisher with Gauss-transformed input variables
389  if (Use["FisherG"])
390  factory->BookMethod( dataloader, TMVA::Types::kFisher, "FisherG", "H:!V:VarTransform=Gauss" );
391 
392  // Composite classifier: ensemble (tree) of boosted Fisher classifiers
393  if (Use["BoostedFisher"])
394  factory->BookMethod( dataloader, TMVA::Types::kFisher, "BoostedFisher",
395  "H:!V:Boost_Num=20:Boost_Transform=log:Boost_Type=AdaBoost:Boost_AdaBoostBeta=0.2:!Boost_DetailedMonitoring" );
396 
397  // Function discrimination analysis (FDA) -- test of various fitters - the recommended one is Minuit (or GA or SA)
398  if (Use["FDA_MC"])
399  factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_MC",
400  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=MC:SampleSize=100000:Sigma=0.1" );
401 
402  if (Use["FDA_GA"]) // can also use Simulated Annealing (SA) algorithm (see Cuts_SA options])
403  factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_GA",
404  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=GA:PopSize=300:Cycles=3:Steps=20:Trim=True:SaveBestGen=1" );
405 
406  if (Use["FDA_SA"]) // can also use Simulated Annealing (SA) algorithm (see Cuts_SA options])
407  factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_SA",
408  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=SA:MaxCalls=15000:KernelTemp=IncAdaptive:InitialTemp=1e+6:MinTemp=1e-6:Eps=1e-10:UseDefaultScale" );
409 
410  if (Use["FDA_MT"])
411  factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_MT",
412  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=MINUIT:ErrorLevel=1:PrintLevel=-1:FitStrategy=2:UseImprove:UseMinos:SetBatch" );
413 
414  if (Use["FDA_GAMT"])
415  factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_GAMT",
416  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=GA:Converger=MINUIT:ErrorLevel=1:PrintLevel=-1:FitStrategy=0:!UseImprove:!UseMinos:SetBatch:Cycles=1:PopSize=5:Steps=5:Trim" );
417 
418  if (Use["FDA_MCMT"])
419  factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_MCMT",
420  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=MC:Converger=MINUIT:ErrorLevel=1:PrintLevel=-1:FitStrategy=0:!UseImprove:!UseMinos:SetBatch:SampleSize=20" );
421 
422  // TMVA ANN: MLP (recommended ANN) -- all ANNs in TMVA are Multilayer Perceptrons
423  if (Use["MLP"])
424  factory->BookMethod( dataloader, TMVA::Types::kMLP, "MLP", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:!UseRegulator" );
425 
426  if (Use["MLPBFGS"])
427  factory->BookMethod( dataloader, TMVA::Types::kMLP, "MLPBFGS", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:TrainingMethod=BFGS:!UseRegulator" );
428 
429  if (Use["MLPBNN"])
430  factory->BookMethod( dataloader, TMVA::Types::kMLP, "MLPBNN", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:TrainingMethod=BFGS:UseRegulator" ); // BFGS training with bayesian regulators
431 
432 
433  // Multi-architecture DNN implementation.
434  if (Use["DNN"])
435  {
436  // General layout.
437  TString layoutString ("Layout=TANH|128,TANH|128,TANH|128,LINEAR");
438 
439  // Training strategies.
440  TString training0("LearningRate=1e-1,Momentum=0.9,Repetitions=1,"
441  "ConvergenceSteps=20,BatchSize=256,TestRepetitions=10,"
442  "WeightDecay=1e-4,Regularization=L2,"
443  "DropConfig=0.0+0.5+0.5+0.5, Multithreading=True");
444  TString training1("LearningRate=1e-2,Momentum=0.9,Repetitions=1,"
445  "ConvergenceSteps=20,BatchSize=256,TestRepetitions=10,"
446  "WeightDecay=1e-4,Regularization=L2,"
447  "DropConfig=0.0+0.0+0.0+0.0, Multithreading=True");
448  TString training2("LearningRate=1e-3,Momentum=0.0,Repetitions=1,"
449  "ConvergenceSteps=20,BatchSize=256,TestRepetitions=10,"
450  "WeightDecay=1e-4,Regularization=L2,"
451  "DropConfig=0.0+0.0+0.0+0.0, Multithreading=True");
452  TString trainingStrategyString ("TrainingStrategy=");
453  trainingStrategyString += training0 + "|" + training1 + "|" + training2;
454 
455  // General Options.
456  TString dnnOptions ("!H:V:ErrorStrategy=CROSSENTROPY:VarTransform=N:"
457  "WeightInitialization=XAVIERUNIFORM");
458  dnnOptions.Append (":"); dnnOptions.Append (layoutString);
459  dnnOptions.Append (":"); dnnOptions.Append (trainingStrategyString);
460 
461  // Standard implementation, no dependencies.
462  TString stdOptions = dnnOptions + ":Architecture=STANDARD";
463  factory->BookMethod(dataloader, TMVA::Types::kDNN, "DNN", stdOptions);
464 
465  // Cuda implementation.
466  if (Use["DNN_GPU"]) {
467  TString gpuOptions = dnnOptions + ":Architecture=GPU";
468  factory->BookMethod(dataloader, TMVA::Types::kDNN, "DNN GPU", gpuOptions);
469  }
470  // Multi-core CPU implementation.
471  if (Use["DNN_CPU"]) {
472  TString cpuOptions = dnnOptions + ":Architecture=CPU";
473  factory->BookMethod(dataloader, TMVA::Types::kDNN, "DNN CPU", cpuOptions);
474  }
475  }
476 
477  // CF(Clermont-Ferrand)ANN
478  if (Use["CFMlpANN"])
479  factory->BookMethod( dataloader, TMVA::Types::kCFMlpANN, "CFMlpANN", "!H:!V:NCycles=2000:HiddenLayers=N+1,N" ); // n_cycles:#nodes:#nodes:...
480 
481  // Tmlp(Root)ANN
482  if (Use["TMlpANN"])
483  factory->BookMethod( dataloader, TMVA::Types::kTMlpANN, "TMlpANN", "!H:!V:NCycles=200:HiddenLayers=N+1,N:LearningMethod=BFGS:ValidationFraction=0.3" ); // n_cycles:#nodes:#nodes:...
484 
485  // Support Vector Machine
486  if (Use["SVM"])
487  factory->BookMethod( dataloader, TMVA::Types::kSVM, "SVM", "Gamma=0.25:Tol=0.001:VarTransform=Norm" );
488 
489  // Boosted Decision Trees
490  if (Use["BDTG"]) // Gradient Boost
491  factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDTG",
492  "!H:!V:NTrees=1000:MinNodeSize=2.5%:BoostType=Grad:Shrinkage=0.10:UseBaggedBoost:BaggedSampleFraction=0.5:nCuts=20:MaxDepth=2" );
493 
494  if (Use["BDT"]) // Adaptive Boost
495  factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDT",
496  "!H:!V:NTrees=850:MinNodeSize=2.5%:MaxDepth=3:BoostType=AdaBoost:AdaBoostBeta=0.5:UseBaggedBoost:BaggedSampleFraction=0.5:SeparationType=GiniIndex:nCuts=20" );
497 
498  if (Use["BDTB"]) // Bagging
499  factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDTB",
500  "!H:!V:NTrees=400:BoostType=Bagging:SeparationType=GiniIndex:nCuts=20" );
501 
502  if (Use["BDTD"]) // Decorrelation + Adaptive Boost
503  factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDTD",
504  "!H:!V:NTrees=400:MinNodeSize=5%:MaxDepth=3:BoostType=AdaBoost:SeparationType=GiniIndex:nCuts=20:VarTransform=Decorrelate" );
505 
506  if (Use["BDTF"]) // Allow Using Fisher discriminant in node splitting for (strong) linearly correlated variables
507  factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDTF",
508  "!H:!V:NTrees=50:MinNodeSize=2.5%:UseFisherCuts:MaxDepth=3:BoostType=AdaBoost:AdaBoostBeta=0.5:SeparationType=GiniIndex:nCuts=20" );
509 
510  // RuleFit -- TMVA implementation of Friedman's method
511  if (Use["RuleFit"])
512  factory->BookMethod( dataloader, TMVA::Types::kRuleFit, "RuleFit",
513  "H:!V:RuleFitModule=RFTMVA:Model=ModRuleLinear:MinImp=0.001:RuleMinDist=0.001:NTrees=20:fEventsMin=0.01:fEventsMax=0.5:GDTau=-1.0:GDTauPrec=0.01:GDStep=0.01:GDNSteps=10000:GDErrScale=1.02" );
514 
515  // For an example of the category classifier usage, see: TMVAClassificationCategory
516  //
517  // --------------------------------------------------------------------------------------------------
518  // Now you can optimize the setting (configuration) of the MVAs using the set of training events
519  // STILL EXPERIMENTAL and only implemented for BDT's !
520  //
521  // factory->OptimizeAllMethods("SigEffAt001","Scan");
522  // factory->OptimizeAllMethods("ROCIntegral","FitGA");
523  //
524  // --------------------------------------------------------------------------------------------------
525 
526  // Now you can tell the factory to train, test, and evaluate the MVAs
527  //
528  // Train MVAs using the set of training events
529  factory->TrainAllMethods();
530 
531  // Evaluate all MVAs using the set of test events
532  factory->TestAllMethods();
533 
534  // Evaluate and compare performance of all configured MVAs
535  factory->EvaluateAllMethods();
536 
537  // --------------------------------------------------------------
538 
539  // Save the output
540  outputFile->Close();
541 
542  std::cout << "==> Wrote root file: " << outputFile->GetName() << std::endl;
543  std::cout << "==> TMVAClassification is done!" << std::endl;
544 
545  delete factory;
546  delete dataloader;
547  // Launch the GUI for the root macros
548  if (!gROOT->IsBatch()) TMVA::TMVAGui( outfileName );
549 
550  return 0;
551 }
552 
553 int main( int argc, char** argv )
554 {
555  // Select methods (don't look at this code - not of interest)
556  TString methodList;
557  for (int i=1; i<argc; i++) {
558  TString regMethod(argv[i]);
559  if(regMethod=="-b" || regMethod=="--batch") continue;
560  if (!methodList.IsNull()) methodList += TString(",");
561  methodList += regMethod;
562  }
563  return TMVAClassification(methodList);
564 }
void AddBackgroundTree(TTree *background, Double_t weight=1.0, Types::ETreeType treetype=Types::kMaxTreeType)
Definition: DataLoader.cxx:382
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
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
MethodBase * BookMethod(DataLoader *loader, TString theMethodName, TString methodTitle, TString theOption="")
Definition: Factory.cxx:337
void TMVAGui(const char *fName="TMVA.root", TString dataset="")
Definition: TMVAGui.cxx:52
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:50
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
#define gROOT
Definition: TROOT.h:364
Basic string class.
Definition: TString.h:137
void TrainAllMethods()
iterates through all booked methods and calls training
Definition: Factory.cxx:822
void AddVariable(const TString &expression, const TString &title, const TString &unit, char type='F', Double_t min=0, Double_t max=0)
Definition: DataLoader.cxx:456
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
Tools & gTools()
Definition: Tools.cxx:79
A specialized string object used for TTree selections.
Definition: TCut.h:27
R__EXTERN TSystem * gSystem
Definition: TSystem.h:549
void EvaluateAllMethods(void)
iterates over all MVAs that have been booked, and calls their evaluation methods
Definition: Factory.cxx:1059
void TestAllMethods()
Definition: Factory.cxx:958
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual Int_t Exec(const char *shellcmd)
Execute a command.
Definition: TSystem.cxx:658
void PrepareTrainingAndTestTree(const TCut &cut, const TString &splitOpt)
Definition: DataLoader.cxx:580
double Double_t
Definition: RtypesCore.h:55
void SetBackgroundWeightExpression(const TString &variable)
Definition: DataLoader.cxx:513
Bool_t IsNull() const
Definition: TString.h:387
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
A TTree object has a header with a name and a title.
Definition: TTree.h:98
void AddSignalTree(TTree *signal, Double_t weight=1.0, Types::ETreeType treetype=Types::kMaxTreeType)
Definition: DataLoader.cxx:353
int main(int argc, char **argv)
void AddSpectator(const TString &expression, const TString &title="", const TString &unit="", Double_t min=0, Double_t max=0)
Definition: DataLoader.cxx:484
virtual void Close(Option_t *option="")
Close a file.
Definition: TFile.cxx:904