1 |
benhoob |
1.1 |
// @(#)root/tmva $Id: TMVAClassification.C 37399 2010-12-08 15:22:07Z evt $
|
2 |
|
|
/**********************************************************************************
|
3 |
|
|
* Project : TMVA - a ROOT-integrated toolkit for multivariate data analysis *
|
4 |
|
|
* Package : TMVA *
|
5 |
|
|
* Root Macro: TMVAClassification *
|
6 |
|
|
* *
|
7 |
|
|
* This macro provides examples for the training and testing of the *
|
8 |
|
|
* TMVA classifiers. *
|
9 |
|
|
* *
|
10 |
|
|
* As input data is used a toy-MC sample consisting of four Gaussian-distributed *
|
11 |
|
|
* and linearly correlated input variables. *
|
12 |
|
|
* *
|
13 |
|
|
* The methods to be used can be switched on and off by means of booleans, or *
|
14 |
|
|
* via the prompt command, for example: *
|
15 |
|
|
* *
|
16 |
|
|
* root -l ./TMVAClassification.C\(\"Fisher,Likelihood\"\) *
|
17 |
|
|
* *
|
18 |
|
|
* (note that the backslashes are mandatory) *
|
19 |
|
|
* If no method given, a default set of classifiers is used. *
|
20 |
|
|
* *
|
21 |
|
|
* The output file "TMVA.root" can be analysed with the use of dedicated *
|
22 |
|
|
* macros (simply say: root -l <macro.C>), which can be conveniently *
|
23 |
|
|
* invoked through a GUI that will appear at the end of the run of this macro. *
|
24 |
|
|
* Launch the GUI via the command: *
|
25 |
|
|
* *
|
26 |
|
|
* root -l ./TMVAGui.C *
|
27 |
|
|
* *
|
28 |
|
|
**********************************************************************************/
|
29 |
|
|
|
30 |
|
|
#include <cstdlib>
|
31 |
|
|
#include <iostream>
|
32 |
|
|
#include <map>
|
33 |
|
|
#include <string>
|
34 |
|
|
|
35 |
|
|
#include "TChain.h"
|
36 |
|
|
#include "TFile.h"
|
37 |
|
|
#include "TTree.h"
|
38 |
|
|
#include "TString.h"
|
39 |
|
|
#include "TObjString.h"
|
40 |
|
|
#include "TSystem.h"
|
41 |
|
|
#include "TROOT.h"
|
42 |
|
|
|
43 |
|
|
#include "TMVAGui.C"
|
44 |
|
|
|
45 |
|
|
#if not defined(__CINT__) || defined(__MAKECINT__)
|
46 |
|
|
// needs to be included when makecint runs (ACLIC)
|
47 |
|
|
#include "TMVA/Factory.h"
|
48 |
|
|
#include "TMVA/Tools.h"
|
49 |
|
|
#endif
|
50 |
|
|
|
51 |
|
|
void TMVAClassification( TString myMethodList = "" )
|
52 |
|
|
{
|
53 |
|
|
// The explicit loading of the shared libTMVA is done in TMVAlogon.C, defined in .rootrc
|
54 |
|
|
// if you use your private .rootrc, or run from a different directory, please copy the
|
55 |
|
|
// corresponding lines from .rootrc
|
56 |
|
|
|
57 |
|
|
// methods to be processed can be given as an argument; use format:
|
58 |
|
|
//
|
59 |
|
|
// mylinux~> root -l TMVAClassification.C\(\"myMethod1,myMethod2,myMethod3\"\)
|
60 |
|
|
//
|
61 |
|
|
// if you like to use a method via the plugin mechanism, we recommend using
|
62 |
|
|
//
|
63 |
|
|
// mylinux~> root -l TMVAClassification.C\(\"P_myMethod\"\)
|
64 |
|
|
// (an example is given for using the BDT as plugin (see below),
|
65 |
|
|
// but of course the real application is when you write your own
|
66 |
|
|
// method based)
|
67 |
|
|
|
68 |
|
|
//---------------------------------------------------------------
|
69 |
|
|
// This loads the library
|
70 |
|
|
TMVA::Tools::Instance();
|
71 |
|
|
|
72 |
|
|
// Default MVA methods to be trained + tested
|
73 |
|
|
std::map<std::string,int> Use;
|
74 |
|
|
|
75 |
|
|
// --- Cut optimisation
|
76 |
|
|
Use["Cuts"] = 1;
|
77 |
|
|
Use["CutsD"] = 1;
|
78 |
|
|
Use["CutsPCA"] = 0;
|
79 |
|
|
Use["CutsGA"] = 0;
|
80 |
|
|
Use["CutsSA"] = 0;
|
81 |
|
|
//
|
82 |
|
|
// --- 1-dimensional likelihood ("naive Bayes estimator")
|
83 |
|
|
Use["Likelihood"] = 1;
|
84 |
|
|
Use["LikelihoodD"] = 0; // the "D" extension indicates decorrelated input variables (see option strings)
|
85 |
|
|
Use["LikelihoodPCA"] = 1; // the "PCA" extension indicates PCA-transformed input variables (see option strings)
|
86 |
|
|
Use["LikelihoodKDE"] = 0;
|
87 |
|
|
Use["LikelihoodMIX"] = 0;
|
88 |
|
|
//
|
89 |
|
|
// --- Mutidimensional likelihood and Nearest-Neighbour methods
|
90 |
|
|
Use["PDERS"] = 1;
|
91 |
|
|
Use["PDERSD"] = 0;
|
92 |
|
|
Use["PDERSPCA"] = 0;
|
93 |
|
|
Use["PDEFoam"] = 1;
|
94 |
|
|
Use["PDEFoamBoost"] = 0; // uses generalised MVA method boosting
|
95 |
|
|
Use["KNN"] = 1; // k-nearest neighbour method
|
96 |
|
|
//
|
97 |
|
|
// --- Linear Discriminant Analysis
|
98 |
|
|
Use["LD"] = 1; // Linear Discriminant identical to Fisher
|
99 |
|
|
Use["Fisher"] = 0;
|
100 |
|
|
Use["FisherG"] = 0;
|
101 |
|
|
Use["BoostedFisher"] = 0; // uses generalised MVA method boosting
|
102 |
|
|
Use["HMatrix"] = 0;
|
103 |
|
|
//
|
104 |
|
|
// --- Function Discriminant analysis
|
105 |
|
|
Use["FDA_GA"] = 1; // minimisation of user-defined function using Genetics Algorithm
|
106 |
|
|
Use["FDA_SA"] = 0;
|
107 |
|
|
Use["FDA_MC"] = 0;
|
108 |
|
|
Use["FDA_MT"] = 0;
|
109 |
|
|
Use["FDA_GAMT"] = 0;
|
110 |
|
|
Use["FDA_MCMT"] = 0;
|
111 |
|
|
//
|
112 |
|
|
// --- Neural Networks (all are feed-forward Multilayer Perceptrons)
|
113 |
|
|
Use["MLP"] = 0; // Recommended ANN
|
114 |
|
|
Use["MLPBFGS"] = 0; // Recommended ANN with optional training method
|
115 |
|
|
Use["MLPBNN"] = 1; // Recommended ANN with BFGS training method and bayesian regulator
|
116 |
|
|
Use["CFMlpANN"] = 0; // Depreciated ANN from ALEPH
|
117 |
|
|
Use["TMlpANN"] = 0; // ROOT's own ANN
|
118 |
|
|
//
|
119 |
|
|
// --- Support Vector Machine
|
120 |
|
|
Use["SVM"] = 1;
|
121 |
|
|
//
|
122 |
|
|
// --- Boosted Decision Trees
|
123 |
|
|
Use["BDT"] = 1; // uses Adaptive Boost
|
124 |
|
|
Use["BDTG"] = 0; // uses Gradient Boost
|
125 |
|
|
Use["BDTB"] = 0; // uses Bagging
|
126 |
|
|
Use["BDTD"] = 0; // decorrelation + Adaptive Boost
|
127 |
|
|
//
|
128 |
|
|
// --- Friedman's RuleFit method, ie, an optimised series of cuts ("rules")
|
129 |
|
|
Use["RuleFit"] = 1;
|
130 |
|
|
// ---------------------------------------------------------------
|
131 |
|
|
|
132 |
|
|
std::cout << std::endl;
|
133 |
|
|
std::cout << "==> Start TMVAClassification" << std::endl;
|
134 |
|
|
|
135 |
|
|
// Select methods (don't look at this code - not of interest)
|
136 |
|
|
if (myMethodList != "") {
|
137 |
|
|
for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) it->second = 0;
|
138 |
|
|
|
139 |
|
|
std::vector<TString> mlist = TMVA::gTools().SplitString( myMethodList, ',' );
|
140 |
|
|
for (UInt_t i=0; i<mlist.size(); i++) {
|
141 |
|
|
std::string regMethod(mlist[i]);
|
142 |
|
|
|
143 |
|
|
if (Use.find(regMethod) == Use.end()) {
|
144 |
|
|
std::cout << "Method \"" << regMethod << "\" not known in TMVA under this name. Choose among the following:" << std::endl;
|
145 |
|
|
for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) std::cout << it->first << " ";
|
146 |
|
|
std::cout << std::endl;
|
147 |
|
|
return;
|
148 |
|
|
}
|
149 |
|
|
Use[regMethod] = 1;
|
150 |
|
|
}
|
151 |
|
|
}
|
152 |
|
|
|
153 |
|
|
// --------------------------------------------------------------------------------------------------
|
154 |
|
|
|
155 |
|
|
// --- Here the preparation phase begins
|
156 |
|
|
|
157 |
|
|
// Create a ROOT output file where TMVA will store ntuples, histograms, etc.
|
158 |
|
|
TString outfileName( "TMVA.root" );
|
159 |
|
|
TFile* outputFile = TFile::Open( outfileName, "RECREATE" );
|
160 |
|
|
|
161 |
|
|
// Create the factory object. Later you can choose the methods
|
162 |
|
|
// whose performance you'd like to investigate. The factory is
|
163 |
|
|
// the only TMVA object you have to interact with
|
164 |
|
|
//
|
165 |
|
|
// The first argument is the base of the name of all the
|
166 |
|
|
// weightfiles in the directory weight/
|
167 |
|
|
//
|
168 |
|
|
// The second argument is the output file for the training results
|
169 |
|
|
// All TMVA output can be suppressed by removing the "!" (not) in
|
170 |
|
|
// front of the "Silent" argument in the option string
|
171 |
|
|
TMVA::Factory *factory = new TMVA::Factory( "TMVAClassification", outputFile,
|
172 |
|
|
"!V:!Silent:Color:DrawProgressBar:Transformations=I;D;P;G,D:AnalysisType=Classification" );
|
173 |
|
|
|
174 |
|
|
// If you wish to modify default settings
|
175 |
|
|
// (please check "src/Config.h" to see all available global options)
|
176 |
|
|
// (TMVA::gConfig().GetVariablePlotting()).fTimesRMS = 8.0;
|
177 |
|
|
// (TMVA::gConfig().GetIONames()).fWeightFileDir = "myWeightDirectory";
|
178 |
|
|
|
179 |
|
|
// Define the input variables that shall be used for the MVA training
|
180 |
|
|
// note that you may also use variable expressions, such as: "3*var1/var2*abs(var3)"
|
181 |
|
|
// [all types of expressions that can also be parsed by TTree::Draw( "expression" )]
|
182 |
|
|
factory->AddVariable( "myvar1 := var1+var2", 'F' );
|
183 |
|
|
factory->AddVariable( "myvar2 := var1-var2", "Expression 2", "", 'F' );
|
184 |
|
|
factory->AddVariable( "var3", "Variable 3", "units", 'F' );
|
185 |
|
|
factory->AddVariable( "var4", "Variable 4", "units", 'F' );
|
186 |
|
|
|
187 |
|
|
// You can add so-called "Spectator variables", which are not used in the MVA training,
|
188 |
|
|
// but will appear in the final "TestTree" produced by TMVA. This TestTree will contain the
|
189 |
|
|
// input variables, the response values of all trained MVAs, and the spectator variables
|
190 |
|
|
factory->AddSpectator( "spec1 := var1*2", "Spectator 1", "units", 'F' );
|
191 |
|
|
factory->AddSpectator( "spec2 := var1*3", "Spectator 2", "units", 'F' );
|
192 |
|
|
|
193 |
|
|
// Read training and test data
|
194 |
|
|
// (it is also possible to use ASCII format as input -> see TMVA Users Guide)
|
195 |
|
|
TString fname = "./tmva_class_example.root";
|
196 |
|
|
|
197 |
|
|
if (gSystem->AccessPathName( fname )) // file does not exist in local directory
|
198 |
|
|
gSystem->Exec("wget http://root.cern.ch/files/tmva_class_example.root");
|
199 |
|
|
|
200 |
|
|
TFile *input = TFile::Open( fname );
|
201 |
|
|
|
202 |
|
|
std::cout << "--- TMVAClassification : Using input file: " << input->GetName() << std::endl;
|
203 |
|
|
|
204 |
|
|
// --- Register the training and test trees
|
205 |
|
|
|
206 |
|
|
TTree *signal = (TTree*)input->Get("TreeS");
|
207 |
|
|
TTree *background = (TTree*)input->Get("TreeB");
|
208 |
|
|
|
209 |
|
|
// global event weights per tree (see below for setting event-wise weights)
|
210 |
|
|
Double_t signalWeight = 1.0;
|
211 |
|
|
Double_t backgroundWeight = 1.0;
|
212 |
|
|
|
213 |
|
|
// You can add an arbitrary number of signal or background trees
|
214 |
|
|
factory->AddSignalTree ( signal, signalWeight );
|
215 |
|
|
factory->AddBackgroundTree( background, backgroundWeight );
|
216 |
|
|
|
217 |
|
|
// To give different trees for training and testing, do as follows:
|
218 |
|
|
// factory->AddSignalTree( signalTrainingTree, signalTrainWeight, "Training" );
|
219 |
|
|
// factory->AddSignalTree( signalTestTree, signalTestWeight, "Test" );
|
220 |
|
|
|
221 |
|
|
// Use the following code instead of the above two or four lines to add signal and background
|
222 |
|
|
// training and test events "by hand"
|
223 |
|
|
// NOTE that in this case one should not give expressions (such as "var1+var2") in the input
|
224 |
|
|
// variable definition, but simply compute the expression before adding the event
|
225 |
|
|
//
|
226 |
|
|
// // --- begin ----------------------------------------------------------
|
227 |
|
|
// std::vector<Double_t> vars( 4 ); // vector has size of number of input variables
|
228 |
|
|
// Float_t treevars[4], weight;
|
229 |
|
|
//
|
230 |
|
|
// // Signal
|
231 |
|
|
// for (UInt_t ivar=0; ivar<4; ivar++) signal->SetBranchAddress( Form( "var%i", ivar+1 ), &(treevars[ivar]) );
|
232 |
|
|
// for (UInt_t i=0; i<signal->GetEntries(); i++) {
|
233 |
|
|
// signal->GetEntry(i);
|
234 |
|
|
// for (UInt_t ivar=0; ivar<4; ivar++) vars[ivar] = treevars[ivar];
|
235 |
|
|
// // add training and test events; here: first half is training, second is testing
|
236 |
|
|
// // note that the weight can also be event-wise
|
237 |
|
|
// if (i < signal->GetEntries()/2.0) factory->AddSignalTrainingEvent( vars, signalWeight );
|
238 |
|
|
// else factory->AddSignalTestEvent ( vars, signalWeight );
|
239 |
|
|
// }
|
240 |
|
|
//
|
241 |
|
|
// // Background (has event weights)
|
242 |
|
|
// background->SetBranchAddress( "weight", &weight );
|
243 |
|
|
// for (UInt_t ivar=0; ivar<4; ivar++) background->SetBranchAddress( Form( "var%i", ivar+1 ), &(treevars[ivar]) );
|
244 |
|
|
// for (UInt_t i=0; i<background->GetEntries(); i++) {
|
245 |
|
|
// background->GetEntry(i);
|
246 |
|
|
// for (UInt_t ivar=0; ivar<4; ivar++) vars[ivar] = treevars[ivar];
|
247 |
|
|
// // add training and test events; here: first half is training, second is testing
|
248 |
|
|
// // note that the weight can also be event-wise
|
249 |
|
|
// if (i < background->GetEntries()/2) factory->AddBackgroundTrainingEvent( vars, backgroundWeight*weight );
|
250 |
|
|
// else factory->AddBackgroundTestEvent ( vars, backgroundWeight*weight );
|
251 |
|
|
// }
|
252 |
|
|
// --- end ------------------------------------------------------------
|
253 |
|
|
//
|
254 |
|
|
// --- end of tree registration
|
255 |
|
|
|
256 |
|
|
// Set individual event weights (the variables must exist in the original TTree)
|
257 |
|
|
// for signal : factory->SetSignalWeightExpression ("weight1*weight2");
|
258 |
|
|
// for background: factory->SetBackgroundWeightExpression("weight1*weight2");
|
259 |
|
|
factory->SetBackgroundWeightExpression( "weight" );
|
260 |
|
|
|
261 |
|
|
// Apply additional cuts on the signal and background samples (can be different)
|
262 |
|
|
TCut mycuts = ""; // for example: TCut mycuts = "abs(var1)<0.5 && abs(var2-0.5)<1";
|
263 |
|
|
TCut mycutb = ""; // for example: TCut mycutb = "abs(var1)<0.5";
|
264 |
|
|
|
265 |
|
|
// Tell the factory how to use the training and testing events
|
266 |
|
|
//
|
267 |
|
|
// If no numbers of events are given, half of the events in the tree are used
|
268 |
|
|
// for training, and the other half for testing:
|
269 |
|
|
// factory->PrepareTrainingAndTestTree( mycut, "SplitMode=random:!V" );
|
270 |
|
|
// To also specify the number of testing events, use:
|
271 |
|
|
// factory->PrepareTrainingAndTestTree( mycut,
|
272 |
|
|
// "NSigTrain=3000:NBkgTrain=3000:NSigTest=3000:NBkgTest=3000:SplitMode=Random:!V" );
|
273 |
|
|
factory->PrepareTrainingAndTestTree( mycuts, mycutb,
|
274 |
|
|
"nTrain_Signal=0:nTrain_Background=0:SplitMode=Random:NormMode=NumEvents:!V" );
|
275 |
|
|
|
276 |
|
|
// ---- Book MVA methods
|
277 |
|
|
//
|
278 |
|
|
// Please lookup the various method configuration options in the corresponding cxx files, eg:
|
279 |
|
|
// src/MethoCuts.cxx, etc, or here: http://tmva.sourceforge.net/optionRef.html
|
280 |
|
|
// it is possible to preset ranges in the option string in which the cut optimisation should be done:
|
281 |
|
|
// "...:CutRangeMin[2]=-1:CutRangeMax[2]=1"...", where [2] is the third input variable
|
282 |
|
|
|
283 |
|
|
// Cut optimisation
|
284 |
|
|
if (Use["Cuts"])
|
285 |
|
|
factory->BookMethod( TMVA::Types::kCuts, "Cuts",
|
286 |
|
|
"!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart" );
|
287 |
|
|
|
288 |
|
|
if (Use["CutsD"])
|
289 |
|
|
factory->BookMethod( TMVA::Types::kCuts, "CutsD",
|
290 |
|
|
"!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart:VarTransform=Decorrelate" );
|
291 |
|
|
|
292 |
|
|
if (Use["CutsPCA"])
|
293 |
|
|
factory->BookMethod( TMVA::Types::kCuts, "CutsPCA",
|
294 |
|
|
"!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart:VarTransform=PCA" );
|
295 |
|
|
|
296 |
|
|
if (Use["CutsGA"])
|
297 |
|
|
factory->BookMethod( TMVA::Types::kCuts, "CutsGA",
|
298 |
|
|
"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" );
|
299 |
|
|
|
300 |
|
|
if (Use["CutsSA"])
|
301 |
|
|
factory->BookMethod( TMVA::Types::kCuts, "CutsSA",
|
302 |
|
|
"!H:!V:FitMethod=SA:EffSel:MaxCalls=150000:KernelTemp=IncAdaptive:InitialTemp=1e+6:MinTemp=1e-6:Eps=1e-10:UseDefaultScale" );
|
303 |
|
|
|
304 |
|
|
// Likelihood ("naive Bayes estimator")
|
305 |
|
|
if (Use["Likelihood"])
|
306 |
|
|
factory->BookMethod( TMVA::Types::kLikelihood, "Likelihood",
|
307 |
|
|
"H:!V:!TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmoothBkg[1]=10:NSmooth=1:NAvEvtPerBin=50" );
|
308 |
|
|
|
309 |
|
|
// Decorrelated likelihood
|
310 |
|
|
if (Use["LikelihoodD"])
|
311 |
|
|
factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodD",
|
312 |
|
|
"!H:!V:TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmooth=5:NAvEvtPerBin=50:VarTransform=Decorrelate" );
|
313 |
|
|
|
314 |
|
|
// PCA-transformed likelihood
|
315 |
|
|
if (Use["LikelihoodPCA"])
|
316 |
|
|
factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodPCA",
|
317 |
|
|
"!H:!V:!TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmooth=5:NAvEvtPerBin=50:VarTransform=PCA" );
|
318 |
|
|
|
319 |
|
|
// Use a kernel density estimator to approximate the PDFs
|
320 |
|
|
if (Use["LikelihoodKDE"])
|
321 |
|
|
factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodKDE",
|
322 |
|
|
"!H:!V:!TransformOutput:PDFInterpol=KDE:KDEtype=Gauss:KDEiter=Adaptive:KDEFineFactor=0.3:KDEborder=None:NAvEvtPerBin=50" );
|
323 |
|
|
|
324 |
|
|
// Use a variable-dependent mix of splines and kernel density estimator
|
325 |
|
|
if (Use["LikelihoodMIX"])
|
326 |
|
|
factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodMIX",
|
327 |
|
|
"!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" );
|
328 |
|
|
|
329 |
|
|
// Test the multi-dimensional probability density estimator
|
330 |
|
|
// here are the options strings for the MinMax and RMS methods, respectively:
|
331 |
|
|
// "!H:!V:VolumeRangeMode=MinMax:DeltaFrac=0.2:KernelEstimator=Gauss:GaussSigma=0.3" );
|
332 |
|
|
// "!H:!V:VolumeRangeMode=RMS:DeltaFrac=3:KernelEstimator=Gauss:GaussSigma=0.3" );
|
333 |
|
|
if (Use["PDERS"])
|
334 |
|
|
factory->BookMethod( TMVA::Types::kPDERS, "PDERS",
|
335 |
|
|
"!H:!V:NormTree=T:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600" );
|
336 |
|
|
|
337 |
|
|
if (Use["PDERSD"])
|
338 |
|
|
factory->BookMethod( TMVA::Types::kPDERS, "PDERSD",
|
339 |
|
|
"!H:!V:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600:VarTransform=Decorrelate" );
|
340 |
|
|
|
341 |
|
|
if (Use["PDERSPCA"])
|
342 |
|
|
factory->BookMethod( TMVA::Types::kPDERS, "PDERSPCA",
|
343 |
|
|
"!H:!V:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600:VarTransform=PCA" );
|
344 |
|
|
|
345 |
|
|
// Multi-dimensional likelihood estimator using self-adapting phase-space binning
|
346 |
|
|
if (Use["PDEFoam"])
|
347 |
|
|
factory->BookMethod( TMVA::Types::kPDEFoam, "PDEFoam",
|
348 |
|
|
"H:!V:SigBgSeparate=F:TailCut=0.001:VolFrac=0.0333:nActiveCells=500:nSampl=2000:nBin=5:Nmin=100:Kernel=None:Compress=T" );
|
349 |
|
|
|
350 |
|
|
if (Use["PDEFoamBoost"])
|
351 |
|
|
factory->BookMethod( TMVA::Types::kPDEFoam, "PDEFoamBoost",
|
352 |
|
|
"!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" );
|
353 |
|
|
|
354 |
|
|
// K-Nearest Neighbour classifier (KNN)
|
355 |
|
|
if (Use["KNN"])
|
356 |
|
|
factory->BookMethod( TMVA::Types::kKNN, "KNN",
|
357 |
|
|
"H:nkNN=20:ScaleFrac=0.8:SigmaFact=1.0:Kernel=Gaus:UseKernel=F:UseWeight=T:!Trim" );
|
358 |
|
|
|
359 |
|
|
// H-Matrix (chi2-squared) method
|
360 |
|
|
if (Use["HMatrix"])
|
361 |
|
|
factory->BookMethod( TMVA::Types::kHMatrix, "HMatrix", "!H:!V" );
|
362 |
|
|
|
363 |
|
|
// Linear discriminant (same as Fisher discriminant)
|
364 |
|
|
if (Use["LD"])
|
365 |
|
|
factory->BookMethod( TMVA::Types::kLD, "LD", "H:!V:VarTransform=None:CreateMVAPdfs:PDFInterpolMVAPdf=Spline2:NbinsMVAPdf=50:NsmoothMVAPdf=10" );
|
366 |
|
|
|
367 |
|
|
// Fisher discriminant (same as LD)
|
368 |
|
|
if (Use["Fisher"])
|
369 |
|
|
factory->BookMethod( TMVA::Types::kFisher, "Fisher", "H:!V:Fisher:CreateMVAPdfs:PDFInterpolMVAPdf=Spline2:NbinsMVAPdf=50:NsmoothMVAPdf=10" );
|
370 |
|
|
|
371 |
|
|
// Fisher with Gauss-transformed input variables
|
372 |
|
|
if (Use["FisherG"])
|
373 |
|
|
factory->BookMethod( TMVA::Types::kFisher, "FisherG", "H:!V:VarTransform=Gauss" );
|
374 |
|
|
|
375 |
|
|
// Composite classifier: ensemble (tree) of boosted Fisher classifiers
|
376 |
|
|
if (Use["BoostedFisher"])
|
377 |
|
|
factory->BookMethod( TMVA::Types::kFisher, "BoostedFisher",
|
378 |
|
|
"H:!V:Boost_Num=20:Boost_Transform=log:Boost_Type=AdaBoost:Boost_AdaBoostBeta=0.2" );
|
379 |
|
|
|
380 |
|
|
// Function discrimination analysis (FDA) -- test of various fitters - the recommended one is Minuit (or GA or SA)
|
381 |
|
|
if (Use["FDA_MC"])
|
382 |
|
|
factory->BookMethod( TMVA::Types::kFDA, "FDA_MC",
|
383 |
|
|
"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" );
|
384 |
|
|
|
385 |
|
|
if (Use["FDA_GA"]) // can also use Simulated Annealing (SA) algorithm (see Cuts_SA options])
|
386 |
|
|
factory->BookMethod( TMVA::Types::kFDA, "FDA_GA",
|
387 |
|
|
"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" );
|
388 |
|
|
|
389 |
|
|
if (Use["FDA_SA"]) // can also use Simulated Annealing (SA) algorithm (see Cuts_SA options])
|
390 |
|
|
factory->BookMethod( TMVA::Types::kFDA, "FDA_SA",
|
391 |
|
|
"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" );
|
392 |
|
|
|
393 |
|
|
if (Use["FDA_MT"])
|
394 |
|
|
factory->BookMethod( TMVA::Types::kFDA, "FDA_MT",
|
395 |
|
|
"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" );
|
396 |
|
|
|
397 |
|
|
if (Use["FDA_GAMT"])
|
398 |
|
|
factory->BookMethod( TMVA::Types::kFDA, "FDA_GAMT",
|
399 |
|
|
"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" );
|
400 |
|
|
|
401 |
|
|
if (Use["FDA_MCMT"])
|
402 |
|
|
factory->BookMethod( TMVA::Types::kFDA, "FDA_MCMT",
|
403 |
|
|
"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" );
|
404 |
|
|
|
405 |
|
|
// TMVA ANN: MLP (recommended ANN) -- all ANNs in TMVA are Multilayer Perceptrons
|
406 |
|
|
if (Use["MLP"])
|
407 |
|
|
factory->BookMethod( TMVA::Types::kMLP, "MLP", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:!UseRegulator" );
|
408 |
|
|
|
409 |
|
|
if (Use["MLPBFGS"])
|
410 |
|
|
factory->BookMethod( TMVA::Types::kMLP, "MLPBFGS", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:TrainingMethod=BFGS:!UseRegulator" );
|
411 |
|
|
|
412 |
|
|
if (Use["MLPBNN"])
|
413 |
|
|
factory->BookMethod( 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
|
414 |
|
|
|
415 |
|
|
// CF(Clermont-Ferrand)ANN
|
416 |
|
|
if (Use["CFMlpANN"])
|
417 |
|
|
factory->BookMethod( TMVA::Types::kCFMlpANN, "CFMlpANN", "!H:!V:NCycles=2000:HiddenLayers=N+1,N" ); // n_cycles:#nodes:#nodes:...
|
418 |
|
|
|
419 |
|
|
// Tmlp(Root)ANN
|
420 |
|
|
if (Use["TMlpANN"])
|
421 |
|
|
factory->BookMethod( TMVA::Types::kTMlpANN, "TMlpANN", "!H:!V:NCycles=200:HiddenLayers=N+1,N:LearningMethod=BFGS:ValidationFraction=0.3" ); // n_cycles:#nodes:#nodes:...
|
422 |
|
|
|
423 |
|
|
// Support Vector Machine
|
424 |
|
|
if (Use["SVM"])
|
425 |
|
|
factory->BookMethod( TMVA::Types::kSVM, "SVM", "Gamma=0.25:Tol=0.001:VarTransform=Norm" );
|
426 |
|
|
|
427 |
|
|
// Boosted Decision Trees
|
428 |
|
|
if (Use["BDTG"]) // Gradient Boost
|
429 |
|
|
factory->BookMethod( TMVA::Types::kBDT, "BDTG",
|
430 |
|
|
"!H:!V:NTrees=1000:BoostType=Grad:Shrinkage=0.10:UseBaggedGrad:GradBaggingFraction=0.5:nCuts=20:NNodesMax=5" );
|
431 |
|
|
|
432 |
|
|
if (Use["BDT"]) // Adaptive Boost
|
433 |
|
|
factory->BookMethod( TMVA::Types::kBDT, "BDT",
|
434 |
|
|
"!H:!V:NTrees=850:nEventsMin=150:MaxDepth=3:BoostType=AdaBoost:AdaBoostBeta=0.5:SeparationType=GiniIndex:nCuts=20:PruneMethod=NoPruning" );
|
435 |
|
|
|
436 |
|
|
if (Use["BDTB"]) // Bagging
|
437 |
|
|
factory->BookMethod( TMVA::Types::kBDT, "BDTB",
|
438 |
|
|
"!H:!V:NTrees=400:BoostType=Bagging:SeparationType=GiniIndex:nCuts=20:PruneMethod=NoPruning" );
|
439 |
|
|
|
440 |
|
|
if (Use["BDTD"]) // Decorrelation + Adaptive Boost
|
441 |
|
|
factory->BookMethod( TMVA::Types::kBDT, "BDTD",
|
442 |
|
|
"!H:!V:NTrees=400:nEventsMin=400:MaxDepth=3:BoostType=AdaBoost:SeparationType=GiniIndex:nCuts=20:PruneMethod=NoPruning:VarTransform=Decorrelate" );
|
443 |
|
|
|
444 |
|
|
// RuleFit -- TMVA implementation of Friedman's method
|
445 |
|
|
if (Use["RuleFit"])
|
446 |
|
|
factory->BookMethod( TMVA::Types::kRuleFit, "RuleFit",
|
447 |
|
|
"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" );
|
448 |
|
|
|
449 |
|
|
// For an example of the category classifier usage, see: TMVAClassificationCategory
|
450 |
|
|
|
451 |
|
|
// --------------------------------------------------------------------------------------------------
|
452 |
|
|
|
453 |
|
|
// ---- Now you can optimize the setting (configuration) of the MVAs using the set of training events
|
454 |
|
|
|
455 |
|
|
// factory->OptimizeAllMethods("SigEffAt001","Scan");
|
456 |
|
|
// factory->OptimizeAllMethods("ROCIntegral","GA");
|
457 |
|
|
|
458 |
|
|
// --------------------------------------------------------------------------------------------------
|
459 |
|
|
|
460 |
|
|
// ---- Now you can tell the factory to train, test, and evaluate the MVAs
|
461 |
|
|
|
462 |
|
|
// Train MVAs using the set of training events
|
463 |
|
|
factory->TrainAllMethods();
|
464 |
|
|
|
465 |
|
|
// ---- Evaluate all MVAs using the set of test events
|
466 |
|
|
factory->TestAllMethods();
|
467 |
|
|
|
468 |
|
|
// ----- Evaluate and compare performance of all configured MVAs
|
469 |
|
|
factory->EvaluateAllMethods();
|
470 |
|
|
|
471 |
|
|
// --------------------------------------------------------------
|
472 |
|
|
|
473 |
|
|
// Save the output
|
474 |
|
|
outputFile->Close();
|
475 |
|
|
|
476 |
|
|
std::cout << "==> Wrote root file: " << outputFile->GetName() << std::endl;
|
477 |
|
|
std::cout << "==> TMVAClassification is done!" << std::endl;
|
478 |
|
|
|
479 |
|
|
delete factory;
|
480 |
|
|
|
481 |
|
|
// Launch the GUI for the root macros
|
482 |
|
|
if (!gROOT->IsBatch()) TMVAGui( outfileName );
|
483 |
|
|
}
|