Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
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 "TMVAC.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
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#ifdef R__HAS_TMVAGPU
123 Use["DNN_GPU"] = 1; // CUDA-accelerated DNN training.
124#else
125 Use["DNN_GPU"] = 0;
126#endif
127
128#ifdef R__HAS_TMVACPU
129 Use["DNN_CPU"] = 1; // Multi-core accelerated DNN.
130#else
131 Use["DNN_CPU"] = 0;
132#endif
133 //
134 // Support Vector Machine
135 Use["SVM"] = 1;
136 //
137 // Boosted Decision Trees
138 Use["BDT"] = 1; // uses Adaptive Boost
139 Use["BDTG"] = 0; // uses Gradient Boost
140 Use["BDTB"] = 0; // uses Bagging
141 Use["BDTD"] = 0; // decorrelation + Adaptive Boost
142 Use["BDTF"] = 0; // allow usage of fisher discriminant for node splitting
143 //
144 // Friedman's RuleFit method, ie, an optimised series of cuts ("rules")
145 Use["RuleFit"] = 1;
146 // ---------------------------------------------------------------
147
148 std::cout << std::endl;
149 std::cout << "==> Start TMVAClassification" << std::endl;
150
151 // Select methods (don't look at this code - not of interest)
152 if (myMethodList != "") {
153 for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) it->second = 0;
154
155 std::vector<TString> mlist = TMVA::gTools().SplitString( myMethodList, ',' );
156 for (UInt_t i=0; i<mlist.size(); i++) {
157 std::string regMethod(mlist[i]);
158
159 if (Use.find(regMethod) == Use.end()) {
160 std::cout << "Method \"" << regMethod << "\" not known in TMVA under this name. Choose among the following:" << std::endl;
161 for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) std::cout << it->first << " ";
162 std::cout << std::endl;
163 return 1;
164 }
165 Use[regMethod] = 1;
166 }
167 }
168
169 // --------------------------------------------------------------------------------------------------
170
171 // Here the preparation phase begins
172
173 // Read training and test data
174 // (it is also possible to use ASCII format as input -> see TMVA Users Guide)
175 TString inputFile = gROOT->GetTutorialDir() + "/tmva/data/tmva_class_example.root";
176 std::unique_ptr<TFile> input{TFile::Open(inputFile)};
177 if (!input || input->IsZombie()) {
178 throw std::runtime_error("ERROR: could not open data file");
179 }
180 std::cout << "--- TMVAClassification : Using input file: " << input->GetName() << std::endl;
181
182 // Register the training and test trees
183
184 TTree *signalTree = (TTree*)input->Get("TreeS");
185 TTree *background = (TTree*)input->Get("TreeB");
186
187 // Create a ROOT output file where TMVA will store ntuples, histograms, etc.
188 TString outfileName("TMVAC.root");
189 std::unique_ptr<TFile> outputFile{TFile::Open(outfileName, "RECREATE")};
190 if (!outputFile || outputFile->IsZombie()) {
191 throw std::runtime_error("ERROR: could not open output file");
192 }
193
194 // Create the factory object. Later you can choose the methods
195 // whose performance you'd like to investigate. The factory is
196 // the only TMVA object you have to interact with
197 //
198 // The first argument is the base of the name of all the
199 // weightfiles in the directory weight/
200 //
201 // The second argument is the output file for the training results
202 // All TMVA output can be suppressed by removing the "!" (not) in
203 // front of the "Silent" argument in the option string
204 auto factory = std::make_unique<TMVA::Factory>(
205 "TMVAClassification", outputFile.get(),
206 "!V:!Silent:Color:DrawProgressBar:Transformations=I;D;P;G,D:AnalysisType=Classification");
207 auto dataloader_raii = std::make_unique<TMVA::DataLoader>("dataset");
208 auto *dataloader = dataloader_raii.get();
209 // If you wish to modify default settings
210 // (please check "src/Config.h" to see all available global options)
211 //
212 // (TMVA::gConfig().GetVariablePlotting()).fTimesRMS = 8.0;
213 // (TMVA::gConfig().GetIONames()).fWeightFileDir = "myWeightDirectory";
214
215 // Define the input variables that shall be used for the MVA training
216 // note that you may also use variable expressions, such as: "3*var1/var2*abs(var3)"
217 // [all types of expressions that can also be parsed by TTree::Draw( "expression" )]
218 dataloader->AddVariable( "myvar1 := var1+var2", 'F' );
219 dataloader->AddVariable( "myvar2 := var1-var2", "Expression 2", "", 'F' );
220 dataloader->AddVariable( "var3", "Variable 3", "units", 'F' );
221 dataloader->AddVariable( "var4", "Variable 4", "units", 'F' );
222
223 // You can add so-called "Spectator variables", which are not used in the MVA training,
224 // but will appear in the final "TestTree" produced by TMVA. This TestTree will contain the
225 // input variables, the response values of all trained MVAs, and the spectator variables
226
227 dataloader->AddSpectator( "spec1 := var1*2", "Spectator 1", "units", 'F' );
228 dataloader->AddSpectator( "spec2 := var1*3", "Spectator 2", "units", 'F' );
229
230
231 // global event weights per tree (see below for setting event-wise weights)
234
235 // You can add an arbitrary number of signal or background trees
236 dataloader->AddSignalTree ( signalTree, signalWeight );
237 dataloader->AddBackgroundTree( background, backgroundWeight );
238
239 // To give different trees for training and testing, do as follows:
240 //
241 // dataloader->AddSignalTree( signalTrainingTree, signalTrainWeight, "Training" );
242 // dataloader->AddSignalTree( signalTestTree, signalTestWeight, "Test" );
243
244 // Use the following code instead of the above two or four lines to add signal and background
245 // training and test events "by hand"
246 // NOTE that in this case one should not give expressions (such as "var1+var2") in the input
247 // variable definition, but simply compute the expression before adding the event
248 // ```cpp
249 // // --- begin ----------------------------------------------------------
250 // std::vector<Double_t> vars( 4 ); // vector has size of number of input variables
251 // Float_t treevars[4], weight;
252 //
253 // // Signal
254 // for (UInt_t ivar=0; ivar<4; ivar++) signalTree->SetBranchAddress( Form( "var%i", ivar+1 ), &(treevars[ivar]) );
255 // for (UInt_t i=0; i<signalTree->GetEntries(); i++) {
256 // signalTree->GetEntry(i);
257 // for (UInt_t ivar=0; ivar<4; ivar++) vars[ivar] = treevars[ivar];
258 // // add training and test events; here: first half is training, second is testing
259 // // note that the weight can also be event-wise
260 // if (i < signalTree->GetEntries()/2.0) dataloader->AddSignalTrainingEvent( vars, signalWeight );
261 // else dataloader->AddSignalTestEvent ( vars, signalWeight );
262 // }
263 //
264 // // Background (has event weights)
265 // background->SetBranchAddress( "weight", &weight );
266 // for (UInt_t ivar=0; ivar<4; ivar++) background->SetBranchAddress( Form( "var%i", ivar+1 ), &(treevars[ivar]) );
267 // for (UInt_t i=0; i<background->GetEntries(); i++) {
268 // background->GetEntry(i);
269 // for (UInt_t ivar=0; ivar<4; ivar++) vars[ivar] = treevars[ivar];
270 // // add training and test events; here: first half is training, second is testing
271 // // note that the weight can also be event-wise
272 // if (i < background->GetEntries()/2) dataloader->AddBackgroundTrainingEvent( vars, backgroundWeight*weight );
273 // else dataloader->AddBackgroundTestEvent ( vars, backgroundWeight*weight );
274 // }
275 // // --- end ------------------------------------------------------------
276 // ```
277 // End of tree registration
278
279 // Set individual event weights (the variables must exist in the original TTree)
280 // - for signal : `dataloader->SetSignalWeightExpression ("weight1*weight2");`
281 // - for background: `dataloader->SetBackgroundWeightExpression("weight1*weight2");`
282 dataloader->SetBackgroundWeightExpression( "weight" );
283
284 // Apply additional cuts on the signal and background samples (can be different)
285 TCut mycuts = ""; // for example: TCut mycuts = "abs(var1)<0.5 && abs(var2-0.5)<1";
286 TCut mycutb = ""; // for example: TCut mycutb = "abs(var1)<0.5";
287
288 // Tell the dataloader how to use the training and testing events
289 //
290 // If no numbers of events are given, half of the events in the tree are used
291 // for training, and the other half for testing:
292 //
293 // dataloader->PrepareTrainingAndTestTree( mycut, "SplitMode=random:!V" );
294 //
295 // To also specify the number of testing events, use:
296 //
297 // dataloader->PrepareTrainingAndTestTree( mycut,
298 // "NSigTrain=3000:NBkgTrain=3000:NSigTest=3000:NBkgTest=3000:SplitMode=Random:!V" );
299 dataloader->PrepareTrainingAndTestTree( mycuts, mycutb,
300 "nTrain_Signal=1000:nTrain_Background=1000:SplitMode=Random:NormMode=NumEvents:!V" );
301
302 // ### Book MVA methods
303 //
304 // Please lookup the various method configuration options in the corresponding cxx files, eg:
305 // src/MethoCuts.cxx, etc, or here: http://tmva.sourceforge.net/old_site/optionRef.html
306 // it is possible to preset ranges in the option string in which the cut optimisation should be done:
307 // "...:CutRangeMin[2]=-1:CutRangeMax[2]=1"...", where [2] is the third input variable
308
309 // Cut optimisation
310 if (Use["Cuts"])
311 factory->BookMethod( dataloader, TMVA::Types::kCuts, "Cuts",
312 "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart" );
313
314 if (Use["CutsD"])
315 factory->BookMethod( dataloader, TMVA::Types::kCuts, "CutsD",
316 "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart:VarTransform=Decorrelate" );
317
318 if (Use["CutsPCA"])
319 factory->BookMethod( dataloader, TMVA::Types::kCuts, "CutsPCA",
320 "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart:VarTransform=PCA" );
321
322 if (Use["CutsGA"])
323 factory->BookMethod( dataloader, TMVA::Types::kCuts, "CutsGA",
324 "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" );
325
326 if (Use["CutsSA"])
327 factory->BookMethod( dataloader, TMVA::Types::kCuts, "CutsSA",
328 "!H:!V:FitMethod=SA:EffSel:MaxCalls=150000:KernelTemp=IncAdaptive:InitialTemp=1e+6:MinTemp=1e-6:Eps=1e-10:UseDefaultScale" );
329
330 // Likelihood ("naive Bayes estimator")
331 if (Use["Likelihood"])
332 factory->BookMethod( dataloader, TMVA::Types::kLikelihood, "Likelihood",
333 "H:!V:TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmoothBkg[1]=10:NSmooth=1:NAvEvtPerBin=50" );
334
335 // Decorrelated likelihood
336 if (Use["LikelihoodD"])
337 factory->BookMethod( dataloader, TMVA::Types::kLikelihood, "LikelihoodD",
338 "!H:!V:TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmooth=5:NAvEvtPerBin=50:VarTransform=Decorrelate" );
339
340 // PCA-transformed likelihood
341 if (Use["LikelihoodPCA"])
342 factory->BookMethod( dataloader, TMVA::Types::kLikelihood, "LikelihoodPCA",
343 "!H:!V:!TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmooth=5:NAvEvtPerBin=50:VarTransform=PCA" );
344
345 // Use a kernel density estimator to approximate the PDFs
346 if (Use["LikelihoodKDE"])
347 factory->BookMethod( dataloader, TMVA::Types::kLikelihood, "LikelihoodKDE",
348 "!H:!V:!TransformOutput:PDFInterpol=KDE:KDEtype=Gauss:KDEiter=Adaptive:KDEFineFactor=0.3:KDEborder=None:NAvEvtPerBin=50" );
349
350 // Use a variable-dependent mix of splines and kernel density estimator
351 if (Use["LikelihoodMIX"])
352 factory->BookMethod( dataloader, TMVA::Types::kLikelihood, "LikelihoodMIX",
353 "!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" );
354
355 // Test the multi-dimensional probability density estimator
356 // here are the options strings for the MinMax and RMS methods, respectively:
357 //
358 // "!H:!V:VolumeRangeMode=MinMax:DeltaFrac=0.2:KernelEstimator=Gauss:GaussSigma=0.3" );
359 // "!H:!V:VolumeRangeMode=RMS:DeltaFrac=3:KernelEstimator=Gauss:GaussSigma=0.3" );
360 if (Use["PDERS"])
361 factory->BookMethod( dataloader, TMVA::Types::kPDERS, "PDERS",
362 "!H:!V:NormTree=T:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600" );
363
364 if (Use["PDERSD"])
365 factory->BookMethod( dataloader, TMVA::Types::kPDERS, "PDERSD",
366 "!H:!V:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600:VarTransform=Decorrelate" );
367
368 if (Use["PDERSPCA"])
369 factory->BookMethod( dataloader, TMVA::Types::kPDERS, "PDERSPCA",
370 "!H:!V:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600:VarTransform=PCA" );
371
372 // Multi-dimensional likelihood estimator using self-adapting phase-space binning
373 if (Use["PDEFoam"])
374 factory->BookMethod( dataloader, TMVA::Types::kPDEFoam, "PDEFoam",
375 "!H:!V:SigBgSeparate=F:TailCut=0.001:VolFrac=0.0666:nActiveCells=500:nSampl=2000:nBin=5:Nmin=100:Kernel=None:Compress=T" );
376
377 if (Use["PDEFoamBoost"])
378 factory->BookMethod( dataloader, TMVA::Types::kPDEFoam, "PDEFoamBoost",
379 "!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" );
380
381 // K-Nearest Neighbour classifier (KNN)
382 if (Use["KNN"])
383 factory->BookMethod( dataloader, TMVA::Types::kKNN, "KNN",
384 "H:nkNN=20:ScaleFrac=0.8:SigmaFact=1.0:Kernel=Gaus:UseKernel=F:UseWeight=T:!Trim" );
385
386 // H-Matrix (chi2-squared) method
387 if (Use["HMatrix"])
388 factory->BookMethod( dataloader, TMVA::Types::kHMatrix, "HMatrix", "!H:!V:VarTransform=None" );
389
390 // Linear discriminant (same as Fisher discriminant)
391 if (Use["LD"])
392 factory->BookMethod( dataloader, TMVA::Types::kLD, "LD", "H:!V:VarTransform=None:CreateMVAPdfs:PDFInterpolMVAPdf=Spline2:NbinsMVAPdf=50:NsmoothMVAPdf=10" );
393
394 // Fisher discriminant (same as LD)
395 if (Use["Fisher"])
396 factory->BookMethod( dataloader, TMVA::Types::kFisher, "Fisher", "H:!V:Fisher:VarTransform=None:CreateMVAPdfs:PDFInterpolMVAPdf=Spline2:NbinsMVAPdf=50:NsmoothMVAPdf=10" );
397
398 // Fisher with Gauss-transformed input variables
399 if (Use["FisherG"])
400 factory->BookMethod( dataloader, TMVA::Types::kFisher, "FisherG", "H:!V:VarTransform=Gauss" );
401
402 // Composite classifier: ensemble (tree) of boosted Fisher classifiers
403 if (Use["BoostedFisher"])
404 factory->BookMethod( dataloader, TMVA::Types::kFisher, "BoostedFisher",
405 "H:!V:Boost_Num=20:Boost_Transform=log:Boost_Type=AdaBoost:Boost_AdaBoostBeta=0.2:!Boost_DetailedMonitoring" );
406
407 // Function discrimination analysis (FDA) -- test of various fitters - the recommended one is Minuit (or GA or SA)
408 if (Use["FDA_MC"])
409 factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_MC",
410 "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" );
411
412 if (Use["FDA_GA"]) // can also use Simulated Annealing (SA) algorithm (see Cuts_SA options])
413 factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_GA",
414 "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=100:Cycles=2:Steps=5:Trim=True:SaveBestGen=1" );
415
416 if (Use["FDA_SA"]) // can also use Simulated Annealing (SA) algorithm (see Cuts_SA options])
417 factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_SA",
418 "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" );
419
420 if (Use["FDA_MT"])
421 factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_MT",
422 "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" );
423
424 if (Use["FDA_GAMT"])
425 factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_GAMT",
426 "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" );
427
428 if (Use["FDA_MCMT"])
429 factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_MCMT",
430 "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" );
431
432 // TMVA ANN: MLP (recommended ANN) -- all ANNs in TMVA are Multilayer Perceptrons
433 if (Use["MLP"])
434 factory->BookMethod( dataloader, TMVA::Types::kMLP, "MLP", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:!UseRegulator" );
435
436 if (Use["MLPBFGS"])
437 factory->BookMethod( dataloader, TMVA::Types::kMLP, "MLPBFGS", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:TrainingMethod=BFGS:!UseRegulator" );
438
439 if (Use["MLPBNN"])
440 factory->BookMethod( dataloader, TMVA::Types::kMLP, "MLPBNN", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=60:HiddenLayers=N+5:TestRate=5:TrainingMethod=BFGS:UseRegulator" ); // BFGS training with bayesian regulators
441
442
443 // Multi-architecture DNN implementation.
444 if (Use["DNN_CPU"] or Use["DNN_GPU"]) {
445 // General layout.
446 TString layoutString ("Layout=TANH|128,TANH|128,TANH|128,LINEAR");
447
448 // Define Training strategy. One could define multiple strategy string separated by the "|" delimiter
449
450 TString trainingStrategyString = ("TrainingStrategy=LearningRate=1e-2,Momentum=0.9,"
451 "ConvergenceSteps=20,BatchSize=100,TestRepetitions=1,"
452 "WeightDecay=1e-4,Regularization=None,"
453 "DropConfig=0.0+0.5+0.5+0.5");
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 // Cuda implementation.
462 if (Use["DNN_GPU"]) {
463 TString gpuOptions = dnnOptions + ":Architecture=GPU";
464 factory->BookMethod(dataloader, TMVA::Types::kDL, "DNN_GPU", gpuOptions);
465 }
466 // Multi-core CPU implementation.
467 if (Use["DNN_CPU"]) {
468 TString cpuOptions = dnnOptions + ":Architecture=CPU";
469 factory->BookMethod(dataloader, TMVA::Types::kDL, "DNN_CPU", cpuOptions);
470 }
471 }
472
473 // CF(Clermont-Ferrand)ANN
474 if (Use["CFMlpANN"])
475 factory->BookMethod( dataloader, TMVA::Types::kCFMlpANN, "CFMlpANN", "!H:!V:NCycles=200:HiddenLayers=N+1,N" ); // n_cycles:#nodes:#nodes:...
476
477 // Tmlp(Root)ANN
478 if (Use["TMlpANN"])
479 factory->BookMethod( dataloader, TMVA::Types::kTMlpANN, "TMlpANN", "!H:!V:NCycles=200:HiddenLayers=N+1,N:LearningMethod=BFGS:ValidationFraction=0.3" ); // n_cycles:#nodes:#nodes:...
480
481 // Support Vector Machine
482 if (Use["SVM"])
483 factory->BookMethod( dataloader, TMVA::Types::kSVM, "SVM", "Gamma=0.25:Tol=0.001:VarTransform=Norm" );
484
485 // Boosted Decision Trees
486 if (Use["BDTG"]) // Gradient Boost
487 factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDTG",
488 "!H:!V:NTrees=1000:MinNodeSize=2.5%:BoostType=Grad:Shrinkage=0.10:UseBaggedBoost:BaggedSampleFraction=0.5:nCuts=20:MaxDepth=2" );
489
490 if (Use["BDT"]) // Adaptive Boost
491 factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDT",
492 "!H:!V:NTrees=850:MinNodeSize=2.5%:MaxDepth=3:BoostType=AdaBoost:AdaBoostBeta=0.5:UseBaggedBoost:BaggedSampleFraction=0.5:SeparationType=GiniIndex:nCuts=20" );
493
494 if (Use["BDTB"]) // Bagging
495 factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDTB",
496 "!H:!V:NTrees=400:BoostType=Bagging:SeparationType=GiniIndex:nCuts=20" );
497
498 if (Use["BDTD"]) // Decorrelation + Adaptive Boost
499 factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDTD",
500 "!H:!V:NTrees=400:MinNodeSize=5%:MaxDepth=3:BoostType=AdaBoost:SeparationType=GiniIndex:nCuts=20:VarTransform=Decorrelate" );
501
502 if (Use["BDTF"]) // Allow Using Fisher discriminant in node splitting for (strong) linearly correlated variables
503 factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDTF",
504 "!H:!V:NTrees=50:MinNodeSize=2.5%:UseFisherCuts:MaxDepth=3:BoostType=AdaBoost:AdaBoostBeta=0.5:SeparationType=GiniIndex:nCuts=20" );
505
506 // RuleFit -- TMVA implementation of Friedman's method
507 if (Use["RuleFit"])
508 factory->BookMethod( dataloader, TMVA::Types::kRuleFit, "RuleFit",
509 "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" );
510
511 // For an example of the category classifier usage, see: TMVAClassificationCategory
512 //
513 // --------------------------------------------------------------------------------------------------
514 // Now you can optimize the setting (configuration) of the MVAs using the set of training events
515 // STILL EXPERIMENTAL and only implemented for BDT's !
516 //
517 // factory->OptimizeAllMethods("SigEffAtBkg0.01","Scan");
518 // factory->OptimizeAllMethods("ROCIntegral","FitGA");
519 //
520 // --------------------------------------------------------------------------------------------------
521
522 // Now you can tell the factory to train, test, and evaluate the MVAs
523 //
524 // Train MVAs using the set of training events
525 factory->TrainAllMethods();
526
527 // Evaluate all MVAs using the set of test events
528 factory->TestAllMethods();
529
530 // Evaluate and compare performance of all configured MVAs
531 factory->EvaluateAllMethods();
532
533 // --------------------------------------------------------------
534
535 // Save the output
536 outputFile->Write();
537
538 std::cout << "==> Wrote root file: " << outputFile->GetName() << std::endl;
539 std::cout << "==> TMVAClassification is done!" << std::endl;
540
541 // Launch the GUI for the root macros
542 if (!gROOT->IsBatch()) TMVA::TMVAGui( outfileName );
543
544 return 0;
545}
546
547int main( int argc, char** argv )
548{
549 // Select methods (don't look at this code - not of interest)
551 for (int i=1; i<argc; i++) {
553 if(regMethod=="-b" || regMethod=="--batch") continue;
554 if (!methodList.IsNull()) methodList += TString(",");
556 }
558}
int main()
Definition Prototype.cxx:12
unsigned int UInt_t
Definition RtypesCore.h:46
double Double_t
Definition RtypesCore.h:59
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void input
#define gROOT
Definition TROOT.h:406
A specialized string object used for TTree selections.
Definition TCut.h:25
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4094
static Tools & Instance()
Definition Tools.cxx:71
std::vector< TString > SplitString(const TString &theOpt, const char separator) const
splits the option string at 'separator' and fills the list 'splitV' with the primitive strings
Definition Tools.cxx:1199
@ kFisher
Definition Types.h:82
@ kTMlpANN
Definition Types.h:85
@ kPDEFoam
Definition Types.h:94
@ kLikelihood
Definition Types.h:79
@ kHMatrix
Definition Types.h:81
@ kRuleFit
Definition Types.h:88
@ kCFMlpANN
Definition Types.h:84
Basic string class.
Definition TString.h:139
A TTree represents a columnar dataset.
Definition TTree.h:79
Tools & gTools()
void TMVAGui(const char *fName="TMVA.root", TString dataset="")