TMVAClassification.C

Go to the documentation of this file.
00001 // @(#)root/tmva $Id: TMVAClassification.C 37399 2010-12-08 15:22:07Z evt $
00002 /**********************************************************************************
00003  * Project   : TMVA - a ROOT-integrated toolkit for multivariate data analysis    *
00004  * Package   : TMVA                                                               *
00005  * Root Macro: TMVAClassification                                                 *
00006  *                                                                                *
00007  * This macro provides examples for the training and testing of the               *
00008  * TMVA classifiers.                                                              *
00009  *                                                                                *
00010  * As input data is used a toy-MC sample consisting of four Gaussian-distributed  *
00011  * and linearly correlated input variables.                                       *
00012  *                                                                                *
00013  * The methods to be used can be switched on and off by means of booleans, or     *
00014  * via the prompt command, for example:                                           *
00015  *                                                                                *
00016  *    root -l ./TMVAClassification.C\(\"Fisher,Likelihood\"\)                     *
00017  *                                                                                *
00018  * (note that the backslashes are mandatory)                                      *
00019  * If no method given, a default set of classifiers is used.                      *
00020  *                                                                                *
00021  * The output file "TMVA.root" can be analysed with the use of dedicated          *
00022  * macros (simply say: root -l <macro.C>), which can be conveniently              *
00023  * invoked through a GUI that will appear at the end of the run of this macro.    *
00024  * Launch the GUI via the command:                                                *
00025  *                                                                                *
00026  *    root -l ./TMVAGui.C                                                         *
00027  *                                                                                *
00028  **********************************************************************************/
00029 
00030 #include <cstdlib>
00031 #include <iostream>
00032 #include <map>
00033 #include <string>
00034 
00035 #include "TChain.h"
00036 #include "TFile.h"
00037 #include "TTree.h"
00038 #include "TString.h"
00039 #include "TObjString.h"
00040 #include "TSystem.h"
00041 #include "TROOT.h"
00042 
00043 #include "TMVAGui.C"
00044 
00045 #if not defined(__CINT__) || defined(__MAKECINT__)
00046 // needs to be included when makecint runs (ACLIC)
00047 #include "TMVA/Factory.h"
00048 #include "TMVA/Tools.h"
00049 #endif
00050 
00051 void TMVAClassification( TString myMethodList = "" )
00052 {
00053    // The explicit loading of the shared libTMVA is done in TMVAlogon.C, defined in .rootrc
00054    // if you use your private .rootrc, or run from a different directory, please copy the
00055    // corresponding lines from .rootrc
00056 
00057    // methods to be processed can be given as an argument; use format:
00058    //
00059    // mylinux~> root -l TMVAClassification.C\(\"myMethod1,myMethod2,myMethod3\"\)
00060    //
00061    // if you like to use a method via the plugin mechanism, we recommend using
00062    //
00063    // mylinux~> root -l TMVAClassification.C\(\"P_myMethod\"\)
00064    // (an example is given for using the BDT as plugin (see below),
00065    // but of course the real application is when you write your own
00066    // method based)
00067 
00068    //---------------------------------------------------------------
00069    // This loads the library
00070    TMVA::Tools::Instance();
00071 
00072    // Default MVA methods to be trained + tested
00073    std::map<std::string,int> Use;
00074 
00075    // --- Cut optimisation
00076    Use["Cuts"]            = 1;
00077    Use["CutsD"]           = 1;
00078    Use["CutsPCA"]         = 0;
00079    Use["CutsGA"]          = 0;
00080    Use["CutsSA"]          = 0;
00081    // 
00082    // --- 1-dimensional likelihood ("naive Bayes estimator")
00083    Use["Likelihood"]      = 1;
00084    Use["LikelihoodD"]     = 0; // the "D" extension indicates decorrelated input variables (see option strings)
00085    Use["LikelihoodPCA"]   = 1; // the "PCA" extension indicates PCA-transformed input variables (see option strings)
00086    Use["LikelihoodKDE"]   = 0;
00087    Use["LikelihoodMIX"]   = 0;
00088    //
00089    // --- Mutidimensional likelihood and Nearest-Neighbour methods
00090    Use["PDERS"]           = 1;
00091    Use["PDERSD"]          = 0;
00092    Use["PDERSPCA"]        = 0;
00093    Use["PDEFoam"]         = 1;
00094    Use["PDEFoamBoost"]    = 0; // uses generalised MVA method boosting
00095    Use["KNN"]             = 1; // k-nearest neighbour method
00096    //
00097    // --- Linear Discriminant Analysis
00098    Use["LD"]              = 1; // Linear Discriminant identical to Fisher
00099    Use["Fisher"]          = 0;
00100    Use["FisherG"]         = 0;
00101    Use["BoostedFisher"]   = 0; // uses generalised MVA method boosting
00102    Use["HMatrix"]         = 0;
00103    //
00104    // --- Function Discriminant analysis
00105    Use["FDA_GA"]          = 1; // minimisation of user-defined function using Genetics Algorithm
00106    Use["FDA_SA"]          = 0;
00107    Use["FDA_MC"]          = 0;
00108    Use["FDA_MT"]          = 0;
00109    Use["FDA_GAMT"]        = 0;
00110    Use["FDA_MCMT"]        = 0;
00111    //
00112    // --- Neural Networks (all are feed-forward Multilayer Perceptrons)
00113    Use["MLP"]             = 0; // Recommended ANN
00114    Use["MLPBFGS"]         = 0; // Recommended ANN with optional training method
00115    Use["MLPBNN"]          = 1; // Recommended ANN with BFGS training method and bayesian regulator
00116    Use["CFMlpANN"]        = 0; // Depreciated ANN from ALEPH
00117    Use["TMlpANN"]         = 0; // ROOT's own ANN
00118    //
00119    // --- Support Vector Machine 
00120    Use["SVM"]             = 1;
00121    // 
00122    // --- Boosted Decision Trees
00123    Use["BDT"]             = 1; // uses Adaptive Boost
00124    Use["BDTG"]            = 0; // uses Gradient Boost
00125    Use["BDTB"]            = 0; // uses Bagging
00126    Use["BDTD"]            = 0; // decorrelation + Adaptive Boost
00127    // 
00128    // --- Friedman's RuleFit method, ie, an optimised series of cuts ("rules")
00129    Use["RuleFit"]         = 1;
00130    // ---------------------------------------------------------------
00131 
00132    std::cout << std::endl;
00133    std::cout << "==> Start TMVAClassification" << std::endl;
00134 
00135    // Select methods (don't look at this code - not of interest)
00136    if (myMethodList != "") {
00137       for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) it->second = 0;
00138 
00139       std::vector<TString> mlist = TMVA::gTools().SplitString( myMethodList, ',' );
00140       for (UInt_t i=0; i<mlist.size(); i++) {
00141          std::string regMethod(mlist[i]);
00142 
00143          if (Use.find(regMethod) == Use.end()) {
00144             std::cout << "Method \"" << regMethod << "\" not known in TMVA under this name. Choose among the following:" << std::endl;
00145             for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) std::cout << it->first << " ";
00146             std::cout << std::endl;
00147             return;
00148          }
00149          Use[regMethod] = 1;
00150       }
00151    }
00152 
00153    // --------------------------------------------------------------------------------------------------
00154 
00155    // --- Here the preparation phase begins
00156 
00157    // Create a ROOT output file where TMVA will store ntuples, histograms, etc.
00158    TString outfileName( "TMVA.root" );
00159    TFile* outputFile = TFile::Open( outfileName, "RECREATE" );
00160 
00161    // Create the factory object. Later you can choose the methods
00162    // whose performance you'd like to investigate. The factory is 
00163    // the only TMVA object you have to interact with
00164    //
00165    // The first argument is the base of the name of all the
00166    // weightfiles in the directory weight/
00167    //
00168    // The second argument is the output file for the training results
00169    // All TMVA output can be suppressed by removing the "!" (not) in
00170    // front of the "Silent" argument in the option string
00171    TMVA::Factory *factory = new TMVA::Factory( "TMVAClassification", outputFile,
00172                                                "!V:!Silent:Color:DrawProgressBar:Transformations=I;D;P;G,D:AnalysisType=Classification" );
00173 
00174    // If you wish to modify default settings
00175    // (please check "src/Config.h" to see all available global options)
00176    //    (TMVA::gConfig().GetVariablePlotting()).fTimesRMS = 8.0;
00177    //    (TMVA::gConfig().GetIONames()).fWeightFileDir = "myWeightDirectory";
00178 
00179    // Define the input variables that shall be used for the MVA training
00180    // note that you may also use variable expressions, such as: "3*var1/var2*abs(var3)"
00181    // [all types of expressions that can also be parsed by TTree::Draw( "expression" )]
00182    factory->AddVariable( "myvar1 := var1+var2", 'F' );
00183    factory->AddVariable( "myvar2 := var1-var2", "Expression 2", "", 'F' );
00184    factory->AddVariable( "var3",                "Variable 3", "units", 'F' );
00185    factory->AddVariable( "var4",                "Variable 4", "units", 'F' );
00186 
00187    // You can add so-called "Spectator variables", which are not used in the MVA training,
00188    // but will appear in the final "TestTree" produced by TMVA. This TestTree will contain the
00189    // input variables, the response values of all trained MVAs, and the spectator variables
00190    factory->AddSpectator( "spec1 := var1*2",  "Spectator 1", "units", 'F' );
00191    factory->AddSpectator( "spec2 := var1*3",  "Spectator 2", "units", 'F' );
00192 
00193    // Read training and test data
00194    // (it is also possible to use ASCII format as input -> see TMVA Users Guide)
00195    TString fname = "./tmva_class_example.root";
00196    
00197    if (gSystem->AccessPathName( fname ))  // file does not exist in local directory
00198       gSystem->Exec("wget http://root.cern.ch/files/tmva_class_example.root");
00199    
00200    TFile *input = TFile::Open( fname );
00201    
00202    std::cout << "--- TMVAClassification       : Using input file: " << input->GetName() << std::endl;
00203    
00204    // --- Register the training and test trees
00205 
00206    TTree *signal     = (TTree*)input->Get("TreeS");
00207    TTree *background = (TTree*)input->Get("TreeB");
00208    
00209    // global event weights per tree (see below for setting event-wise weights)
00210    Double_t signalWeight     = 1.0;
00211    Double_t backgroundWeight = 1.0;
00212    
00213    // You can add an arbitrary number of signal or background trees
00214    factory->AddSignalTree    ( signal,     signalWeight     );
00215    factory->AddBackgroundTree( background, backgroundWeight );
00216    
00217    // To give different trees for training and testing, do as follows:
00218    //    factory->AddSignalTree( signalTrainingTree, signalTrainWeight, "Training" );
00219    //    factory->AddSignalTree( signalTestTree,     signalTestWeight,  "Test" );
00220    
00221    // Use the following code instead of the above two or four lines to add signal and background
00222    // training and test events "by hand"
00223    // NOTE that in this case one should not give expressions (such as "var1+var2") in the input
00224    //      variable definition, but simply compute the expression before adding the event
00225    //
00226    //     // --- begin ----------------------------------------------------------
00227    //     std::vector<Double_t> vars( 4 ); // vector has size of number of input variables
00228    //     Float_t  treevars[4], weight;
00229    //     
00230    //     // Signal
00231    //     for (UInt_t ivar=0; ivar<4; ivar++) signal->SetBranchAddress( Form( "var%i", ivar+1 ), &(treevars[ivar]) );
00232    //     for (UInt_t i=0; i<signal->GetEntries(); i++) {
00233    //        signal->GetEntry(i);
00234    //        for (UInt_t ivar=0; ivar<4; ivar++) vars[ivar] = treevars[ivar];
00235    //        // add training and test events; here: first half is training, second is testing
00236    //        // note that the weight can also be event-wise
00237    //        if (i < signal->GetEntries()/2.0) factory->AddSignalTrainingEvent( vars, signalWeight );
00238    //        else                              factory->AddSignalTestEvent    ( vars, signalWeight );
00239    //     }
00240    //   
00241    //     // Background (has event weights)
00242    //     background->SetBranchAddress( "weight", &weight );
00243    //     for (UInt_t ivar=0; ivar<4; ivar++) background->SetBranchAddress( Form( "var%i", ivar+1 ), &(treevars[ivar]) );
00244    //     for (UInt_t i=0; i<background->GetEntries(); i++) {
00245    //        background->GetEntry(i);
00246    //        for (UInt_t ivar=0; ivar<4; ivar++) vars[ivar] = treevars[ivar];
00247    //        // add training and test events; here: first half is training, second is testing
00248    //        // note that the weight can also be event-wise
00249    //        if (i < background->GetEntries()/2) factory->AddBackgroundTrainingEvent( vars, backgroundWeight*weight );
00250    //        else                                factory->AddBackgroundTestEvent    ( vars, backgroundWeight*weight );
00251    //     }
00252          // --- end ------------------------------------------------------------
00253    //
00254    // --- end of tree registration 
00255 
00256    // Set individual event weights (the variables must exist in the original TTree)
00257    //    for signal    : factory->SetSignalWeightExpression    ("weight1*weight2");
00258    //    for background: factory->SetBackgroundWeightExpression("weight1*weight2");
00259    factory->SetBackgroundWeightExpression( "weight" );
00260 
00261    // Apply additional cuts on the signal and background samples (can be different)
00262    TCut mycuts = ""; // for example: TCut mycuts = "abs(var1)<0.5 && abs(var2-0.5)<1";
00263    TCut mycutb = ""; // for example: TCut mycutb = "abs(var1)<0.5";
00264 
00265    // Tell the factory how to use the training and testing events
00266    //
00267    // If no numbers of events are given, half of the events in the tree are used 
00268    // for training, and the other half for testing:
00269    //    factory->PrepareTrainingAndTestTree( mycut, "SplitMode=random:!V" );
00270    // To also specify the number of testing events, use:
00271    //    factory->PrepareTrainingAndTestTree( mycut,
00272    //                                         "NSigTrain=3000:NBkgTrain=3000:NSigTest=3000:NBkgTest=3000:SplitMode=Random:!V" );
00273    factory->PrepareTrainingAndTestTree( mycuts, mycutb,
00274                                         "nTrain_Signal=0:nTrain_Background=0:SplitMode=Random:NormMode=NumEvents:!V" );
00275 
00276    // ---- Book MVA methods
00277    //
00278    // Please lookup the various method configuration options in the corresponding cxx files, eg:
00279    // src/MethoCuts.cxx, etc, or here: http://tmva.sourceforge.net/optionRef.html
00280    // it is possible to preset ranges in the option string in which the cut optimisation should be done:
00281    // "...:CutRangeMin[2]=-1:CutRangeMax[2]=1"...", where [2] is the third input variable
00282 
00283    // Cut optimisation
00284    if (Use["Cuts"])
00285       factory->BookMethod( TMVA::Types::kCuts, "Cuts",
00286                            "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart" );
00287 
00288    if (Use["CutsD"])
00289       factory->BookMethod( TMVA::Types::kCuts, "CutsD",
00290                            "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart:VarTransform=Decorrelate" );
00291 
00292    if (Use["CutsPCA"])
00293       factory->BookMethod( TMVA::Types::kCuts, "CutsPCA",
00294                            "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart:VarTransform=PCA" );
00295 
00296    if (Use["CutsGA"])
00297       factory->BookMethod( TMVA::Types::kCuts, "CutsGA",
00298                            "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" );
00299 
00300    if (Use["CutsSA"])
00301       factory->BookMethod( TMVA::Types::kCuts, "CutsSA",
00302                            "!H:!V:FitMethod=SA:EffSel:MaxCalls=150000:KernelTemp=IncAdaptive:InitialTemp=1e+6:MinTemp=1e-6:Eps=1e-10:UseDefaultScale" );
00303 
00304    // Likelihood ("naive Bayes estimator")
00305    if (Use["Likelihood"])
00306       factory->BookMethod( TMVA::Types::kLikelihood, "Likelihood",
00307                            "H:!V:!TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmoothBkg[1]=10:NSmooth=1:NAvEvtPerBin=50" );
00308 
00309    // Decorrelated likelihood
00310    if (Use["LikelihoodD"])
00311       factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodD",
00312                            "!H:!V:TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmooth=5:NAvEvtPerBin=50:VarTransform=Decorrelate" );
00313 
00314    // PCA-transformed likelihood
00315    if (Use["LikelihoodPCA"])
00316       factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodPCA",
00317                            "!H:!V:!TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmooth=5:NAvEvtPerBin=50:VarTransform=PCA" ); 
00318 
00319    // Use a kernel density estimator to approximate the PDFs
00320    if (Use["LikelihoodKDE"])
00321       factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodKDE",
00322                            "!H:!V:!TransformOutput:PDFInterpol=KDE:KDEtype=Gauss:KDEiter=Adaptive:KDEFineFactor=0.3:KDEborder=None:NAvEvtPerBin=50" ); 
00323 
00324    // Use a variable-dependent mix of splines and kernel density estimator
00325    if (Use["LikelihoodMIX"])
00326       factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodMIX",
00327                            "!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" ); 
00328 
00329    // Test the multi-dimensional probability density estimator
00330    // here are the options strings for the MinMax and RMS methods, respectively:
00331    //      "!H:!V:VolumeRangeMode=MinMax:DeltaFrac=0.2:KernelEstimator=Gauss:GaussSigma=0.3" );
00332    //      "!H:!V:VolumeRangeMode=RMS:DeltaFrac=3:KernelEstimator=Gauss:GaussSigma=0.3" );
00333    if (Use["PDERS"])
00334       factory->BookMethod( TMVA::Types::kPDERS, "PDERS",
00335                            "!H:!V:NormTree=T:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600" );
00336 
00337    if (Use["PDERSD"])
00338       factory->BookMethod( TMVA::Types::kPDERS, "PDERSD",
00339                            "!H:!V:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600:VarTransform=Decorrelate" );
00340 
00341    if (Use["PDERSPCA"])
00342       factory->BookMethod( TMVA::Types::kPDERS, "PDERSPCA",
00343                            "!H:!V:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600:VarTransform=PCA" );
00344 
00345    // Multi-dimensional likelihood estimator using self-adapting phase-space binning
00346    if (Use["PDEFoam"])
00347       factory->BookMethod( TMVA::Types::kPDEFoam, "PDEFoam",
00348                            "H:!V:SigBgSeparate=F:TailCut=0.001:VolFrac=0.0333:nActiveCells=500:nSampl=2000:nBin=5:Nmin=100:Kernel=None:Compress=T" );
00349 
00350    if (Use["PDEFoamBoost"])
00351       factory->BookMethod( TMVA::Types::kPDEFoam, "PDEFoamBoost",
00352                            "!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" );
00353 
00354    // K-Nearest Neighbour classifier (KNN)
00355    if (Use["KNN"])
00356       factory->BookMethod( TMVA::Types::kKNN, "KNN",
00357                            "H:nkNN=20:ScaleFrac=0.8:SigmaFact=1.0:Kernel=Gaus:UseKernel=F:UseWeight=T:!Trim" );
00358 
00359    // H-Matrix (chi2-squared) method
00360    if (Use["HMatrix"])
00361       factory->BookMethod( TMVA::Types::kHMatrix, "HMatrix", "!H:!V" );
00362 
00363    // Linear discriminant (same as Fisher discriminant)
00364    if (Use["LD"])
00365       factory->BookMethod( TMVA::Types::kLD, "LD", "H:!V:VarTransform=None:CreateMVAPdfs:PDFInterpolMVAPdf=Spline2:NbinsMVAPdf=50:NsmoothMVAPdf=10" );
00366 
00367    // Fisher discriminant (same as LD)
00368    if (Use["Fisher"])
00369       factory->BookMethod( TMVA::Types::kFisher, "Fisher", "H:!V:Fisher:CreateMVAPdfs:PDFInterpolMVAPdf=Spline2:NbinsMVAPdf=50:NsmoothMVAPdf=10" );
00370 
00371    // Fisher with Gauss-transformed input variables
00372    if (Use["FisherG"])
00373       factory->BookMethod( TMVA::Types::kFisher, "FisherG", "H:!V:VarTransform=Gauss" );
00374 
00375    // Composite classifier: ensemble (tree) of boosted Fisher classifiers
00376    if (Use["BoostedFisher"])
00377       factory->BookMethod( TMVA::Types::kFisher, "BoostedFisher", 
00378                            "H:!V:Boost_Num=20:Boost_Transform=log:Boost_Type=AdaBoost:Boost_AdaBoostBeta=0.2" );
00379 
00380    // Function discrimination analysis (FDA) -- test of various fitters - the recommended one is Minuit (or GA or SA)
00381    if (Use["FDA_MC"])
00382       factory->BookMethod( TMVA::Types::kFDA, "FDA_MC",
00383                            "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" );
00384 
00385    if (Use["FDA_GA"]) // can also use Simulated Annealing (SA) algorithm (see Cuts_SA options])
00386       factory->BookMethod( TMVA::Types::kFDA, "FDA_GA",
00387                            "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" );
00388 
00389    if (Use["FDA_SA"]) // can also use Simulated Annealing (SA) algorithm (see Cuts_SA options])
00390       factory->BookMethod( TMVA::Types::kFDA, "FDA_SA",
00391                            "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" );
00392 
00393    if (Use["FDA_MT"])
00394       factory->BookMethod( TMVA::Types::kFDA, "FDA_MT",
00395                            "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" );
00396 
00397    if (Use["FDA_GAMT"])
00398       factory->BookMethod( TMVA::Types::kFDA, "FDA_GAMT",
00399                            "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" );
00400 
00401    if (Use["FDA_MCMT"])
00402       factory->BookMethod( TMVA::Types::kFDA, "FDA_MCMT",
00403                            "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" );
00404 
00405    // TMVA ANN: MLP (recommended ANN) -- all ANNs in TMVA are Multilayer Perceptrons
00406    if (Use["MLP"])
00407       factory->BookMethod( TMVA::Types::kMLP, "MLP", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:!UseRegulator" );
00408 
00409    if (Use["MLPBFGS"])
00410       factory->BookMethod( TMVA::Types::kMLP, "MLPBFGS", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:TrainingMethod=BFGS:!UseRegulator" );
00411 
00412    if (Use["MLPBNN"])
00413       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
00414 
00415    // CF(Clermont-Ferrand)ANN
00416    if (Use["CFMlpANN"])
00417       factory->BookMethod( TMVA::Types::kCFMlpANN, "CFMlpANN", "!H:!V:NCycles=2000:HiddenLayers=N+1,N"  ); // n_cycles:#nodes:#nodes:...  
00418 
00419    // Tmlp(Root)ANN
00420    if (Use["TMlpANN"])
00421       factory->BookMethod( TMVA::Types::kTMlpANN, "TMlpANN", "!H:!V:NCycles=200:HiddenLayers=N+1,N:LearningMethod=BFGS:ValidationFraction=0.3"  ); // n_cycles:#nodes:#nodes:...
00422 
00423    // Support Vector Machine
00424    if (Use["SVM"])
00425       factory->BookMethod( TMVA::Types::kSVM, "SVM", "Gamma=0.25:Tol=0.001:VarTransform=Norm" );
00426 
00427    // Boosted Decision Trees
00428    if (Use["BDTG"]) // Gradient Boost
00429       factory->BookMethod( TMVA::Types::kBDT, "BDTG",
00430                            "!H:!V:NTrees=1000:BoostType=Grad:Shrinkage=0.10:UseBaggedGrad:GradBaggingFraction=0.5:nCuts=20:NNodesMax=5" );
00431 
00432    if (Use["BDT"])  // Adaptive Boost
00433       factory->BookMethod( TMVA::Types::kBDT, "BDT",
00434                            "!H:!V:NTrees=850:nEventsMin=150:MaxDepth=3:BoostType=AdaBoost:AdaBoostBeta=0.5:SeparationType=GiniIndex:nCuts=20:PruneMethod=NoPruning" );
00435 
00436    if (Use["BDTB"]) // Bagging
00437       factory->BookMethod( TMVA::Types::kBDT, "BDTB",
00438                            "!H:!V:NTrees=400:BoostType=Bagging:SeparationType=GiniIndex:nCuts=20:PruneMethod=NoPruning" );
00439 
00440    if (Use["BDTD"]) // Decorrelation + Adaptive Boost
00441       factory->BookMethod( TMVA::Types::kBDT, "BDTD",
00442                            "!H:!V:NTrees=400:nEventsMin=400:MaxDepth=3:BoostType=AdaBoost:SeparationType=GiniIndex:nCuts=20:PruneMethod=NoPruning:VarTransform=Decorrelate" );
00443 
00444    // RuleFit -- TMVA implementation of Friedman's method
00445    if (Use["RuleFit"])
00446       factory->BookMethod( TMVA::Types::kRuleFit, "RuleFit",
00447                            "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" );
00448 
00449    // For an example of the category classifier usage, see: TMVAClassificationCategory
00450 
00451    // --------------------------------------------------------------------------------------------------
00452 
00453    // ---- Now you can optimize the setting (configuration) of the MVAs using the set of training events
00454 
00455    // factory->OptimizeAllMethods("SigEffAt001","Scan");
00456    // factory->OptimizeAllMethods("ROCIntegral","GA");
00457 
00458    // --------------------------------------------------------------------------------------------------
00459 
00460    // ---- Now you can tell the factory to train, test, and evaluate the MVAs
00461 
00462    // Train MVAs using the set of training events
00463    factory->TrainAllMethods();
00464 
00465    // ---- Evaluate all MVAs using the set of test events
00466    factory->TestAllMethods();
00467 
00468    // ----- Evaluate and compare performance of all configured MVAs
00469    factory->EvaluateAllMethods();
00470 
00471    // --------------------------------------------------------------
00472 
00473    // Save the output
00474    outputFile->Close();
00475 
00476    std::cout << "==> Wrote root file: " << outputFile->GetName() << std::endl;
00477    std::cout << "==> TMVAClassification is done!" << std::endl;
00478 
00479    delete factory;
00480 
00481    // Launch the GUI for the root macros
00482    if (!gROOT->IsBatch()) TMVAGui( outfileName );
00483 }

Generated on Tue Jul 5 15:26:36 2011 for ROOT_528-00b_version by  doxygen 1.5.1