This macro provides a simple example on how to use the trained regression MVAs within an analysis module
- Project : TMVA - a Root-integrated toolkit for multivariate data analysis
- Package : TMVA
- Exectuable: TMVARegressionApplication
0.037791967392
8.90558886528
Processing /mnt/build/workspace/root-makedoc-v612/rootspi/rdoc/src/v6-12-00-patches/tutorials/tmva/TMVARegressionApplication.C...
==> Start TMVARegressionApplication
: Booking "BDTG method" of type "BDT" from dataset/weights/TMVARegression_BDTG.weights.xml.
<HEADER> DataSetInfo : [Default] : Added class "Regression"
: Booked classifier "BDTG" of type: "BDT"
: Booking "FDA_GA method" of type "FDA" from dataset/weights/TMVARegression_FDA_GA.weights.xml.
: User-defined formula string : "(0)+(1)*x0+(2)*x1"
: TFormula-compatible formula string: "[0]+[1]*[3]+[2]*[4]"
: Booked classifier "FDA_GA" of type: "FDA"
: Booking "KNN method" of type "KNN" from dataset/weights/TMVARegression_KNN.weights.xml.
: Creating kd-tree with 1000 events
: Computing scale factor for 1d distributions: (ifrac, bottom, top) = (80%, 10%, 90%)
<HEADER> ModulekNN : Optimizing tree for 2 variables with 1000 values
: <Fill> Class 1 has 1000 events
: Booked classifier "KNN" of type: "KNN"
: Booking "LD method" of type "LD" from dataset/weights/TMVARegression_LD.weights.xml.
: Booked classifier "LD" of type: "LD"
: Booking "MLP method" of type "MLP" from dataset/weights/TMVARegression_MLP.weights.xml.
<HEADER> MLP : Building Network.
: Initializing weights
: Booked classifier "MLP" of type: "MLP"
: Booking "PDEFoam method" of type "PDEFoam" from dataset/weights/TMVARegression_PDEFoam.weights.xml.
: Read foams from file: dataset/weights/TMVARegression_PDEFoam.weights_foams.root
: Booked classifier "PDEFoam" of type: "PDEFoam"
--- TMVARegressionApp : Using input file: ./files/tmva_reg_example.root
--- Select signal sample
--- Processing: 10000 events
--- ... Processing event: 0
--- ... Processing event: 1000
--- ... Processing event: 2000
--- ... Processing event: 3000
--- ... Processing event: 4000
--- ... Processing event: 5000
--- ... Processing event: 6000
--- ... Processing event: 7000
--- ... Processing event: 8000
--- ... Processing event: 9000
--- End of event loop: Real time 0:00:03, CP time 3.370
--- Created root file: "TMVARegApp.root" containing the MVA output histograms
==> TMVARegressionApplication is done!
#include <cstdlib>
#include <vector>
#include <iostream>
#include <map>
#include <string>
void TMVARegressionApplication( TString myMethodList = "" )
{
std::map<std::string,int> Use;
Use["PDERS"] = 0;
Use["PDEFoam"] = 1;
Use["KNN"] = 1;
Use["LD"] = 1;
Use["FDA_GA"] = 1;
Use["FDA_MC"] = 0;
Use["FDA_MT"] = 0;
Use["FDA_GAMT"] = 0;
Use["MLP"] = 1;
Use["DNN_CPU"] = 0;
Use["SVM"] = 0;
Use["BDT"] = 0;
Use["BDTG"] = 1;
std::cout << std::endl;
std::cout << "==> Start TMVARegressionApplication" << std::endl;
if (myMethodList != "") {
for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) it->second = 0;
for (
UInt_t i=0; i<mlist.size(); i++) {
std::string regMethod(mlist[i]);
if (Use.find(regMethod) == Use.end()) {
std::cout << "Method \"" << regMethod << "\" not known in TMVA under this name. Choose among the following:" << std::endl;
for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) std::cout << it->first << " ";
std::cout << std::endl;
return;
}
Use[regMethod] = 1;
}
}
TString dir = "dataset/weights/";
TString prefix = "TMVARegression";
for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) {
if (it->second) {
TString methodName = it->first + " method";
TString weightfile = dir + prefix + "_" + TString(it->first) + ".weights.xml";
reader->
BookMVA( methodName, weightfile );
}
}
for (std::map<std::string,int>::iterator it = Use.begin(); it != Use.end(); it++) {
TH1* h =
new TH1F( it->first.c_str(), TString(it->first) +
" method", 100, -100, 600 );
if (it->second) hists[++nhists] =
h;
}
nhists++;
TFile *input(0);
TString fname = "./tmva_reg_example.root";
}
else {
input =
TFile::Open(
"http://root.cern.ch/files/tmva_reg_example.root",
"CACHEREAD");
}
if (!input) {
std::cout << "ERROR: could not open data file" << std::endl;
exit(1);
}
std::cout <<
"--- TMVARegressionApp : Using input file: " << input->
GetName() << std::endl;
TTree* theTree = (TTree*)input->
Get(
"TreeR");
std::cout << "--- Select signal sample" << std::endl;
theTree->SetBranchAddress( "var1", &var1 );
theTree->SetBranchAddress( "var2", &var2 );
std::cout << "--- Processing: " << theTree->GetEntries() << " events" << std::endl;
for (
Long64_t ievt=0; ievt<theTree->GetEntries();ievt++) {
if (ievt%1000 == 0) {
std::cout << "--- ... Processing event: " << ievt << std::endl;
}
theTree->GetEntry(ievt);
for (
Int_t ih=0; ih<nhists; ih++) {
}
}
std::cout <<
"--- End of event loop: "; sw.
Print();
TFile *target = new TFile( "TMVARegApp.root","RECREATE" );
for (
Int_t ih=0; ih<nhists; ih++) hists[ih]->Write();
target->Close();
std::cout <<
"--- Created root file: \"" << target->
GetName()
<< "\" containing the MVA output histograms" << std::endl;
std::cout << "==> TMVARegressionApplication is done!" << std::endl << std::endl;
}
int main(
int argc,
char** argv )
{
TString methodList;
for (int i=1; i<argc; i++) {
TString regMethod(argv[i]);
if(regMethod=="-b" || regMethod=="--batch") continue;
if (!methodList.IsNull()) methodList += TString(",");
methodList += regMethod;
}
TMVARegressionApplication(methodList);
return 0;
}
- Author
- Andreas Hoecker
Definition in file TMVARegressionApplication.C.