Logo ROOT   6.08/07
Reference Guide
TMVARegression.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 ///
10 /// The methods to be used can be switched on and off by means of booleans, or
11 /// via the prompt command, for example:
12 ///
13 /// root -l TMVARegression.C\(\"LD,MLP\"\)
14 ///
15 /// (note that the backslashes are mandatory)
16 /// If no method given, a default set is used.
17 ///
18 /// The output file "TMVAReg.root" can be analysed with the use of dedicated
19 /// macros (simply say: root -l <macro.C>), which can be conveniently
20 /// invoked through a GUI that will appear at the end of the run of this macro.
21 /// - Project : TMVA - a Root-integrated toolkit for multivariate data analysis
22 /// - Package : TMVA
23 /// - Root Macro: TMVARegression
24 ///
25 /// \macro_output
26 /// \macro_code
27 /// \author Andreas Hoecker
28 
29 #include <cstdlib>
30 #include <iostream>
31 #include <map>
32 #include <string>
33 
34 #include "TChain.h"
35 #include "TFile.h"
36 #include "TTree.h"
37 #include "TString.h"
38 #include "TObjString.h"
39 #include "TSystem.h"
40 #include "TROOT.h"
41 
42 #include "TMVA/Tools.h"
43 #include "TMVA/Factory.h"
44 #include "TMVA/DataLoader.h"
45 #include "TMVA/TMVARegGui.h"
46 
47 
48 using namespace TMVA;
49 
50 void TMVARegression( TString myMethodList = "" )
51 {
52  // The explicit loading of the shared libTMVA is done in TMVAlogon.C, defined in .rootrc
53  // if you use your private .rootrc, or run from a different directory, please copy the
54  // corresponding lines from .rootrc
55 
56  // methods to be processed can be given as an argument; use format:
57  //
58  // mylinux~> root -l TMVARegression.C\(\"myMethod1,myMethod2,myMethod3\"\)
59  //
60 
61  //---------------------------------------------------------------
62  // This loads the library
64 
65 
66 
67  // Default MVA methods to be trained + tested
68  std::map<std::string,int> Use;
69 
70  // Mutidimensional likelihood and Nearest-Neighbour methods
71  Use["PDERS"] = 0;
72  Use["PDEFoam"] = 1;
73  Use["KNN"] = 1;
74  //
75  // Linear Discriminant Analysis
76  Use["LD"] = 1;
77  //
78  // Function Discriminant analysis
79  Use["FDA_GA"] = 1;
80  Use["FDA_MC"] = 0;
81  Use["FDA_MT"] = 0;
82  Use["FDA_GAMT"] = 0;
83  //
84  // Neural Network
85  Use["MLP"] = 1;
86  Use["DNN"] = 0;
87  //
88  // Support Vector Machine
89  Use["SVM"] = 0;
90  //
91  // Boosted Decision Trees
92  Use["BDT"] = 0;
93  Use["BDTG"] = 1;
94  // ---------------------------------------------------------------
95 
96  std::cout << std::endl;
97  std::cout << "==> Start TMVARegression" << std::endl;
98 
99  // Select methods (don't look at this code - not of interest)
100  if (myMethodList != "") {
101  for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) it->second = 0;
102 
103  std::vector<TString> mlist = gTools().SplitString( myMethodList, ',' );
104  for (UInt_t i=0; i<mlist.size(); i++) {
105  std::string regMethod(mlist[i]);
106 
107  if (Use.find(regMethod) == Use.end()) {
108  std::cout << "Method \"" << regMethod << "\" not known in TMVA under this name. Choose among the following:" << std::endl;
109  for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) std::cout << it->first << " ";
110  std::cout << std::endl;
111  return;
112  }
113  Use[regMethod] = 1;
114  }
115  }
116 
117  // --------------------------------------------------------------------------------------------------
118 
119  // Here the preparation phase begins
120 
121  // Create a new root output file
122  TString outfileName( "TMVAReg.root" );
123  TFile* outputFile = TFile::Open( outfileName, "RECREATE" );
124 
125  // Create the factory object. Later you can choose the methods
126  // whose performance you'd like to investigate. The factory will
127  // then run the performance analysis for you.
128  //
129  // The first argument is the base of the name of all the
130  // weightfiles in the directory weight/
131  //
132  // The second argument is the output file for the training results
133  // All TMVA output can be suppressed by removing the "!" (not) in
134  // front of the "Silent" argument in the option string
135  TMVA::Factory *factory = new TMVA::Factory( "TMVARegression", outputFile,
136  "!V:!Silent:Color:DrawProgressBar:AnalysisType=Regression" );
137 
138 
140  // If you wish to modify default settings
141  // (please check "src/Config.h" to see all available global options)
142  //
143  // (TMVA::gConfig().GetVariablePlotting()).fTimesRMS = 8.0;
144  // (TMVA::gConfig().GetIONames()).fWeightFileDir = "myWeightDirectory";
145 
146  // Define the input variables that shall be used for the MVA training
147  // note that you may also use variable expressions, such as: "3*var1/var2*abs(var3)"
148  // [all types of expressions that can also be parsed by TTree::Draw( "expression" )]
149  dataloader->AddVariable( "var1", "Variable 1", "units", 'F' );
150  dataloader->AddVariable( "var2", "Variable 2", "units", 'F' );
151 
152  // You can add so-called "Spectator variables", which are not used in the MVA training,
153  // but will appear in the final "TestTree" produced by TMVA. This TestTree will contain the
154  // input variables, the response values of all trained MVAs, and the spectator variables
155  dataloader->AddSpectator( "spec1:=var1*2", "Spectator 1", "units", 'F' );
156  dataloader->AddSpectator( "spec2:=var1*3", "Spectator 2", "units", 'F' );
157 
158  // Add the variable carrying the regression target
159  dataloader->AddTarget( "fvalue" );
160 
161  // It is also possible to declare additional targets for multi-dimensional regression, ie:
162  // factory->AddTarget( "fvalue2" );
163  // BUT: this is currently ONLY implemented for MLP
164 
165  // Read training and test data (see TMVAClassification for reading ASCII files)
166  // load the signal and background event samples from ROOT trees
167  TFile *input(0);
168  TString fname = "./tmva_reg_example.root";
169  if (!gSystem->AccessPathName( fname ))
170  input = TFile::Open( fname ); // check if file in local directory exists
171  else
172  input = TFile::Open( "http://root.cern.ch/files/tmva_reg_example.root" ); // if not: download from ROOT server
173 
174  if (!input) {
175  std::cout << "ERROR: could not open data file" << std::endl;
176  exit(1);
177  }
178  std::cout << "--- TMVARegression : Using input file: " << input->GetName() << std::endl;
179 
180  // Register the regression tree
181 
182  TTree *regTree = (TTree*)input->Get("TreeR");
183 
184  // global event weights per tree (see below for setting event-wise weights)
185  Double_t regWeight = 1.0;
186 
187  // You can add an arbitrary number of regression trees
188  dataloader->AddRegressionTree( regTree, regWeight );
189 
190  // This would set individual event weights (the variables defined in the
191  // expression need to exist in the original TTree)
192  dataloader->SetWeightExpression( "var1", "Regression" );
193 
194  // Apply additional cuts on the signal and background samples (can be different)
195  TCut mycut = ""; // for example: TCut mycut = "abs(var1)<0.5 && abs(var2-0.5)<1";
196 
197  // tell the DataLoader to use all remaining events in the trees after training for testing:
198  dataloader->PrepareTrainingAndTestTree( mycut,
199  "nTrain_Regression=1000:nTest_Regression=0:SplitMode=Random:NormMode=NumEvents:!V" );
200  //
201  // dataloader->PrepareTrainingAndTestTree( mycut,
202  // "nTrain_Regression=0:nTest_Regression=0:SplitMode=Random:NormMode=NumEvents:!V" );
203 
204  // If no numbers of events are given, half of the events in the tree are used
205  // for training, and the other half for testing:
206  //
207  // dataloader->PrepareTrainingAndTestTree( mycut, "SplitMode=random:!V" );
208 
209  // Book MVA methods
210  //
211  // Please lookup the various method configuration options in the corresponding cxx files, eg:
212  // src/MethoCuts.cxx, etc, or here: http://tmva.sourceforge.net/optionRef.html
213  // it is possible to preset ranges in the option string in which the cut optimisation should be done:
214  // "...:CutRangeMin[2]=-1:CutRangeMax[2]=1"...", where [2] is the third input variable
215 
216  // PDE - RS method
217  if (Use["PDERS"])
218  factory->BookMethod( dataloader, TMVA::Types::kPDERS, "PDERS",
219  "!H:!V:NormTree=T:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=40:NEventsMax=60:VarTransform=None" );
220  // And the options strings for the MinMax and RMS methods, respectively:
221  //
222  // "!H:!V:VolumeRangeMode=MinMax:DeltaFrac=0.2:KernelEstimator=Gauss:GaussSigma=0.3" );
223  // "!H:!V:VolumeRangeMode=RMS:DeltaFrac=3:KernelEstimator=Gauss:GaussSigma=0.3" );
224 
225  if (Use["PDEFoam"])
226  factory->BookMethod( dataloader, TMVA::Types::kPDEFoam, "PDEFoam",
227  "!H:!V:MultiTargetRegression=F:TargetSelection=Mpv:TailCut=0.001:VolFrac=0.0666:nActiveCells=500:nSampl=2000:nBin=5:Compress=T:Kernel=None:Nmin=10:VarTransform=None" );
228 
229  // K-Nearest Neighbour classifier (KNN)
230  if (Use["KNN"])
231  factory->BookMethod( dataloader, TMVA::Types::kKNN, "KNN",
232  "nkNN=20:ScaleFrac=0.8:SigmaFact=1.0:Kernel=Gaus:UseKernel=F:UseWeight=T:!Trim" );
233 
234  // Linear discriminant
235  if (Use["LD"])
236  factory->BookMethod( dataloader, TMVA::Types::kLD, "LD",
237  "!H:!V:VarTransform=None" );
238 
239  // Function discrimination analysis (FDA) -- test of various fitters - the recommended one is Minuit (or GA or SA)
240  if (Use["FDA_MC"])
241  factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_MC",
242  "!H:!V:Formula=(0)+(1)*x0+(2)*x1:ParRanges=(-100,100);(-100,100);(-100,100):FitMethod=MC:SampleSize=100000:Sigma=0.1:VarTransform=D" );
243 
244  if (Use["FDA_GA"]) // can also use Simulated Annealing (SA) algorithm (see Cuts_SA options) .. the formula of this example is good for parabolas
245  factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_GA",
246  "!H:!V:Formula=(0)+(1)*x0+(2)*x1:ParRanges=(-100,100);(-100,100);(-100,100):FitMethod=GA:PopSize=100:Cycles=3:Steps=30:Trim=True:SaveBestGen=1:VarTransform=Norm" );
247 
248  if (Use["FDA_MT"])
249  factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_MT",
250  "!H:!V:Formula=(0)+(1)*x0+(2)*x1:ParRanges=(-100,100);(-100,100);(-100,100);(-10,10):FitMethod=MINUIT:ErrorLevel=1:PrintLevel=-1:FitStrategy=2:UseImprove:UseMinos:SetBatch" );
251 
252  if (Use["FDA_GAMT"])
253  factory->BookMethod( dataloader, TMVA::Types::kFDA, "FDA_GAMT",
254  "!H:!V:Formula=(0)+(1)*x0+(2)*x1:ParRanges=(-100,100);(-100,100);(-100,100):FitMethod=GA:Converger=MINUIT:ErrorLevel=1:PrintLevel=-1:FitStrategy=0:!UseImprove:!UseMinos:SetBatch:Cycles=1:PopSize=5:Steps=5:Trim" );
255 
256  // Neural network (MLP)
257  if (Use["MLP"])
258  factory->BookMethod( dataloader, TMVA::Types::kMLP, "MLP", "!H:!V:VarTransform=Norm:NeuronType=tanh:NCycles=20000:HiddenLayers=N+20:TestRate=6:TrainingMethod=BFGS:Sampling=0.3:SamplingEpoch=0.8:ConvergenceImprove=1e-6:ConvergenceTests=15:!UseRegulator" );
259 
260  if (Use["DNN"])
261  {
262  /*
263  TString layoutString ("Layout=TANH|(N+100)*2,LINEAR");
264  TString layoutString ("Layout=SOFTSIGN|100,SOFTSIGN|50,SOFTSIGN|20,LINEAR");
265  TString layoutString ("Layout=RELU|300,RELU|100,RELU|30,RELU|10,LINEAR");
266  TString layoutString ("Layout=SOFTSIGN|50,SOFTSIGN|30,SOFTSIGN|20,SOFTSIGN|10,LINEAR");
267  TString layoutString ("Layout=TANH|50,TANH|30,TANH|20,TANH|10,LINEAR");
268  TString layoutString ("Layout=SOFTSIGN|50,SOFTSIGN|20,LINEAR");
269  TString layoutString ("Layout=TANH|100,TANH|30,LINEAR");
270  */
271  TString layoutString ("Layout=TANH|100,LINEAR");
272 
273  TString training0 ("LearningRate=1e-5,Momentum=0.5,Repetitions=1,ConvergenceSteps=500,BatchSize=50,TestRepetitions=7,WeightDecay=0.01,Regularization=NONE,DropConfig=0.5+0.5+0.5+0.5,DropRepetitions=2");
274  TString training1 ("LearningRate=1e-5,Momentum=0.9,Repetitions=1,ConvergenceSteps=170,BatchSize=30,TestRepetitions=7,WeightDecay=0.01,Regularization=L2,DropConfig=0.1+0.1+0.1,DropRepetitions=1");
275  TString training2 ("LearningRate=1e-5,Momentum=0.3,Repetitions=1,ConvergenceSteps=150,BatchSize=40,TestRepetitions=7,WeightDecay=0.01,Regularization=NONE");
276  TString training3 ("LearningRate=1e-6,Momentum=0.1,Repetitions=1,ConvergenceSteps=500,BatchSize=100,TestRepetitions=7,WeightDecay=0.0001,Regularization=NONE");
277 
278  TString trainingStrategyString ("TrainingStrategy=");
279  trainingStrategyString += training0 + "|" + training1 + "|" + training2 + "|" + training3;
280 
281 
282  // TString trainingStrategyString ("TrainingStrategy=LearningRate=1e-1,Momentum=0.3,Repetitions=3,ConvergenceSteps=20,BatchSize=30,TestRepetitions=7,WeightDecay=0.0,L1=false,DropFraction=0.0,DropRepetitions=5");
283 
284  TString nnOptions ("!H:V:ErrorStrategy=SUMOFSQUARES:VarTransform=G:WeightInitialization=XAVIERUNIFORM");
285  // TString nnOptions ("!H:V:VarTransform=Normalize:ErrorStrategy=CHECKGRADIENTS");
286  nnOptions.Append (":"); nnOptions.Append (layoutString);
287  nnOptions.Append (":"); nnOptions.Append (trainingStrategyString);
288 
289  factory->BookMethod(dataloader, TMVA::Types::kDNN, "DNN", nnOptions ); // NN
290  }
291 
292 
293 
294  // Support Vector Machine
295  if (Use["SVM"])
296  factory->BookMethod( dataloader, TMVA::Types::kSVM, "SVM", "Gamma=0.25:Tol=0.001:VarTransform=Norm" );
297 
298  // Boosted Decision Trees
299  if (Use["BDT"])
300  factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDT",
301  "!H:!V:NTrees=100:MinNodeSize=1.0%:BoostType=AdaBoostR2:SeparationType=RegressionVariance:nCuts=20:PruneMethod=CostComplexity:PruneStrength=30" );
302 
303  if (Use["BDTG"])
304  factory->BookMethod( dataloader, TMVA::Types::kBDT, "BDTG",
305  "!H:!V:NTrees=2000::BoostType=Grad:Shrinkage=0.1:UseBaggedBoost:BaggedSampleFraction=0.5:nCuts=20:MaxDepth=3:MaxDepth=4" );
306  // --------------------------------------------------------------------------------------------------
307 
308  // Now you can tell the factory to train, test, and evaluate the MVAs
309 
310  // Train MVAs using the set of training events
311  factory->TrainAllMethods();
312 
313  // Evaluate all MVAs using the set of test events
314  factory->TestAllMethods();
315 
316  // Evaluate and compare performance of all configured MVAs
317  factory->EvaluateAllMethods();
318 
319  // --------------------------------------------------------------
320 
321  // Save the output
322  outputFile->Close();
323 
324  std::cout << "==> Wrote root file: " << outputFile->GetName() << std::endl;
325  std::cout << "==> TMVARegression is done!" << std::endl;
326 
327  delete factory;
328  delete dataloader;
329 
330  // Launch the GUI for the root macros
331  if (!gROOT->IsBatch()) TMVA::TMVARegGui( outfileName );
332 }
333 
334 int main( int argc, char** argv )
335 {
336  // Select methods (don't look at this code - not of interest)
337  TString methodList;
338  for (int i=1; i<argc; i++) {
339  TString regMethod(argv[i]);
340  if(regMethod=="-b" || regMethod=="--batch") continue;
341  if (!methodList.IsNull()) methodList += TString(",");
342  methodList += regMethod;
343  }
344  TMVARegression(methodList);
345  return 0;
346 }
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
#define gROOT
Definition: TROOT.h:364
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
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
void AddRegressionTree(TTree *tree, Double_t weight=1.0, Types::ETreeType treetype=Types::kMaxTreeType)
Definition: DataLoader.h:120
void PrepareTrainingAndTestTree(const TCut &cut, const TString &splitOpt)
Definition: DataLoader.cxx:580
double Double_t
Definition: RtypesCore.h:55
void AddTarget(const TString &expression, const TString &title="", const TString &unit="", Double_t min=0, Double_t max=0)
Definition: DataLoader.cxx:472
void SetWeightExpression(const TString &variable, const TString &className="")
Definition: DataLoader.cxx:519
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
void TMVARegGui(const char *fName="TMVAReg.root", TString dataset="")
Definition: TMVARegGui.cxx:46
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