LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
TrainMVA.C File Reference
#include "TFile.h"
#include "TTree.h"
#include "TMVA/Factory.h"
#include <string>
#include <vector>
#include <map>
#include <iostream>
#include "MVAPIDResult.h"

Go to the source code of this file.

Functions

void BuildTree (std::string inFile, std::string outFile)
 
void TrainMVA (std::vector< std::string > signalFiles, std::vector< std::string > backgroundFiles, std::string outputFile, std::string jobName)
 
void PrintRes (std::string inFile)
 

Function Documentation

void BuildTree ( std::string  inFile,
std::string  outFile 
)

Definition at line 12 of file TrainMVA.C.

References anab::MVAPIDResult::concentration, anab::MVAPIDResult::conicalness, anab::MVAPIDResult::coreHaloRatio, anab::MVAPIDResult::dEdxEnd, anab::MVAPIDResult::dEdxEndRatio, anab::MVAPIDResult::dEdxStart, anab::MVAPIDResult::evalRatio, anab::MVAPIDResult::isStoppingReco, anab::MVAPIDResult::isTrack, anab::MVAPIDResult::length, and anab::MVAPIDResult::nSpacePoints.

12  {
13  TFile* fIn=new TFile(inFile.c_str());
14  TTree* tr=(TTree*)fIn->Get("pid/MVAPID");
15  std::vector<anab::MVAPIDResult>* mvares=0;
16  tr->SetBranchAddress("MVAResult",&mvares);
17  TFile* fOut=new TFile(outFile.c_str(),"RECREATE");
18  TTree* mvaTree = new TTree("mvaTree","mvaTree");
19 
20  float evalRatio, concentration, coreHaloRatio, conicalness;
21  float dEdxStart, dEdxEnd, dEdxEndRatio;
22  float length;
23  int isTrack, isStoppingReco;
24 
25  mvaTree->Branch("evalRatio",&evalRatio,"evalRatio/F");
26  mvaTree->Branch("concentration",&concentration,"concentration/F");
27  mvaTree->Branch("coreHaloRatio",&coreHaloRatio,"coreHaloRatio/F");
28  mvaTree->Branch("conicalness",&conicalness,"conicalness/F");
29  mvaTree->Branch("dEdxStart",&dEdxStart,"dEdxStart/F");
30  mvaTree->Branch("dEdxEnd",&dEdxEnd,"dEdxEnd/F");
31  mvaTree->Branch("dEdxEndRatio",&dEdxEndRatio,"dEdxEndRatio/F");
32  mvaTree->Branch("length",&length,"length/F");
33  mvaTree->Branch("isTrack",&isTrack,"isTrack/I");
34  mvaTree->Branch("isStoppingReco",&isStoppingReco,"isStoppingReco/I");
35 
36  for(int iEntry=0;iEntry<tr->GetEntries();++iEntry){
37  tr->GetEntry(iEntry);
38  if(!mvares->size()) continue;
39 
40  anab::MVAPIDResult* biggestTrack=&((*mvares)[0]);
41  for(unsigned int iRes=0;iRes!=mvares->size();++iRes){
42  if((((*mvares)[iRes]).nSpacePoints)>(biggestTrack->nSpacePoints)){
43  biggestTrack=&((*mvares)[iRes]);
44  }
45  }
46 
47  evalRatio=biggestTrack->evalRatio;
48  concentration=biggestTrack->concentration;
49  coreHaloRatio=biggestTrack->coreHaloRatio;
50  conicalness=biggestTrack->conicalness;
51  dEdxStart=biggestTrack->dEdxStart;
52  dEdxEnd=biggestTrack->dEdxEnd;
53  dEdxEndRatio=biggestTrack->dEdxEndRatio;
54  length=biggestTrack->length;
55  isTrack=biggestTrack->isTrack;
56  isStoppingReco=biggestTrack->isStoppingReco;
57 
58  mvaTree->Fill();
59  }
60 
61  fIn->Close();
62  fOut->Write();
63  fOut->Close();
64 }
void PrintRes ( std::string  inFile)

Definition at line 101 of file TrainMVA.C.

References anab::MVAPIDResult::mvaOutput, and anab::MVAPIDResult::nSpacePoints.

101  {
102 
103  TFile* fIn=new TFile(inFile.c_str());
104  TTree* tr=(TTree*)fIn->Get("pid/ANAB");
105  std::vector<anab::MVAPIDResult>* mvares=0;
106  tr->SetBranchAddress("MVAResult",&mvares);
107 
108  for(int iEntry=0;iEntry<tr->GetEntries();++iEntry){
109  tr->GetEntry(iEntry);
110  if(!mvares->size()) continue;
111 
112  anab::MVAPIDResult* biggestTrack=&((*mvares)[0]);
113  for(unsigned int iRes=0;iRes!=mvares->size();++iRes){
114  if((((*mvares)[iRes]).nSpacePoints)>(biggestTrack->nSpacePoints)){
115  biggestTrack=&((*mvares)[iRes]);
116  }
117  }
118 
119  std::cout<<biggestTrack->mvaOutput.at(string("ANN"))<<std::endl;
120  }
121 
122  fIn->Close();
123 }
std::map< std::string, double > mvaOutput
Definition: MVAPIDResult.h:27
void TrainMVA ( std::vector< std::string >  signalFiles,
std::vector< std::string >  backgroundFiles,
std::string  outputFile,
std::string  jobName 
)

Definition at line 66 of file TrainMVA.C.

66  {
67 
68  TFile* fOut = new TFile(outputFile.c_str(),"RECREATE");
69  TMVA::Factory* factory = new TMVA::Factory( jobName.c_str(), fOut, "" );
70 
71  std::vector<TTree*> sigTrees;
72 
73  for(std::vector<std::string>::iterator fIter=signalFiles.begin();fIter!=signalFiles.end();++fIter){
74  TFile* fIn=new TFile(fIter->c_str());
75  factory->AddSignalTree((TTree*)fIn->Get("mvaTree"));
76  }
77 
78  for(std::vector<std::string>::iterator fIter=backgroundFiles.begin();fIter!=backgroundFiles.end();++fIter){
79  TFile* fIn=new TFile(fIter->c_str());
80  factory->AddBackgroundTree((TTree*)fIn->Get("mvaTree"));
81  }
82 
83  factory->AddVariable("evalRatio",'F');
84  factory->AddVariable("concentration",'F');
85  factory->AddVariable("coreHaloRatio",'F');
86  factory->AddVariable("conicalness",'F');
87  factory->AddVariable("dEdxStart",'F');
88  factory->AddVariable("dEdxEnd",'F');
89  factory->AddVariable("dEdxEndRatio",'F');
90 
91  factory->BookMethod( TMVA::Types::kTMlpANN, "ANN", "" );
92  factory->BookMethod( TMVA::Types::kBDT, "BDT", "" );
93  factory->TrainAllMethods();
94  factory->TestAllMethods();
95  factory->EvaluateAllMethods();
96  fOut->Write();
97  fOut->Close();
98 }
intermediate_table::iterator iterator