00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include <cstdlib>
00027 #include <iostream>
00028 #include <map>
00029 #include <string>
00030
00031 #include "TChain.h"
00032 #include "TFile.h"
00033 #include "TTree.h"
00034 #include "TString.h"
00035 #include "TObjString.h"
00036 #include "TSystem.h"
00037 #include "TROOT.h"
00038
00039 #include "TMVARegGui.C"
00040
00041 #if not defined(__CINT__) || defined(__MAKECINT__)
00042 #include "TMVA/Tools.h"
00043 #include "TMVA/Factory.h"
00044 #endif
00045
00046 using namespace TMVA;
00047
00048 void TMVARegression( TString myMethodList = "" )
00049 {
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061 TMVA::Tools::Instance();
00062
00063
00064 std::map<std::string,int> Use;
00065
00066
00067 Use["PDERS"] = 0;
00068 Use["PDEFoam"] = 1;
00069 Use["KNN"] = 1;
00070
00071
00072 Use["LD"] = 1;
00073
00074
00075 Use["FDA_GA"] = 1;
00076 Use["FDA_MC"] = 0;
00077 Use["FDA_MT"] = 0;
00078 Use["FDA_GAMT"] = 0;
00079
00080
00081 Use["MLP"] = 1;
00082
00083
00084 Use["SVM"] = 0;
00085
00086
00087 Use["BDT"] = 0;
00088 Use["BDTG"] = 1;
00089
00090
00091 std::cout << std::endl;
00092 std::cout << "==> Start TMVARegression" << std::endl;
00093
00094
00095 if (myMethodList != "") {
00096 for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) it->second = 0;
00097
00098 std::vector<TString> mlist = gTools().SplitString( myMethodList, ',' );
00099 for (UInt_t i=0; i<mlist.size(); i++) {
00100 std::string regMethod(mlist[i]);
00101
00102 if (Use.find(regMethod) == Use.end()) {
00103 std::cout << "Method \"" << regMethod << "\" not known in TMVA under this name. Choose among the following:" << std::endl;
00104 for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) std::cout << it->first << " ";
00105 std::cout << std::endl;
00106 return;
00107 }
00108 Use[regMethod] = 1;
00109 }
00110 }
00111
00112
00113
00114
00115
00116
00117 TString outfileName( "TMVAReg.root" );
00118 TFile* outputFile = TFile::Open( outfileName, "RECREATE" );
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130 TMVA::Factory *factory = new TMVA::Factory( "TMVARegression", outputFile,
00131 "!V:!Silent:Color:DrawProgressBar" );
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141 factory->AddVariable( "var1", "Variable 1", "units", 'F' );
00142 factory->AddVariable( "var2", "Variable 2", "units", 'F' );
00143
00144
00145
00146
00147 factory->AddSpectator( "spec1:=var1*2", "Spectator 1", "units", 'F' );
00148 factory->AddSpectator( "spec2:=var1*3", "Spectator 2", "units", 'F' );
00149
00150
00151 factory->AddTarget( "fvalue" );
00152
00153
00154
00155
00156
00157
00158
00159 TFile *input(0);
00160 TString fname = "./tmva_reg_example.root";
00161 if (!gSystem->AccessPathName( fname ))
00162 input = TFile::Open( fname );
00163 else
00164 input = TFile::Open( "http://root.cern.ch/files/tmva_reg_example.root" );
00165
00166 if (!input) {
00167 std::cout << "ERROR: could not open data file" << std::endl;
00168 exit(1);
00169 }
00170 std::cout << "--- TMVARegression : Using input file: " << input->GetName() << std::endl;
00171
00172
00173
00174 TTree *regTree = (TTree*)input->Get("TreeR");
00175
00176
00177 Double_t regWeight = 1.0;
00178
00179
00180 factory->AddRegressionTree( regTree, regWeight );
00181
00182
00183
00184 factory->SetWeightExpression( "var1", "Regression" );
00185
00186
00187 TCut mycut = "";
00188
00189
00190 factory->PrepareTrainingAndTestTree( mycut,
00191 "nTrain_Regression=0:nTest_Regression=0:SplitMode=Random:NormMode=NumEvents:!V" );
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205 if (Use["PDERS"])
00206 factory->BookMethod( TMVA::Types::kPDERS, "PDERS",
00207 "!H:!V:NormTree=T:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=40:NEventsMax=60:VarTransform=None" );
00208
00209
00210
00211
00212 if (Use["PDEFoam"])
00213 factory->BookMethod( TMVA::Types::kPDEFoam, "PDEFoam",
00214 "!H:!V:MultiTargetRegression=F:TargetSelection=Mpv:TailCut=0.001:VolFrac=0.0333:nActiveCells=500:nSampl=2000:nBin=5:Compress=T:Kernel=None:Nmin=10:VarTransform=None" );
00215
00216
00217 if (Use["KNN"])
00218 factory->BookMethod( TMVA::Types::kKNN, "KNN",
00219 "nkNN=20:ScaleFrac=0.8:SigmaFact=1.0:Kernel=Gaus:UseKernel=F:UseWeight=T:!Trim" );
00220
00221
00222 if (Use["LD"])
00223 factory->BookMethod( TMVA::Types::kLD, "LD",
00224 "!H:!V:VarTransform=None" );
00225
00226
00227 if (Use["FDA_MC"])
00228 factory->BookMethod( TMVA::Types::kFDA, "FDA_MC",
00229 "!H:!V:Formula=(0)+(1)*x0+(2)*x1:ParRanges=(-100,100);(-100,100);(-100,100):FitMethod=MC:SampleSize=100000:Sigma=0.1:VarTransform=D" );
00230
00231 if (Use["FDA_GA"])
00232 factory->BookMethod( TMVA::Types::kFDA, "FDA_GA",
00233 "!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" );
00234
00235 if (Use["FDA_MT"])
00236 factory->BookMethod( TMVA::Types::kFDA, "FDA_MT",
00237 "!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" );
00238
00239 if (Use["FDA_GAMT"])
00240 factory->BookMethod( TMVA::Types::kFDA, "FDA_GAMT",
00241 "!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" );
00242
00243
00244 if (Use["MLP"])
00245 factory->BookMethod( 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" );
00246
00247
00248 if (Use["SVM"])
00249 factory->BookMethod( TMVA::Types::kSVM, "SVM", "Gamma=0.25:Tol=0.001:VarTransform=Norm" );
00250
00251
00252 if (Use["BDT"])
00253 factory->BookMethod( TMVA::Types::kBDT, "BDT",
00254 "!H:!V:NTrees=100:nEventsMin=5:BoostType=AdaBoostR2:SeparationType=RegressionVariance:nCuts=20:PruneMethod=CostComplexity:PruneStrength=30" );
00255
00256 if (Use["BDTG"])
00257 factory->BookMethod( TMVA::Types::kBDT, "BDTG",
00258 "!H:!V:NTrees=2000::BoostType=Grad:Shrinkage=0.1:UseBaggedGrad:GradBaggingFraction=0.5:nCuts=20:MaxDepth=3:NNodesMax=15" );
00259
00260
00261
00262
00263
00264 factory->TrainAllMethods();
00265
00266
00267 factory->TestAllMethods();
00268
00269
00270 factory->EvaluateAllMethods();
00271
00272
00273
00274
00275 outputFile->Close();
00276
00277 std::cout << "==> Wrote root file: " << outputFile->GetName() << std::endl;
00278 std::cout << "==> TMVARegression is done!" << std::endl;
00279
00280 delete factory;
00281
00282
00283 if (!gROOT->IsBatch()) TMVARegGui( outfileName );
00284 }